Actual source code: sieveVec.c
1: #define PETSCVEC_DLL
2: /*
3: Implementation of PETSc Vec using Sieve fields
4: */
5: #include <Mesh.hh>
6: #include private/vecimpl.h
7: #include src/inline/axpy.h
8: #include src/inline/setval.h
9: #include petscblaslapack.h
13: PetscErrorCode VecSet_Sieve(Vec v, PetscScalar alpha)
14: {
15: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) v->data;
16: PetscInt n = v->map.n;
17: PetscScalar *xx = (PetscScalar *) field->restrict(*field->getPatches()->begin());
18: PetscErrorCode ierr;
21: if (alpha == 0.0) {
22: PetscMemzero(xx, n*sizeof(PetscScalar));
23: } else {
24: SET(xx, n, alpha);
25: }
26: return(0);
27: }
31: PetscErrorCode VecScale_Sieve(Vec v, PetscScalar alpha)
32: {
33: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) v->data;
34: PetscBLASInt bn = (PetscBLASInt) v->map.n, one = 1;
35: PetscErrorCode ierr;
38: if (alpha == 0.0) {
39: VecSet_Sieve(v, alpha);
40: } else if (alpha != 1.0) {
41: PetscScalar a = alpha;
43: BLASscal_(&bn, &a, (PetscScalar *) field->restrict(*field->getPatches()->begin()), &one);
44: PetscLogFlops(v->map.n);
45: }
46: return(0);
47: }
51: PetscErrorCode VecCopy_Sieve(Vec x, Vec y)
52: {
53: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) x->data;
54: PetscScalar *yy;
55: PetscErrorCode ierr;
58: if (x != y) {
59: VecGetArray(y, &yy);
60: PetscMemcpy(yy, (PetscScalar *) field->restrict(*field->getPatches()->begin()), x->map.n*sizeof(PetscScalar));
61: VecRestoreArray(y, &yy);
62: }
63: return(0);
64: }
68: PetscErrorCode VecAXPY_Sieve(Vec y, PetscScalar alpha, Vec x)
69: {
70: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) y->data;
71: PetscBLASInt bn = (PetscBLASInt) y->map.n, one = 1;
72: PetscScalar *xarray;
73: PetscErrorCode ierr;
76: /* assume that the BLAS handles alpha == 1.0 efficiently since we have no fast code for it */
77: if (alpha != 0.0) {
78: PetscScalar oalpha = alpha;
80: VecGetArray(x, &xarray);
81: BLASaxpy_(&bn, &oalpha, xarray, &one, (PetscScalar *) field->restrict(*field->getPatches()->begin()), &one);
82: VecRestoreArray(x, &xarray);
83: PetscLogFlops(2*y->map.n);
84: }
85: return(0);
86: }
90: PetscErrorCode VecAYPX_Sieve(Vec y, PetscScalar alpha, Vec x)
91: {
92: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) y->data;
93: PetscScalar *yy = (PetscScalar *) field->restrict(*field->getPatches()->begin());
94: PetscInt n = y->map.n;
95: PetscScalar *xx;
96: PetscErrorCode ierr;
99: if (alpha == 0.0) {
100: VecCopy_Sieve(x, y);
101: } else if (alpha == 1.0) {
102: VecAXPY_Sieve(y, alpha, x);
103: } else {
104: VecGetArray(x, &xx);
105: #if defined(PETSC_USE_FORTRAN_KERNEL_AYPX)
106: {
107: PetscScalar oalpha = alpha;
108: fortranaypx_(&n, &oalpha, xx, yy);
109: }
110: #else
111: {
112: PetscInt i;
113: for (i=0; i<n; i++) {
114: yy[i] = xx[i] + alpha*yy[i];
115: }
116: }
117: #endif
118: VecRestoreArray(x, &xx);
119: PetscLogFlops(2*n);
120: }
121: return(0);
122: }
126: PetscErrorCode VecAXPBY_Sieve(Vec y, PetscScalar alpha, PetscScalar beta, Vec x)
127: {
128: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) y->data;
129: PetscScalar *yy = (PetscScalar *) field->restrict(*field->getPatches()->begin());
130: PetscScalar *xx ,a = alpha,b = beta;
131: PetscInt n = y->map.n, i;
132: PetscErrorCode ierr;
135: if (a == 0.0) {
136: VecScale_Sieve(y, beta);
137: } else if (b == 1.0) {
138: VecAXPY_Sieve(y, alpha, x);
139: } else if (a == 1.0) {
140: VecAYPX_Sieve(y, beta, x);
141: } else if (b == 0.0) {
142: VecGetArray(x, &xx);
143: for(i = 0; i < n; i++) {
144: yy[i] = a*xx[i];
145: }
146: VecRestoreArray(x, &xx);
147: PetscLogFlops(x->map.n);
148: } else {
149: VecGetArray(x, &xx);
150: for(i = 0; i < n; i++) {
151: yy[i] = a*xx[i] + b*yy[i];
152: }
153: VecRestoreArray(x, &xx);
154: PetscLogFlops(3*x->map.n);
155: }
156: return(0);
157: }
161: PetscErrorCode VecMAXPY_Sieve(Vec x, PetscInt nv, const PetscScalar *alpha, Vec *y)
162: {
163: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) x->data;
164: PetscInt n = x->map.n,j,j_rem;
165: PetscScalar *xx = (PetscScalar *) field->restrict(*field->getPatches()->begin());
166: PetscScalar *yy0,*yy1,*yy2,*yy3,alpha0,alpha1,alpha2,alpha3;
167: PetscErrorCode ierr;
169: #if defined(PETSC_HAVE_PRAGMA_DISJOINT)
170: #pragma disjoint(*xx,*yy0,*yy1,*yy2,*yy3,*alpha)
171: #endif
174: PetscLogFlops(nv*2*n);
176: switch (j_rem=nv&0x3) {
177: case 3:
178: VecGetArray(y[0],&yy0);
179: VecGetArray(y[1],&yy1);
180: VecGetArray(y[2],&yy2);
181: alpha0 = alpha[0];
182: alpha1 = alpha[1];
183: alpha2 = alpha[2];
184: alpha += 3;
185: APXY3(xx,alpha0,alpha1,alpha2,yy0,yy1,yy2,n);
186: VecRestoreArray(y[0],&yy0);
187: VecRestoreArray(y[1],&yy1);
188: VecRestoreArray(y[2],&yy2);
189: y += 3;
190: break;
191: case 2:
192: VecGetArray(y[0],&yy0);
193: VecGetArray(y[1],&yy1);
194: alpha0 = alpha[0];
195: alpha1 = alpha[1];
196: alpha +=2;
197: APXY2(xx,alpha0,alpha1,yy0,yy1,n);
198: VecRestoreArray(y[0],&yy0);
199: VecRestoreArray(y[1],&yy1);
200: y +=2;
201: break;
202: case 1:
203: VecGetArray(y[0],&yy0);
204: alpha0 = *alpha++;
205: {PetscBLASInt nn = (PetscBLASInt)n; APXY(xx,alpha0,yy0,nn);}
206: VecRestoreArray(y[0],&yy0);
207: y +=1;
208: break;
209: }
210: for (j=j_rem; j<nv; j+=4) {
211: VecGetArray(y[0],&yy0);
212: VecGetArray(y[1],&yy1);
213: VecGetArray(y[2],&yy2);
214: VecGetArray(y[3],&yy3);
215: alpha0 = alpha[0];
216: alpha1 = alpha[1];
217: alpha2 = alpha[2];
218: alpha3 = alpha[3];
219: alpha += 4;
221: APXY4(xx,alpha0,alpha1,alpha2,alpha3,yy0,yy1,yy2,yy3,n);
222: VecRestoreArray(y[0],&yy0);
223: VecRestoreArray(y[1],&yy1);
224: VecRestoreArray(y[2],&yy2);
225: VecRestoreArray(y[3],&yy3);
226: y += 4;
227: }
228: return(0);
229: }
233: PetscErrorCode VecWAXPY_Sieve(Vec w, PetscScalar alpha, Vec x, Vec y)
234: {
235: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) w->data;
236: PetscScalar *ww = (PetscScalar *) field->restrict(*field->getPatches()->begin());
237: PetscInt n = w->map.n, i;
238: PetscScalar *yy, *xx;
239: PetscErrorCode ierr;
242: VecGetArray(y, &yy);
243: VecGetArray(x, &xx);
244: if (alpha == 1.0) {
245: PetscLogFlops(n);
246: /* could call BLAS axpy after call to memcopy, but may be slower */
247: for (i=0; i<n; i++) ww[i] = yy[i] + xx[i];
248: } else if (alpha == -1.0) {
249: PetscLogFlops(n);
250: for (i=0; i<n; i++) ww[i] = yy[i] - xx[i];
251: } else if (alpha == 0.0) {
252: PetscMemcpy(ww,yy,n*sizeof(PetscScalar));
253: } else {
254: PetscScalar oalpha = alpha;
255: #if defined(PETSC_USE_FORTRAN_KERNEL_WAXPY)
256: fortranwaxpy_(&n,&oalpha,xx,yy,ww);
257: #else
258: for (i=0; i<n; i++) ww[i] = yy[i] + oalpha * xx[i];
259: #endif
260: PetscLogFlops(2*n);
261: }
262: VecRestoreArray(y, &yy);
263: VecRestoreArray(x, &xx);
264: return(0);
265: }
269: PetscErrorCode VecPointwiseMult_Sieve(Vec w, Vec x, Vec y)
270: {
271: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) w->data;
272: PetscScalar *ww = (PetscScalar *) field->restrict(*field->getPatches()->begin());
273: PetscInt n = w->map.n, i;
274: PetscScalar *xx, *yy;
275: PetscErrorCode ierr;
278: VecGetArray(x, &xx);
279: if (x != y) {
280: VecGetArray(y, &yy);
281: } else {
282: yy = xx;
283: }
285: if (ww == xx) {
286: for (i=0; i<n; i++) ww[i] *= yy[i];
287: } else if (ww == yy) {
288: for (i=0; i<n; i++) ww[i] *= xx[i];
289: } else {
290: #if defined(PETSC_USE_FORTRAN_KERNEL_XTIMESY)
291: fortranxtimesy_(xx,yy,ww,&n);
292: #else
293: for (i=0; i<n; i++) ww[i] = xx[i] * yy[i];
294: #endif
295: }
296: VecRestoreArray(x, &xx);
297: if (x != y) {
298: VecRestoreArray(y, &yy);
299: }
300: PetscLogFlops(n);
301: return(0);
302: }
306: PetscErrorCode VecPointwiseDivide_Sieve(Vec w, Vec x, Vec y)
307: {
308: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) w->data;
309: PetscScalar *ww = (PetscScalar *) field->restrict(*field->getPatches()->begin());
310: PetscInt n = w->map.n, i;
311: PetscScalar *xx, *yy;
312: PetscErrorCode ierr;
315: VecGetArray(x, &xx);
316: if (x != y) {
317: VecGetArray(y, &yy);
318: } else {
319: yy = xx;
320: }
321: for (i=0; i<n; i++) {
322: ww[i] = xx[i] / yy[i];
323: }
324: VecRestoreArray(x, &xx);
325: if (x != y) {
326: VecRestoreArray(y, &yy);
327: }
328: PetscLogFlops(n);
329: return(0);
330: }
334: PetscErrorCode VecGetArray_Sieve(Vec v, PetscScalar *a[])
335: {
336: ALE::Mesh::field_type *field = (ALE::Mesh::field_type *) v->data;
340: if (v->array_gotten) {
341: SETERRQ(PETSC_ERR_ORDER,"Array has already been gotten for this vector,you may\n\
342: have forgotten a call to VecRestoreArray()");
343: }
344: v->array_gotten = PETSC_TRUE;
345: *a = (PetscScalar *) field->restrict(*field->getPatches()->begin());
346: PetscObjectTakeAccess(v);
347: return(0);
348: }
352: PetscErrorCode VecRestoreArray_Sieve(Vec v, PetscScalar *a[])
353: {
357: if (!v->array_gotten) {
358: SETERRQ(PETSC_ERR_ORDER,"Array has not been gotten for this vector, you may\n\
359: have forgotten a call to VecGetArray()");
360: }
361: v->array_gotten = PETSC_FALSE;
362: if (a) *a = PETSC_NULL; /* now user cannot accidently use it again */
363: PetscObjectGrantAccess(v);
364: return(0);
365: }
367: static struct _VecOps DvOps = {
368: 0, /* 1 */
369: VecDuplicateVecs_Default,
370: VecDestroyVecs_Default,
371: 0,
372: 0,
373: 0,
374: 0,
375: 0,
376: VecScale_Sieve,
377: VecCopy_Sieve, /* 10 */
378: VecSet_Sieve,
379: 0,
380: VecAXPY_Sieve,
381: VecAXPBY_Sieve,
382: VecMAXPY_Sieve,
383: VecAYPX_Sieve,
384: VecWAXPY_Sieve,
385: VecPointwiseMult_Sieve,
386: VecPointwiseDivide_Sieve,
387: 0, /* 20 */
388: 0,
389: 0,
390: VecGetArray_Sieve,
391: 0,
392: 0,
393: VecRestoreArray_Sieve,
394: 0,
395: 0,
396: 0,
397: 0, /* 30 */
398: 0,
399: 0,
400: 0,
401: 0,
402: 0,
403: 0,
404: 0,
405: 0,
406: 0,
407: 0, /* 40 */
408: 0,
409: 0,
410: 0, /* VecViewNative... */
411: 0,
412: 0,
413: 0,
414: 0,
415: 0,
416: 0,
417: 0, /* 50 */
418: 0,
419: 0,
420: 0,
421: 0};
425: PetscErrorCode VecCreate_Sieve_Private(Vec v, ALE::Mesh::field_type *field)
426: {
427: PetscErrorCode ierr;
430: PetscLogObjectMemory(v, sizeof(ALE::Mesh::field_type) + v->map.n*sizeof(double));
431: PetscMemcpy(v->ops, &DvOps, sizeof(DvOps));
432: v->data = (void *) field;
433: v->mapping = PETSC_NULL;
434: v->bmapping = PETSC_NULL;
435: v->petscnative = PETSC_FALSE;
437: if (v->map.bs == -1) v->map.bs = 1;
438: PetscMapInitialize(v->comm, &v->map);
439: v->stash.insertmode = NOT_SET_VALUES;
440:
441: PetscObjectChangeTypeName((PetscObject) v, VECSIEVE);
442: return(0);
443: }
445: /*MC
446: VECSIEVE - VECSIEVE = "sieve" - The parallel vector based upon Sieve fields
448: Options Database Keys:
449: . -vec_type sieve - sets the vector type to VECSIEVE during a call to VecSetFromOptions()
451: Level: beginner
453: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VECSIEVE, VecType
454: M*/
459: PetscErrorCode VecCreate_Sieve(Vec v)
460: {
461: ALE::Mesh::field_type *field = new ALE::Mesh::field_type(v->comm, 0);
465: VecCreate_Sieve_Private(v, field);
466: return(0);
467: }
472: /*@C
473: VecCreateSieve - Creates a parallel vector implemented by the Sieve Field.
475: Collective on MPI_Comm
477: Input Parameter:
478: . field - The Sieve Field
480: Output Parameter:
481: . v - The vector
483: Level: intermediate
485: Concepts: vectors^creating with array
487: .seealso: VecCreate(), VecDuplicate(), VecDuplicateVecs()
488: @*/
489: PetscErrorCode VecCreateSieve(ALE::Mesh::field_type *field, Vec *v)
490: {
494: VecCreate(field->comm(), v);
495: (*v)->map.n = field->getSize();
496: (*v)->map.N = field->getGlobalOffsets()[field->commSize()];
497: VecCreate_Sieve_Private(*v, field);
498: return(0);
499: }