Actual source code: vscat.c
1: #define PETSCVEC_DLL
2: /*
3: Code for creating scatters between vectors. This file
4: includes the code for scattering between sequential vectors and
5: some special cases for parallel scatters.
6: */
8: #include src/vec/is/isimpl.h
9: #include private/vecimpl.h
11: /* Logging support */
12: PetscCookie PETSCVEC_DLLEXPORT VEC_SCATTER_COOKIE = 0;
14: #if defined(PETSC_USE_DEBUG)
15: /*
16: Checks if any indices are less than zero and generates an error
17: */
20: static PetscErrorCode VecScatterCheckIndices_Private(PetscInt nmax,PetscInt n,PetscInt *idx)
21: {
22: PetscInt i;
25: for (i=0; i<n; i++) {
26: if (idx[i] < 0) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Negative index %D at %D location",idx[i],i);
27: if (idx[i] >= nmax) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Index %D at %D location greater than max %D",idx[i],i,nmax);
28: }
29: return(0);
30: }
31: #endif
33: /*
34: This is special scatter code for when the entire parallel vector is
35: copied to each processor.
37: This code was written by Cameron Cooper, Occidental College, Fall 1995,
38: will working at ANL as a SERS student.
39: */
42: PetscErrorCode VecScatterBegin_MPI_ToAll(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
43: {
44: #if defined(PETSC_USE_64BIT_INDICES)
46: SETERRQ(PETSC_ERR_SUP,"Not implemented due to limited MPI support for extremely long gathers");
47: #else
49: PetscInt yy_n,xx_n;
50: PetscScalar *xv,*yv;
53: VecGetLocalSize(y,&yy_n);
54: VecGetLocalSize(x,&xx_n);
55: VecGetArray(y,&yv);
56: if (x != y) {VecGetArray(x,&xv);}
58: if (mode & SCATTER_REVERSE) {
59: PetscScalar *xvt,*xvt2;
60: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
61: PetscInt i;
63: if (addv == INSERT_VALUES) {
64: PetscInt rstart,rend;
65: /*
66: copy the correct part of the local vector into the local storage of
67: the MPI one. Note: this operation only makes sense if all the local
68: vectors have the same values
69: */
70: VecGetOwnershipRange(y,&rstart,&rend);
71: PetscMemcpy(yv,xv+rstart,yy_n*sizeof(PetscScalar));
72: } else {
73: MPI_Comm comm;
74: PetscMPIInt rank;
75: PetscObjectGetComm((PetscObject)y,&comm);
76: MPI_Comm_rank(comm,&rank);
77: if (scat->work1) xvt = scat->work1;
78: else {
79: PetscMalloc(xx_n*sizeof(PetscScalar),&xvt);
80: scat->work1 = xvt;
81: }
82: if (!rank) { /* I am the zeroth processor, values are accumulated here */
83: if (scat->work2) xvt2 = scat->work2;
84: else {
85: PetscMalloc(xx_n*sizeof(PetscScalar),& xvt2);
86: scat->work2 = xvt2;
87: }
88: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,xvt2,scat->count,y->map.range,MPIU_SCALAR,0,ctx->comm);
89: #if defined(PETSC_USE_COMPLEX)
90: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,ctx->comm);
91: #else
92: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,ctx->comm);
93: #endif
94: if (addv == ADD_VALUES) {
95: for (i=0; i<xx_n; i++) {
96: xvt[i] += xvt2[i];
97: }
98: #if !defined(PETSC_USE_COMPLEX)
99: } else if (addv == MAX_VALUES) {
100: for (i=0; i<xx_n; i++) {
101: xvt[i] = PetscMax(xvt[i],xvt2[i]);
102: }
103: #endif
104: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
105: MPI_Scatterv(xvt,scat->count,y->map.range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
106: } else {
107: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,0, 0,0,MPIU_SCALAR,0,ctx->comm);
108: #if defined(PETSC_USE_COMPLEX)
109: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,ctx->comm);
110: #else
111: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,ctx->comm);
112: #endif
113: MPI_Scatterv(0,scat->count,y->map.range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
114: }
115: }
116: } else {
117: PetscScalar *yvt;
118: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
119: PetscInt i;
121: if (addv == INSERT_VALUES) {
122: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,x->map.range,MPIU_SCALAR,ctx->comm);
123: } else {
124: if (scat->work1) yvt = scat->work1;
125: else {
126: PetscMalloc(yy_n*sizeof(PetscScalar),&yvt);
127: scat->work1 = yvt;
128: }
129: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,x->map.range,MPIU_SCALAR,ctx->comm);
130: if (addv == ADD_VALUES){
131: for (i=0; i<yy_n; i++) {
132: yv[i] += yvt[i];
133: }
134: #if !defined(PETSC_USE_COMPLEX)
135: } else if (addv == MAX_VALUES) {
136: for (i=0; i<yy_n; i++) {
137: yv[i] = PetscMax(yv[i],yvt[i]);
138: }
139: #endif
140: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
141: }
142: }
143: VecRestoreArray(y,&yv);
144: if (x != y) {VecRestoreArray(x,&xv);}
145: #endif
146: return(0);
147: }
149: /*
150: This is special scatter code for when the entire parallel vector is
151: copied to processor 0.
153: */
156: PetscErrorCode VecScatterBegin_MPI_ToOne(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
157: {
158: #if defined(PETSC_USE_64BIT_INDICES)
160: SETERRQ(PETSC_ERR_SUP,"Not implemented due to limited MPI support for extremely long gathers");
161: #else
163: PetscMPIInt rank;
164: PetscInt yy_n,xx_n;
165: PetscScalar *xv,*yv;
166: MPI_Comm comm;
169: VecGetLocalSize(y,&yy_n);
170: VecGetLocalSize(x,&xx_n);
171: VecGetArray(x,&xv);
172: VecGetArray(y,&yv);
174: PetscObjectGetComm((PetscObject)x,&comm);
175: MPI_Comm_rank(comm,&rank);
177: /* -------- Reverse scatter; spread from processor 0 to other processors */
178: if (mode & SCATTER_REVERSE) {
179: PetscScalar *yvt;
180: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
181: PetscInt i;
183: if (addv == INSERT_VALUES) {
184: MPI_Scatterv(xv,scat->count,y->map.range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
185: } else {
186: if (scat->work2) yvt = scat->work2;
187: else {
188: PetscMalloc(xx_n*sizeof(PetscScalar),&yvt);
189: scat->work2 = yvt;
190: }
191: MPI_Scatterv(xv,scat->count,y->map.range,MPIU_SCALAR,yvt,yy_n,MPIU_SCALAR,0,ctx->comm);
192: if (addv == ADD_VALUES) {
193: for (i=0; i<yy_n; i++) {
194: yv[i] += yvt[i];
195: }
196: #if !defined(PETSC_USE_COMPLEX)
197: } else if (addv == MAX_VALUES) {
198: for (i=0; i<yy_n; i++) {
199: yv[i] = PetscMax(yv[i],yvt[i]);
200: }
201: #endif
202: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
203: }
204: /* --------- Forward scatter; gather all values onto processor 0 */
205: } else {
206: PetscScalar *yvt = 0;
207: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
208: PetscInt i;
210: if (addv == INSERT_VALUES) {
211: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,x->map.range,MPIU_SCALAR,0,ctx->comm);
212: } else {
213: if (!rank) {
214: if (scat->work1) yvt = scat->work1;
215: else {
216: PetscMalloc(yy_n*sizeof(PetscScalar),&yvt);
217: scat->work1 = yvt;
218: }
219: }
220: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,x->map.range,MPIU_SCALAR,0,ctx->comm);
221: if (!rank) {
222: if (addv == ADD_VALUES) {
223: for (i=0; i<yy_n; i++) {
224: yv[i] += yvt[i];
225: }
226: #if !defined(PETSC_USE_COMPLEX)
227: } else if (addv == MAX_VALUES) {
228: for (i=0; i<yy_n; i++) {
229: yv[i] = PetscMax(yv[i],yvt[i]);
230: }
231: #endif
232: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
233: }
234: }
235: }
236: VecRestoreArray(x,&xv);
237: VecRestoreArray(y,&yv);
238: #endif
239: return(0);
240: }
242: /*
243: The follow to are used for both VecScatterBegin_MPI_ToAll() and VecScatterBegin_MPI_ToOne()
244: */
247: PetscErrorCode VecScatterDestroy_MPI_ToAll(VecScatter ctx)
248: {
249: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
250: PetscErrorCode ierr;
253: PetscFree(scat->work1);
254: PetscFree(scat->work2);
255: PetscFree2(ctx->todata,scat->count);
256: PetscHeaderDestroy(ctx);
257: return(0);
258: }
262: PetscErrorCode VecScatterDestroy_SGtoSG(VecScatter ctx)
263: {
264: PetscErrorCode ierr;
267: PetscFree4(ctx->todata,((VecScatter_Seq_General*)ctx->todata)->vslots,ctx->fromdata,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
268: PetscHeaderDestroy(ctx);
269: return(0);
270: }
274: PetscErrorCode VecScatterDestroy_SGtoSS(VecScatter ctx)
275: {
276: PetscErrorCode ierr;
279: PetscFree3(ctx->todata,ctx->fromdata,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
280: PetscHeaderDestroy(ctx);
281: return(0);
282: }
286: PetscErrorCode VecScatterDestroy_SStoSG(VecScatter ctx)
287: {
288: PetscErrorCode ierr;
291: PetscFree3(ctx->todata,((VecScatter_Seq_General*)ctx->todata)->vslots,ctx->fromdata);
292: PetscHeaderDestroy(ctx);
293: return(0);
294: }
298: PetscErrorCode VecScatterDestroy_SStoSS(VecScatter ctx)
299: {
303: PetscFree2(ctx->todata,ctx->fromdata);
304: PetscHeaderDestroy(ctx);
305: return(0);
306: }
308: /* -------------------------------------------------------------------------*/
311: PetscErrorCode VecScatterCopy_MPI_ToAll(VecScatter in,VecScatter out)
312: {
313: VecScatter_MPI_ToAll *in_to = (VecScatter_MPI_ToAll*)in->todata,*sto;
314: PetscErrorCode ierr;
315: PetscMPIInt size;
316: PetscInt i;
319: out->postrecvs = 0;
320: out->begin = in->begin;
321: out->end = in->end;
322: out->copy = in->copy;
323: out->destroy = in->destroy;
324: out->view = in->view;
326: MPI_Comm_size(out->comm,&size);
327: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&sto->count);
328: sto->type = in_to->type;
329: for (i=0; i<size; i++) {
330: sto->count[i] = in_to->count[i];
331: }
332: sto->work1 = 0;
333: sto->work2 = 0;
334: out->todata = (void*)sto;
335: out->fromdata = (void*)0;
336: return(0);
337: }
339: /* --------------------------------------------------------------------------------------*/
340: /*
341: Scatter: sequential general to sequential general
342: */
345: PetscErrorCode VecScatterBegin_SGtoSG(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
346: {
347: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
348: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
349: PetscErrorCode ierr;
350: PetscInt i,n = gen_from->n,*fslots,*tslots;
351: PetscScalar *xv,*yv;
352:
354: VecGetArray(x,&xv);
355: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
356: if (mode & SCATTER_REVERSE){
357: gen_to = (VecScatter_Seq_General*)ctx->fromdata;
358: gen_from = (VecScatter_Seq_General*)ctx->todata;
359: }
360: fslots = gen_from->vslots;
361: tslots = gen_to->vslots;
363: if (addv == INSERT_VALUES) {
364: for (i=0; i<n; i++) {yv[tslots[i]] = xv[fslots[i]];}
365: } else if (addv == ADD_VALUES) {
366: for (i=0; i<n; i++) {yv[tslots[i]] += xv[fslots[i]];}
367: #if !defined(PETSC_USE_COMPLEX)
368: } else if (addv == MAX_VALUES) {
369: for (i=0; i<n; i++) {yv[tslots[i]] = PetscMax(yv[tslots[i]],xv[fslots[i]]);}
370: #endif
371: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
372: VecRestoreArray(x,&xv);
373: if (x != y) {VecRestoreArray(y,&yv);}
374: return(0);
375: }
377: /*
378: Scatter: sequential general to sequential stride 1
379: */
382: PetscErrorCode VecScatterBegin_SGtoSS_Stride1(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
383: {
384: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
385: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
386: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
387: PetscErrorCode ierr;
388: PetscInt first = gen_to->first;
389: PetscScalar *xv,*yv;
390:
392: VecGetArray(x,&xv);
393: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
394: if (mode & SCATTER_REVERSE){
395: xv += first;
396: if (addv == INSERT_VALUES) {
397: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
398: } else if (addv == ADD_VALUES) {
399: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
400: #if !defined(PETSC_USE_COMPLEX)
401: } else if (addv == MAX_VALUES) {
402: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
403: #endif
404: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
405: } else {
406: yv += first;
407: if (addv == INSERT_VALUES) {
408: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
409: } else if (addv == ADD_VALUES) {
410: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
411: #if !defined(PETSC_USE_COMPLEX)
412: } else if (addv == MAX_VALUES) {
413: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
414: #endif
415: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
416: }
417: VecRestoreArray(x,&xv);
418: if (x != y) {VecRestoreArray(y,&yv);}
419: return(0);
420: }
422: /*
423: Scatter: sequential general to sequential stride
424: */
427: PetscErrorCode VecScatterBegin_SGtoSS(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
428: {
429: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
430: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
431: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
432: PetscErrorCode ierr;
433: PetscInt first = gen_to->first,step = gen_to->step;
434: PetscScalar *xv,*yv;
435:
437: VecGetArray(x,&xv);
438: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
440: if (mode & SCATTER_REVERSE){
441: if (addv == INSERT_VALUES) {
442: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
443: } else if (addv == ADD_VALUES) {
444: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
445: #if !defined(PETSC_USE_COMPLEX)
446: } else if (addv == MAX_VALUES) {
447: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
448: #endif
449: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
450: } else {
451: if (addv == INSERT_VALUES) {
452: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
453: } else if (addv == ADD_VALUES) {
454: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
455: #if !defined(PETSC_USE_COMPLEX)
456: } else if (addv == MAX_VALUES) {
457: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
458: #endif
459: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
460: }
461: VecRestoreArray(x,&xv);
462: if (x != y) {VecRestoreArray(y,&yv);}
463: return(0);
464: }
466: /*
467: Scatter: sequential stride 1 to sequential general
468: */
471: PetscErrorCode VecScatterBegin_SStoSG_Stride1(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
472: {
473: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
474: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
475: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
476: PetscErrorCode ierr;
477: PetscInt first = gen_from->first;
478: PetscScalar *xv,*yv;
479:
481: VecGetArray(x,&xv);
482: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
484: if (mode & SCATTER_REVERSE){
485: yv += first;
486: if (addv == INSERT_VALUES) {
487: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
488: } else if (addv == ADD_VALUES) {
489: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
490: #if !defined(PETSC_USE_COMPLEX)
491: } else if (addv == MAX_VALUES) {
492: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
493: #endif
494: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
495: } else {
496: xv += first;
497: if (addv == INSERT_VALUES) {
498: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
499: } else if (addv == ADD_VALUES) {
500: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
501: #if !defined(PETSC_USE_COMPLEX)
502: } else if (addv == MAX_VALUES) {
503: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
504: #endif
505: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
506: }
507: VecRestoreArray(x,&xv);
508: if (x != y) {VecRestoreArray(y,&yv);}
509: return(0);
510: }
514: /*
515: Scatter: sequential stride to sequential general
516: */
517: PetscErrorCode VecScatterBegin_SStoSG(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
518: {
519: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
520: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
521: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
522: PetscErrorCode ierr;
523: PetscInt first = gen_from->first,step = gen_from->step;
524: PetscScalar *xv,*yv;
525:
527: VecGetArray(x,&xv);
528: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
530: if (mode & SCATTER_REVERSE){
531: if (addv == INSERT_VALUES) {
532: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
533: } else if (addv == ADD_VALUES) {
534: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
535: #if !defined(PETSC_USE_COMPLEX)
536: } else if (addv == MAX_VALUES) {
537: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
538: #endif
539: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
540: } else {
541: if (addv == INSERT_VALUES) {
542: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
543: } else if (addv == ADD_VALUES) {
544: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
545: #if !defined(PETSC_USE_COMPLEX)
546: } else if (addv == MAX_VALUES) {
547: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
548: #endif
549: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
550: }
551: VecRestoreArray(x,&xv);
552: if (x != y) {VecRestoreArray(y,&yv);}
553: return(0);
554: }
556: /*
557: Scatter: sequential stride to sequential stride
558: */
561: PetscErrorCode VecScatterBegin_SStoSS(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
562: {
563: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
564: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
565: PetscInt i,n = gen_from->n,to_first = gen_to->first,to_step = gen_to->step;
566: PetscErrorCode ierr;
567: PetscInt from_first = gen_from->first,from_step = gen_from->step;
568: PetscScalar *xv,*yv;
569:
571: VecGetArray(x,&xv);
572: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
574: if (mode & SCATTER_REVERSE){
575: from_first = gen_to->first;
576: to_first = gen_from->first;
577: from_step = gen_to->step;
578: to_step = gen_from->step;
579: }
581: if (addv == INSERT_VALUES) {
582: if (to_step == 1 && from_step == 1) {
583: PetscMemcpy(yv+to_first,xv+from_first,n*sizeof(PetscScalar));
584: } else {
585: for (i=0; i<n; i++) {
586: yv[to_first + i*to_step] = xv[from_first+i*from_step];
587: }
588: }
589: } else if (addv == ADD_VALUES) {
590: if (to_step == 1 && from_step == 1) {
591: yv += to_first; xv += from_first;
592: for (i=0; i<n; i++) {
593: yv[i] += xv[i];
594: }
595: } else {
596: for (i=0; i<n; i++) {
597: yv[to_first + i*to_step] += xv[from_first+i*from_step];
598: }
599: }
600: #if !defined(PETSC_USE_COMPLEX)
601: } else if (addv == MAX_VALUES) {
602: if (to_step == 1 && from_step == 1) {
603: yv += to_first; xv += from_first;
604: for (i=0; i<n; i++) {
605: yv[i] = PetscMax(yv[i],xv[i]);
606: }
607: } else {
608: for (i=0; i<n; i++) {
609: yv[to_first + i*to_step] = PetscMax(yv[to_first + i*to_step],xv[from_first+i*from_step]);
610: }
611: }
612: #endif
613: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
614: VecRestoreArray(x,&xv);
615: if (x != y) {VecRestoreArray(y,&yv);}
616: return(0);
617: }
619: /* --------------------------------------------------------------------------*/
624: PetscErrorCode VecScatterCopy_SGToSG(VecScatter in,VecScatter out)
625: {
626: PetscErrorCode ierr;
627: VecScatter_Seq_General *in_to = (VecScatter_Seq_General*)in->todata,*out_to;
628: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from;
629:
631: out->postrecvs = 0;
632: out->begin = in->begin;
633: out->end = in->end;
634: out->copy = in->copy;
635: out->destroy = in->destroy;
636: out->view = in->view;
638: PetscMalloc4(1,VecScatter_Seq_General,&out_to,in_to->n,PetscInt,&out_to->vslots,1,VecScatter_Seq_General,&out_from,in_from->n,PetscInt,&out_from->vslots);
639: out_to->n = in_to->n;
640: out_to->type = in_to->type;
641: out_to->nonmatching_computed = PETSC_FALSE;
642: out_to->slots_nonmatching = 0;
643: out_to->is_copy = PETSC_FALSE;
644: PetscMemcpy(out_to->vslots,in_to->vslots,(out_to->n)*sizeof(PetscInt));
646: out_from->n = in_from->n;
647: out_from->type = in_from->type;
648: out_from->nonmatching_computed = PETSC_FALSE;
649: out_from->slots_nonmatching = 0;
650: out_from->is_copy = PETSC_FALSE;
651: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
653: out->todata = (void*)out_to;
654: out->fromdata = (void*)out_from;
655: return(0);
656: }
661: PetscErrorCode VecScatterCopy_SGToSS(VecScatter in,VecScatter out)
662: {
663: PetscErrorCode ierr;
664: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to;
665: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from;
666:
668: out->postrecvs = 0;
669: out->begin = in->begin;
670: out->end = in->end;
671: out->copy = in->copy;
672: out->destroy = in->destroy;
673: out->view = in->view;
675: PetscMalloc3(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_General,&out_from,in_from->n,PetscInt,&out_from->vslots);
676: out_to->n = in_to->n;
677: out_to->type = in_to->type;
678: out_to->first = in_to->first;
679: out_to->step = in_to->step;
680: out_to->type = in_to->type;
682: out_from->n = in_from->n;
683: out_from->type = in_from->type;
684: out_from->nonmatching_computed = PETSC_FALSE;
685: out_from->slots_nonmatching = 0;
686: out_from->is_copy = PETSC_FALSE;
687: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
689: out->todata = (void*)out_to;
690: out->fromdata = (void*)out_from;
691: return(0);
692: }
694: /* --------------------------------------------------------------------------*/
695: /*
696: Scatter: parallel to sequential vector, sequential strides for both.
697: */
700: PetscErrorCode VecScatterCopy_SStoSS(VecScatter in,VecScatter out)
701: {
702: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to;
703: VecScatter_Seq_Stride *in_from = (VecScatter_Seq_Stride*)in->fromdata,*out_from;
704: PetscErrorCode ierr;
707: out->postrecvs = 0;
708: out->begin = in->begin;
709: out->end = in->end;
710: out->copy = in->copy;
711: out->destroy = in->destroy;
712: out->view = in->view;
714: PetscMalloc2(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_Stride,&out_from);
715: out_to->n = in_to->n;
716: out_to->type = in_to->type;
717: out_to->first = in_to->first;
718: out_to->step = in_to->step;
719: out_to->type = in_to->type;
720: out_from->n = in_from->n;
721: out_from->type = in_from->type;
722: out_from->first = in_from->first;
723: out_from->step = in_from->step;
724: out_from->type = in_from->type;
725: out->todata = (void*)out_to;
726: out->fromdata = (void*)out_from;
727: return(0);
728: }
730: EXTERN PetscErrorCode VecScatterCreate_PtoS(PetscInt,PetscInt *,PetscInt,PetscInt *,Vec,Vec,PetscInt,VecScatter);
731: EXTERN PetscErrorCode VecScatterCreate_PtoP(PetscInt,PetscInt *,PetscInt,PetscInt *,Vec,Vec,VecScatter);
732: EXTERN PetscErrorCode VecScatterCreate_StoP(PetscInt,PetscInt *,PetscInt,PetscInt *,Vec,PetscInt,VecScatter);
734: /* =======================================================================*/
735: #define VEC_SEQ_ID 0
736: #define VEC_MPI_ID 1
738: /*
739: Blocksizes we have optimized scatters for
740: */
742: #define VecScatterOptimizedBS(mbs) ((2 <= mbs && mbs <= 8) || mbs == 12)
746: /*@C
747: VecScatterCreate - Creates a vector scatter context.
749: Collective on Vec
751: Input Parameters:
752: + xin - a vector that defines the shape (parallel data layout of the vector)
753: of vectors from which we scatter
754: . yin - a vector that defines the shape (parallel data layout of the vector)
755: of vectors to which we scatter
756: . ix - the indices of xin to scatter (if PETSC_NULL scatters all values)
757: - iy - the indices of yin to hold results (if PETSC_NULL fills entire vector yin)
759: Output Parameter:
760: . newctx - location to store the new scatter context
762: Options Database Keys:
763: + -vecscatter_merge - VecScatterBegin() handles all of the communication, VecScatterEnd() is a nop
764: eliminates the chance for overlap of computation and communication
765: . -vecscatter_ssend - Uses MPI_Ssend_init() instead of MPI_Send_init()
766: . -vecscatter_sendfirst - Posts sends before receives
767: . -vecscatter_rr - use ready receiver mode for MPI sends
768: - -vecscatter_packtogether - Pack all messages before sending, receive all messages before unpacking
769: ONLY implemented for BLOCK SIZE of 4 and 12! (others easily added)
771: Level: intermediate
773: Notes:
774: In calls to VecScatter() you can use different vectors than the xin and
775: yin you used above; BUT they must have the same parallel data layout, for example,
776: they could be obtained from VecDuplicate().
777: A VecScatter context CANNOT be used in two or more simultaneous scatters;
778: that is you cannot call a second VecScatterBegin() with the same scatter
779: context until the VecScatterEnd() has been called on the first VecScatterBegin().
780: In this case a separate VecScatter is needed for each concurrent scatter.
782: Concepts: scatter^between vectors
783: Concepts: gather^between vectors
785: .seealso: VecScatterDestroy(), VecScatterCreateToAll(), VecScatterCreateToZero()
786: @*/
787: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterCreate(Vec xin,IS ix,Vec yin,IS iy,VecScatter *newctx)
788: {
789: VecScatter ctx;
791: PetscMPIInt size;
792: PetscInt totalv,xin_type = VEC_SEQ_ID,yin_type = VEC_SEQ_ID,*range;
793: MPI_Comm comm,ycomm;
794: PetscTruth ixblock,iyblock,iystride,islocal,cando,flag;
795: IS tix = 0,tiy = 0;
799: /*
800: Determine if the vectors are "parallel", ie. it shares a comm with other processors, or
801: sequential (it does not share a comm). The difference is that parallel vectors treat the
802: index set as providing indices in the global parallel numbering of the vector, with
803: sequential vectors treat the index set as providing indices in the local sequential
804: numbering
805: */
806: PetscObjectGetComm((PetscObject)xin,&comm);
807: MPI_Comm_size(comm,&size);
808: if (size > 1) {xin_type = VEC_MPI_ID;}
810: PetscObjectGetComm((PetscObject)yin,&ycomm);
811: MPI_Comm_size(ycomm,&size);
812: if (size > 1) {comm = ycomm; yin_type = VEC_MPI_ID;}
813:
814: /* generate the Scatter context */
815: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
816: ctx->inuse = PETSC_FALSE;
818: ctx->beginandendtogether = PETSC_FALSE;
819: PetscOptionsHasName(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether);
820: if (ctx->beginandendtogether) {
821: PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
822: }
823: PetscOptionsHasName(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether);
824: if (ctx->packtogether) {
825: PetscInfo(ctx,"Pack all messages before sending\n");
826: }
828: VecGetLocalSize(xin,&ctx->from_n);
829: VecGetLocalSize(yin,&ctx->to_n);
831: /*
832: if ix or iy is not included; assume just grabbing entire vector
833: */
834: if (!ix && xin_type == VEC_SEQ_ID) {
835: ISCreateStride(comm,ctx->from_n,0,1,&ix);
836: tix = ix;
837: } else if (!ix && xin_type == VEC_MPI_ID) {
838: if (yin_type == VEC_MPI_ID) {
839: PetscInt ntmp, low;
840: VecGetLocalSize(xin,&ntmp);
841: VecGetOwnershipRange(xin,&low,PETSC_NULL);
842: ISCreateStride(comm,ntmp,low,1,&ix);
843: } else{
844: PetscInt Ntmp;
845: VecGetSize(xin,&Ntmp);
846: ISCreateStride(comm,Ntmp,0,1,&ix);
847: }
848: tix = ix;
849: } else if (!ix) {
850: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
851: }
853: if (!iy && yin_type == VEC_SEQ_ID) {
854: ISCreateStride(comm,ctx->to_n,0,1,&iy);
855: tiy = iy;
856: } else if (!iy && yin_type == VEC_MPI_ID) {
857: if (xin_type == VEC_MPI_ID) {
858: PetscInt ntmp, low;
859: VecGetLocalSize(yin,&ntmp);
860: VecGetOwnershipRange(yin,&low,PETSC_NULL);
861: ISCreateStride(comm,ntmp,low,1,&iy);
862: } else{
863: PetscInt Ntmp;
864: VecGetSize(yin,&Ntmp);
865: ISCreateStride(comm,Ntmp,0,1,&iy);
866: }
867: tiy = iy;
868: } else if (!iy) {
869: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
870: }
872: /* ===========================================================================================================
873: Check for special cases
874: ==========================================================================================================*/
875: /* ---------------------------------------------------------------------------*/
876: if (xin_type == VEC_SEQ_ID && yin_type == VEC_SEQ_ID) {
877: if (ix->type == IS_GENERAL && iy->type == IS_GENERAL){
878: PetscInt nx,ny,*idx,*idy;
879: VecScatter_Seq_General *to,*from;
881: ISGetLocalSize(ix,&nx);
882: ISGetLocalSize(iy,&ny);
883: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
884: ISGetIndices(ix,&idx);
885: ISGetIndices(iy,&idy);
886: PetscMalloc4(1,VecScatter_Seq_General,&to,nx,PetscInt,&to->vslots,1,VecScatter_Seq_General,&from,nx,PetscInt,&from->vslots);
887: to->n = nx;
888: #if defined(PETSC_USE_DEBUG)
889: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
890: #endif
891: PetscMemcpy(to->vslots,idy,nx*sizeof(PetscInt));
892: from->n = nx;
893: #if defined(PETSC_USE_DEBUG)
894: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
895: #endif
896: PetscMemcpy(from->vslots,idx,nx*sizeof(PetscInt));
897: to->type = VEC_SCATTER_SEQ_GENERAL;
898: from->type = VEC_SCATTER_SEQ_GENERAL;
899: ctx->todata = (void*)to;
900: ctx->fromdata = (void*)from;
901: ctx->postrecvs = 0;
902: ctx->begin = VecScatterBegin_SGtoSG;
903: ctx->end = 0;
904: ctx->destroy = VecScatterDestroy_SGtoSG;
905: ctx->copy = VecScatterCopy_SGToSG;
906: PetscInfo(xin,"Special case: sequential vector general scatter\n");
907: goto functionend;
908: } else if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
909: PetscInt nx,ny,to_first,to_step,from_first,from_step;
910: VecScatter_Seq_Stride *from8,*to8;
912: ISGetLocalSize(ix,&nx);
913: ISGetLocalSize(iy,&ny);
914: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
915: ISStrideGetInfo(iy,&to_first,&to_step);
916: ISStrideGetInfo(ix,&from_first,&from_step);
917: PetscMalloc2(1,VecScatter_Seq_Stride,&to8,1,VecScatter_Seq_Stride,&from8);
918: to8->n = nx;
919: to8->first = to_first;
920: to8->step = to_step;
921: from8->n = nx;
922: from8->first = from_first;
923: from8->step = from_step;
924: to8->type = VEC_SCATTER_SEQ_STRIDE;
925: from8->type = VEC_SCATTER_SEQ_STRIDE;
926: ctx->todata = (void*)to8;
927: ctx->fromdata = (void*)from8;
928: ctx->postrecvs = 0;
929: ctx->begin = VecScatterBegin_SStoSS;
930: ctx->end = 0;
931: ctx->destroy = VecScatterDestroy_SStoSS;
932: ctx->copy = VecScatterCopy_SStoSS;
933: PetscInfo(xin,"Special case: sequential vector stride to stride\n");
934: goto functionend;
935: } else if (ix->type == IS_GENERAL && iy->type == IS_STRIDE){
936: PetscInt nx,ny,*idx,first,step;
937: VecScatter_Seq_General *from9;
938: VecScatter_Seq_Stride *to9;
940: ISGetLocalSize(ix,&nx);
941: ISGetIndices(ix,&idx);
942: ISGetLocalSize(iy,&ny);
943: ISStrideGetInfo(iy,&first,&step);
944: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
945: PetscMalloc3(1,VecScatter_Seq_Stride,&to9,1,VecScatter_Seq_General,&from9,nx,PetscInt,&from9->vslots);
946: to9->n = nx;
947: to9->first = first;
948: to9->step = step;
949: from9->n = nx;
950: #if defined(PETSC_USE_DEBUG)
951: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
952: #endif
953: PetscMemcpy(from9->vslots,idx,nx*sizeof(PetscInt));
954: ctx->todata = (void*)to9; ctx->fromdata = (void*)from9;
955: ctx->postrecvs = 0;
956: if (step == 1) ctx->begin = VecScatterBegin_SGtoSS_Stride1;
957: else ctx->begin = VecScatterBegin_SGtoSS;
958: ctx->destroy = VecScatterDestroy_SGtoSS;
959: ctx->end = 0;
960: ctx->copy = VecScatterCopy_SGToSS;
961: to9->type = VEC_SCATTER_SEQ_STRIDE;
962: from9->type = VEC_SCATTER_SEQ_GENERAL;
963: PetscInfo(xin,"Special case: sequential vector general to stride\n");
964: goto functionend;
965: } else if (ix->type == IS_STRIDE && iy->type == IS_GENERAL){
966: PetscInt nx,ny,*idy,first,step;
967: VecScatter_Seq_General *to10;
968: VecScatter_Seq_Stride *from10;
970: ISGetLocalSize(ix,&nx);
971: ISGetIndices(iy,&idy);
972: ISGetLocalSize(iy,&ny);
973: ISStrideGetInfo(ix,&first,&step);
974: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
975: PetscMalloc3(1,VecScatter_Seq_General,&to10,nx,PetscInt,&to10->vslots,1,VecScatter_Seq_Stride,&from10);
976: from10->n = nx;
977: from10->first = first;
978: from10->step = step;
979: to10->n = nx;
980: #if defined(PETSC_USE_DEBUG)
981: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
982: #endif
983: PetscMemcpy(to10->vslots,idy,nx*sizeof(PetscInt));
984: ctx->todata = (void*)to10;
985: ctx->fromdata = (void*)from10;
986: ctx->postrecvs = 0;
987: if (step == 1) ctx->begin = VecScatterBegin_SStoSG_Stride1;
988: else ctx->begin = VecScatterBegin_SStoSG;
989: ctx->destroy = VecScatterDestroy_SStoSG;
990: ctx->end = 0;
991: ctx->copy = 0;
992: to10->type = VEC_SCATTER_SEQ_GENERAL;
993: from10->type = VEC_SCATTER_SEQ_STRIDE;
994: PetscInfo(xin,"Special case: sequential vector stride to general\n");
995: goto functionend;
996: } else {
997: PetscInt nx,ny,*idx,*idy;
998: VecScatter_Seq_General *to11,*from11;
999: PetscTruth idnx,idny;
1001: ISGetLocalSize(ix,&nx);
1002: ISGetLocalSize(iy,&ny);
1003: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1005: ISIdentity(ix,&idnx);
1006: ISIdentity(iy,&idny);
1007: if (idnx && idny) {
1008: VecScatter_Seq_Stride *to13,*from13;
1009: PetscMalloc2(1,VecScatter_Seq_Stride,&to13,1,VecScatter_Seq_Stride,&from13);
1010: to13->n = nx;
1011: to13->first = 0;
1012: to13->step = 1;
1013: from13->n = nx;
1014: from13->first = 0;
1015: from13->step = 1;
1016: to13->type = VEC_SCATTER_SEQ_STRIDE;
1017: from13->type = VEC_SCATTER_SEQ_STRIDE;
1018: ctx->todata = (void*)to13;
1019: ctx->fromdata = (void*)from13;
1020: ctx->postrecvs = 0;
1021: ctx->begin = VecScatterBegin_SStoSS;
1022: ctx->end = 0;
1023: ctx->destroy = VecScatterDestroy_SStoSS;
1024: ctx->copy = VecScatterCopy_SStoSS;
1025: PetscInfo(xin,"Special case: sequential copy\n");
1026: goto functionend;
1027: }
1029: ISGetIndices(iy,&idy);
1030: ISGetIndices(ix,&idx);
1031: PetscMalloc4(1,VecScatter_Seq_General,&to11,nx,PetscInt,&to11->vslots,1,VecScatter_Seq_General,&from11,nx,PetscInt,&from11->vslots);
1032: to11->n = nx;
1033: #if defined(PETSC_USE_DEBUG)
1034: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1035: #endif
1036: PetscMemcpy(to11->vslots,idy,nx*sizeof(PetscInt));
1037: from11->n = nx;
1038: #if defined(PETSC_USE_DEBUG)
1039: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
1040: #endif
1041: PetscMemcpy(from11->vslots,idx,nx*sizeof(PetscInt));
1042: to11->type = VEC_SCATTER_SEQ_GENERAL;
1043: from11->type = VEC_SCATTER_SEQ_GENERAL;
1044: ctx->todata = (void*)to11;
1045: ctx->fromdata = (void*)from11;
1046: ctx->postrecvs = 0;
1047: ctx->begin = VecScatterBegin_SGtoSG;
1048: ctx->end = 0;
1049: ctx->destroy = VecScatterDestroy_SGtoSG;
1050: ctx->copy = VecScatterCopy_SGToSG;
1051: ISRestoreIndices(ix,&idx);
1052: ISRestoreIndices(iy,&idy);
1053: PetscInfo(xin,"Sequential vector scatter with block indices\n");
1054: goto functionend;
1055: }
1056: }
1057: /* ---------------------------------------------------------------------------*/
1058: if (xin_type == VEC_MPI_ID && yin_type == VEC_SEQ_ID) {
1060: /* ===========================================================================================================
1061: Check for special cases
1062: ==========================================================================================================*/
1063: islocal = PETSC_FALSE;
1064: /* special case extracting (subset of) local portion */
1065: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1066: PetscInt nx,ny,to_first,to_step,from_first,from_step;
1067: PetscInt start,end;
1068: VecScatter_Seq_Stride *from12,*to12;
1070: VecGetOwnershipRange(xin,&start,&end);
1071: ISGetLocalSize(ix,&nx);
1072: ISStrideGetInfo(ix,&from_first,&from_step);
1073: ISGetLocalSize(iy,&ny);
1074: ISStrideGetInfo(iy,&to_first,&to_step);
1075: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1076: if (ix->min >= start && ix->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1077: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1078: if (cando) {
1079: PetscMalloc2(1,VecScatter_Seq_Stride,&to12,1,VecScatter_Seq_Stride,&from12);
1080: to12->n = nx;
1081: to12->first = to_first;
1082: to12->step = to_step;
1083: from12->n = nx;
1084: from12->first = from_first-start;
1085: from12->step = from_step;
1086: to12->type = VEC_SCATTER_SEQ_STRIDE;
1087: from12->type = VEC_SCATTER_SEQ_STRIDE;
1088: ctx->todata = (void*)to12;
1089: ctx->fromdata = (void*)from12;
1090: ctx->postrecvs = 0;
1091: ctx->begin = VecScatterBegin_SStoSS;
1092: ctx->end = 0;
1093: ctx->destroy = VecScatterDestroy_SStoSS;
1094: ctx->copy = VecScatterCopy_SStoSS;
1095: PetscInfo(xin,"Special case: processors only getting local values\n");
1096: goto functionend;
1097: }
1098: } else {
1099: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1100: }
1102: /* test for special case of all processors getting entire vector */
1103: totalv = 0;
1104: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1105: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1106: PetscMPIInt *count;
1107: VecScatter_MPI_ToAll *sto;
1109: ISGetLocalSize(ix,&nx);
1110: ISStrideGetInfo(ix,&from_first,&from_step);
1111: ISGetLocalSize(iy,&ny);
1112: ISStrideGetInfo(iy,&to_first,&to_step);
1113: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1114: VecGetSize(xin,&N);
1115: if (nx != N) {
1116: totalv = 0;
1117: } else if (from_first == 0 && from_step == 1 &&
1118: from_first == to_first && from_step == to_step){
1119: totalv = 1;
1120: } else totalv = 0;
1121: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1123: if (cando) {
1125: MPI_Comm_size(ctx->comm,&size);
1126: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1127: range = xin->map.range;
1128: for (i=0; i<size; i++) {
1129: count[i] = range[i+1] - range[i];
1130: }
1131: sto->count = count;
1132: sto->work1 = 0;
1133: sto->work2 = 0;
1134: sto->type = VEC_SCATTER_MPI_TOALL;
1135: ctx->todata = (void*)sto;
1136: ctx->fromdata = 0;
1137: ctx->postrecvs = 0;
1138: ctx->begin = VecScatterBegin_MPI_ToAll;
1139: ctx->end = 0;
1140: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1141: ctx->copy = VecScatterCopy_MPI_ToAll;
1142: PetscInfo(xin,"Special case: all processors get entire parallel vector\n");
1143: goto functionend;
1144: }
1145: } else {
1146: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1147: }
1149: /* test for special case of processor 0 getting entire vector */
1150: totalv = 0;
1151: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1152: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1153: PetscMPIInt rank,*count;
1154: VecScatter_MPI_ToAll *sto;
1156: PetscObjectGetComm((PetscObject)xin,&comm);
1157: MPI_Comm_rank(comm,&rank);
1158: ISGetLocalSize(ix,&nx);
1159: ISStrideGetInfo(ix,&from_first,&from_step);
1160: ISGetLocalSize(iy,&ny);
1161: ISStrideGetInfo(iy,&to_first,&to_step);
1162: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1163: if (!rank) {
1164: VecGetSize(xin,&N);
1165: if (nx != N) {
1166: totalv = 0;
1167: } else if (from_first == 0 && from_step == 1 &&
1168: from_first == to_first && from_step == to_step){
1169: totalv = 1;
1170: } else totalv = 0;
1171: } else {
1172: if (!nx) totalv = 1;
1173: else totalv = 0;
1174: }
1175: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1177: if (cando) {
1178: MPI_Comm_size(ctx->comm,&size);
1179: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1180: range = xin->map.range;
1181: for (i=0; i<size; i++) {
1182: count[i] = range[i+1] - range[i];
1183: }
1184: sto->count = count;
1185: sto->work1 = 0;
1186: sto->work2 = 0;
1187: sto->type = VEC_SCATTER_MPI_TOONE;
1188: ctx->todata = (void*)sto;
1189: ctx->fromdata = 0;
1190: ctx->postrecvs = 0;
1191: ctx->begin = VecScatterBegin_MPI_ToOne;
1192: ctx->end = 0;
1193: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1194: ctx->copy = VecScatterCopy_MPI_ToAll;
1195: PetscInfo(xin,"Special case: processor zero gets entire parallel vector, rest get none\n");
1196: goto functionend;
1197: }
1198: } else {
1199: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1200: }
1202: ISBlock(ix,&ixblock);
1203: ISBlock(iy,&iyblock);
1204: ISStride(iy,&iystride);
1205: if (ixblock) {
1206: /* special case block to block */
1207: if (iyblock) {
1208: PetscInt nx,ny,*idx,*idy,bsx,bsy;
1209: ISBlockGetBlockSize(iy,&bsy);
1210: ISBlockGetBlockSize(ix,&bsx);
1211: if (bsx == bsy && VecScatterOptimizedBS(bsx)) {
1212: ISBlockGetSize(ix,&nx);
1213: ISBlockGetIndices(ix,&idx);
1214: ISBlockGetSize(iy,&ny);
1215: ISBlockGetIndices(iy,&idy);
1216: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1217: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,bsx,ctx);
1218: ISBlockRestoreIndices(ix,&idx);
1219: ISBlockRestoreIndices(iy,&idy);
1220: PetscInfo(xin,"Special case: blocked indices\n");
1221: goto functionend;
1222: }
1223: /* special case block to stride */
1224: } else if (iystride) {
1225: PetscInt ystart,ystride,ysize,bsx;
1226: ISStrideGetInfo(iy,&ystart,&ystride);
1227: ISGetLocalSize(iy,&ysize);
1228: ISBlockGetBlockSize(ix,&bsx);
1229: /* see if stride index set is equivalent to block index set */
1230: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1231: PetscInt nx,*idx,*idy,il;
1232: ISBlockGetSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1233: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1234: PetscMalloc(nx*sizeof(PetscInt),&idy);
1235: if (nx) {
1236: idy[0] = ystart;
1237: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1238: }
1239: VecScatterCreate_PtoS(nx,idx,nx,idy,xin,yin,bsx,ctx);
1240: PetscFree(idy);
1241: ISBlockRestoreIndices(ix,&idx);
1242: PetscInfo(xin,"Special case: blocked indices to stride\n");
1243: goto functionend;
1244: }
1245: }
1246: }
1247: /* left over general case */
1248: {
1249: PetscInt nx,ny,*idx,*idy;
1250: ISGetLocalSize(ix,&nx);
1251: ISGetIndices(ix,&idx);
1252: ISGetLocalSize(iy,&ny);
1253: ISGetIndices(iy,&idy);
1254: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1255: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,1,ctx);
1256: ISRestoreIndices(ix,&idx);
1257: ISRestoreIndices(iy,&idy);
1258: PetscInfo(xin,"General case: MPI to Seq\n");
1259: goto functionend;
1260: }
1261: }
1262: /* ---------------------------------------------------------------------------*/
1263: if (xin_type == VEC_SEQ_ID && yin_type == VEC_MPI_ID) {
1264: /* ===========================================================================================================
1265: Check for special cases
1266: ==========================================================================================================*/
1267: /* special case local copy portion */
1268: islocal = PETSC_FALSE;
1269: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1270: PetscInt nx,ny,to_first,to_step,from_step,start,end,from_first;
1271: VecScatter_Seq_Stride *from,*to;
1273: VecGetOwnershipRange(yin,&start,&end);
1274: ISGetLocalSize(ix,&nx);
1275: ISStrideGetInfo(ix,&from_first,&from_step);
1276: ISGetLocalSize(iy,&ny);
1277: ISStrideGetInfo(iy,&to_first,&to_step);
1278: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1279: if (iy->min >= start && iy->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1280: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,yin->comm);
1281: if (cando) {
1282: PetscMalloc2(1,VecScatter_Seq_Stride,&to,1,VecScatter_Seq_Stride,&from);
1283: to->n = nx;
1284: to->first = to_first-start;
1285: to->step = to_step;
1286: from->n = nx;
1287: from->first = from_first;
1288: from->step = from_step;
1289: to->type = VEC_SCATTER_SEQ_STRIDE;
1290: from->type = VEC_SCATTER_SEQ_STRIDE;
1291: ctx->todata = (void*)to;
1292: ctx->fromdata = (void*)from;
1293: ctx->postrecvs = 0;
1294: ctx->begin = VecScatterBegin_SStoSS;
1295: ctx->end = 0;
1296: ctx->destroy = VecScatterDestroy_SStoSS;
1297: ctx->copy = VecScatterCopy_SStoSS;
1298: PetscInfo(xin,"Special case: sequential stride to MPI stride\n");
1299: goto functionend;
1300: }
1301: } else {
1302: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,yin->comm);
1303: }
1304: if (ix->type == IS_BLOCK && iy->type == IS_STRIDE){
1305: PetscInt ystart,ystride,ysize,bsx;
1306: ISStrideGetInfo(iy,&ystart,&ystride);
1307: ISGetLocalSize(iy,&ysize);
1308: ISBlockGetBlockSize(ix,&bsx);
1309: /* see if stride index set is equivalent to block index set */
1310: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1311: PetscInt nx,*idx,*idy,il;
1312: ISBlockGetSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1313: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1314: PetscMalloc(nx*sizeof(PetscInt),&idy);
1315: if (nx) {
1316: idy[0] = ystart;
1317: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1318: }
1319: VecScatterCreate_StoP(nx,idx,nx,idy,yin,bsx,ctx);
1320: PetscFree(idy);
1321: ISBlockRestoreIndices(ix,&idx);
1322: PetscInfo(xin,"Special case: Blocked indices to stride\n");
1323: goto functionend;
1324: }
1325: }
1327: /* general case */
1328: {
1329: PetscInt nx,ny,*idx,*idy;
1330: ISGetLocalSize(ix,&nx);
1331: ISGetIndices(ix,&idx);
1332: ISGetLocalSize(iy,&ny);
1333: ISGetIndices(iy,&idy);
1334: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1335: VecScatterCreate_StoP(nx,idx,ny,idy,yin,1,ctx);
1336: ISRestoreIndices(ix,&idx);
1337: ISRestoreIndices(iy,&idy);
1338: PetscInfo(xin,"General case: Seq to MPI\n");
1339: goto functionend;
1340: }
1341: }
1342: /* ---------------------------------------------------------------------------*/
1343: if (xin_type == VEC_MPI_ID && yin_type == VEC_MPI_ID) {
1344: /* no special cases for now */
1345: PetscInt nx,ny,*idx,*idy;
1346: ISGetLocalSize(ix,&nx);
1347: ISGetIndices(ix,&idx);
1348: ISGetLocalSize(iy,&ny);
1349: ISGetIndices(iy,&idy);
1350: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1351: VecScatterCreate_PtoP(nx,idx,ny,idy,xin,yin,ctx);
1352: ISRestoreIndices(ix,&idx);
1353: ISRestoreIndices(iy,&idy);
1354: PetscInfo(xin,"General case: MPI to MPI\n");
1355: goto functionend;
1356: }
1358: functionend:
1359: *newctx = ctx;
1360: if (tix) {ISDestroy(tix);}
1361: if (tiy) {ISDestroy(tiy);}
1362: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view_info",&flag);
1363: if (flag) {
1364: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm),PETSC_VIEWER_ASCII_INFO);
1365: VecScatterView(ctx,PETSC_VIEWER_STDOUT_(comm));
1366: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm));
1367: }
1368: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view",&flag);
1369: if (flag) {
1370: VecScatterView(ctx,PETSC_VIEWER_STDOUT_(comm));
1371: }
1372: return(0);
1373: }
1375: /* ------------------------------------------------------------------*/
1378: /*@
1379: VecScatterPostRecvs - Posts the receives required for the ready-receiver
1380: version of the VecScatter routines.
1382: Collective on VecScatter
1384: Input Parameters:
1385: + x - the vector from which we scatter (not needed,can be null)
1386: . y - the vector to which we scatter
1387: . addv - either ADD_VALUES or INSERT_VALUES
1388: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1389: SCATTER_FORWARD, SCATTER_REVERSE
1390: - inctx - scatter context generated by VecScatterCreate()
1392: Output Parameter:
1393: . y - the vector to which we scatter
1395: Level: advanced
1397: Notes:
1398: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1399: the SCATTER_FORWARD.
1400: The vectors x and y cannot be the same. y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1402: This scatter is far more general than the conventional
1403: scatter, since it can be a gather or a scatter or a combination,
1404: depending on the indices ix and iy. If x is a parallel vector and y
1405: is sequential, VecScatterBegin() can serve to gather values to a
1406: single processor. Similarly, if y is parallel and x sequential, the
1407: routine can scatter from one processor to many processors.
1409: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1410: @*/
1411: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterPostRecvs(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter inctx)
1412: {
1419: if (inctx->postrecvs) {
1420: (*inctx->postrecvs)(x,y,addv,mode,inctx);
1421: }
1422: return(0);
1423: }
1425: /* ------------------------------------------------------------------*/
1428: /*@
1429: VecScatterGetMerged - Returns true if the scatter is completed in the VecScatterBegin()
1430: and the VecScatterEnd() does nothing
1432: Not Collective
1434: Input Parameter:
1435: . ctx - scatter context created with VecScatterCreate()
1437: Output Parameter:
1438: . flg - PETSC_TRUE if the VecScatterBegin/End() are all done during the VecScatterBegin()
1440: Level: developer
1442: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1443: @*/
1444: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterGetMerged(VecScatter ctx,PetscTruth *flg)
1445: {
1448: *flg = ctx->beginandendtogether;
1449: return(0);
1450: }
1454: /*@
1455: VecScatterBegin - Begins a generalized scatter from one vector to
1456: another. Complete the scattering phase with VecScatterEnd().
1458: Collective on VecScatter and Vec
1460: Input Parameters:
1461: + x - the vector from which we scatter
1462: . y - the vector to which we scatter
1463: . addv - either ADD_VALUES or INSERT_VALUES, with INSERT_VALUES mode any location
1464: not scattered to retains its old value; i.e. the vector is NOT first zeroed.
1465: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1466: SCATTER_FORWARD or SCATTER_REVERSE
1467: - inctx - scatter context generated by VecScatterCreate()
1469: Level: intermediate
1471: Notes:
1472: The vectors x and y need not be the same vectors used in the call
1473: to VecScatterCreate(), but x must have the same parallel data layout
1474: as that passed in as the x to VecScatterCreate(), similarly for the y.
1475: Most likely they have been obtained from VecDuplicate().
1477: You cannot change the values in the input vector between the calls to VecScatterBegin()
1478: and VecScatterEnd().
1480: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1481: the SCATTER_FORWARD.
1482:
1483: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1485: This scatter is far more general than the conventional
1486: scatter, since it can be a gather or a scatter or a combination,
1487: depending on the indices ix and iy. If x is a parallel vector and y
1488: is sequential, VecScatterBegin() can serve to gather values to a
1489: single processor. Similarly, if y is parallel and x sequential, the
1490: routine can scatter from one processor to many processors.
1492: Concepts: scatter^between vectors
1493: Concepts: gather^between vectors
1495: .seealso: VecScatterCreate(), VecScatterEnd()
1496: @*/
1497: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterBegin(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter inctx)
1498: {
1500: #if defined(PETSC_USE_DEBUG)
1501: PetscInt to_n,from_n;
1502: #endif
1508: if (inctx->inuse) SETERRQ(PETSC_ERR_ARG_WRONGSTATE," Scatter ctx already in use");
1509: #if defined(PETSC_USE_DEBUG)
1510: /*
1511: Error checking to make sure these vectors match the vectors used
1512: to create the vector scatter context. -1 in the from_n and to_n indicate the
1513: vector lengths are unknown (for example with mapped scatters) and thus
1514: no error checking is performed.
1515: */
1516: if (inctx->from_n >= 0 && inctx->to_n >= 0) {
1517: VecGetLocalSize(x,&from_n);
1518: VecGetLocalSize(y,&to_n);
1519: if (mode & SCATTER_REVERSE) {
1520: if (to_n != inctx->from_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter (scatter reverse and vector to != ctx from size)");
1521: if (from_n != inctx->to_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter (scatter reverse and vector from != ctx to size)");
1522: } else {
1523: if (to_n != inctx->to_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter (scatter forward and vector to != ctx to size)");
1524: if (from_n != inctx->from_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter (scatter forward and vector from != ctx from size)");
1525: }
1526: }
1527: #endif
1529: inctx->inuse = PETSC_TRUE;
1530: PetscLogEventBarrierBegin(VEC_ScatterBarrier,0,0,0,0,inctx->comm);
1531: (*inctx->begin)(x,y,addv,mode,inctx);
1532: if (inctx->beginandendtogether && inctx->end) {
1533: inctx->inuse = PETSC_FALSE;
1534: (*inctx->end)(x,y,addv,mode,inctx);
1535: }
1536: PetscLogEventBarrierEnd(VEC_ScatterBarrier,0,0,0,0,inctx->comm);
1537: return(0);
1538: }
1540: /* --------------------------------------------------------------------*/
1543: /*@
1544: VecScatterEnd - Ends a generalized scatter from one vector to another. Call
1545: after first calling VecScatterBegin().
1547: Collective on VecScatter and Vec
1549: Input Parameters:
1550: + x - the vector from which we scatter
1551: . y - the vector to which we scatter
1552: . addv - either ADD_VALUES or INSERT_VALUES.
1553: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1554: SCATTER_FORWARD, SCATTER_REVERSE
1555: - ctx - scatter context generated by VecScatterCreate()
1557: Level: intermediate
1559: Notes:
1560: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1561: the SCATTER_FORWARD.
1562: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1564: .seealso: VecScatterBegin(), VecScatterCreate()
1565: @*/
1566: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterEnd(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
1567: {
1574: ctx->inuse = PETSC_FALSE;
1575: if (!ctx->end) return(0);
1576: if (!ctx->beginandendtogether) {
1577: PetscLogEventBegin(VEC_ScatterEnd,ctx,x,y,0);
1578: (*(ctx)->end)(x,y,addv,mode,ctx);
1579: PetscLogEventEnd(VEC_ScatterEnd,ctx,x,y,0);
1580: }
1581: return(0);
1582: }
1586: /*@
1587: VecScatterDestroy - Destroys a scatter context created by
1588: VecScatterCreate().
1590: Collective on VecScatter
1592: Input Parameter:
1593: . ctx - the scatter context
1595: Level: intermediate
1597: .seealso: VecScatterCreate(), VecScatterCopy()
1598: @*/
1599: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterDestroy(VecScatter ctx)
1600: {
1605: if (--ctx->refct > 0) return(0);
1607: /* if memory was published with AMS then destroy it */
1608: PetscObjectDepublish(ctx);
1610: (*ctx->destroy)(ctx);
1611: return(0);
1612: }
1616: /*@
1617: VecScatterCopy - Makes a copy of a scatter context.
1619: Collective on VecScatter
1621: Input Parameter:
1622: . sctx - the scatter context
1624: Output Parameter:
1625: . ctx - the context copy
1627: Level: advanced
1629: .seealso: VecScatterCreate(), VecScatterDestroy()
1630: @*/
1631: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterCopy(VecScatter sctx,VecScatter *ctx)
1632: {
1638: if (!sctx->copy) SETERRQ(PETSC_ERR_SUP,"Cannot copy this type");
1639: PetscHeaderCreate(*ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",sctx->comm,VecScatterDestroy,VecScatterView);
1640: (*ctx)->to_n = sctx->to_n;
1641: (*ctx)->from_n = sctx->from_n;
1642: (*sctx->copy)(sctx,*ctx);
1643: return(0);
1644: }
1647: /* ------------------------------------------------------------------*/
1650: /*@
1651: VecScatterView - Views a vector scatter context.
1653: Collective on VecScatter
1655: Input Parameters:
1656: + ctx - the scatter context
1657: - viewer - the viewer for displaying the context
1659: Level: intermediate
1661: @*/
1662: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterView(VecScatter ctx,PetscViewer viewer)
1663: {
1668: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ctx->comm);
1670: if (!ctx->view) SETERRQ(PETSC_ERR_SUP,"Cannot view this type of scatter context yet");
1672: (*ctx->view)(ctx,viewer);
1673: return(0);
1674: }
1678: /*@
1679: VecScatterRemap - Remaps the "from" and "to" indices in a
1680: vector scatter context. FOR EXPERTS ONLY!
1682: Collective on VecScatter
1684: Input Parameters:
1685: + scat - vector scatter context
1686: . from - remapping for "from" indices (may be PETSC_NULL)
1687: - to - remapping for "to" indices (may be PETSC_NULL)
1689: Level: developer
1691: Notes: In the parallel case the todata is actually the indices
1692: from which the data is TAKEN! The from stuff is where the
1693: data is finally put. This is VERY VERY confusing!
1695: In the sequential case the todata is the indices where the
1696: data is put and the fromdata is where it is taken from.
1697: This is backwards from the paralllel case! CRY! CRY! CRY!
1699: @*/
1700: PetscErrorCode PETSCVEC_DLLEXPORT VecScatterRemap(VecScatter scat,PetscInt *rto,PetscInt *rfrom)
1701: {
1702: VecScatter_Seq_General *to,*from;
1703: VecScatter_MPI_General *mto;
1704: PetscInt i;
1711: from = (VecScatter_Seq_General *)scat->fromdata;
1712: mto = (VecScatter_MPI_General *)scat->todata;
1714: if (mto->type == VEC_SCATTER_MPI_TOALL) SETERRQ(PETSC_ERR_ARG_SIZ,"Not for to all scatter");
1716: if (rto) {
1717: if (mto->type == VEC_SCATTER_MPI_GENERAL) {
1718: /* handle off processor parts */
1719: for (i=0; i<mto->starts[mto->n]; i++) {
1720: mto->indices[i] = rto[mto->indices[i]];
1721: }
1722: /* handle local part */
1723: to = &mto->local;
1724: for (i=0; i<to->n; i++) {
1725: to->vslots[i] = rto[to->vslots[i]];
1726: }
1727: } else if (from->type == VEC_SCATTER_SEQ_GENERAL) {
1728: for (i=0; i<from->n; i++) {
1729: from->vslots[i] = rto[from->vslots[i]];
1730: }
1731: } else if (from->type == VEC_SCATTER_SEQ_STRIDE) {
1732: VecScatter_Seq_Stride *sto = (VecScatter_Seq_Stride*)from;
1733:
1734: /* if the remapping is the identity and stride is identity then skip remap */
1735: if (sto->step == 1 && sto->first == 0) {
1736: for (i=0; i<sto->n; i++) {
1737: if (rto[i] != i) {
1738: SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1739: }
1740: }
1741: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1742: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1743: }
1745: if (rfrom) {
1746: SETERRQ(PETSC_ERR_SUP,"Unable to remap the FROM in scatters yet");
1747: }
1749: /*
1750: Mark then vector lengths as unknown because we do not know the
1751: lengths of the remapped vectors
1752: */
1753: scat->from_n = -1;
1754: scat->to_n = -1;
1756: return(0);
1757: }