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 private/isimpl.h
9: #include private/vecimpl.h
11: /* Logging support */
12: PetscCookie VEC_SCATTER_COOKIE;
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,const 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 copied to each processor.
36: This code was written by Cameron Cooper, Occidental College, Fall 1995,
37: will working at ANL as a SERS student.
38: */
41: PetscErrorCode VecScatterBegin_MPI_ToAll(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
42: {
44: PetscInt yy_n,xx_n;
45: PetscScalar *xv,*yv;
48: VecGetLocalSize(y,&yy_n);
49: VecGetLocalSize(x,&xx_n);
50: VecGetArray(y,&yv);
51: if (x != y) {VecGetArray(x,&xv);}
53: if (mode & SCATTER_REVERSE) {
54: PetscScalar *xvt,*xvt2;
55: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
56: PetscInt i;
57: PetscMPIInt *disply = scat->displx;
59: if (addv == INSERT_VALUES) {
60: PetscInt rstart,rend;
61: /*
62: copy the correct part of the local vector into the local storage of
63: the MPI one. Note: this operation only makes sense if all the local
64: vectors have the same values
65: */
66: VecGetOwnershipRange(y,&rstart,&rend);
67: PetscMemcpy(yv,xv+rstart,yy_n*sizeof(PetscScalar));
68: } else {
69: MPI_Comm comm;
70: PetscMPIInt rank;
71: PetscObjectGetComm((PetscObject)y,&comm);
72: MPI_Comm_rank(comm,&rank);
73: if (scat->work1) xvt = scat->work1;
74: else {
75: PetscMalloc(xx_n*sizeof(PetscScalar),&xvt);
76: scat->work1 = xvt;
77: }
78: if (!rank) { /* I am the zeroth processor, values are accumulated here */
79: if (scat->work2) xvt2 = scat->work2;
80: else {
81: PetscMalloc(xx_n*sizeof(PetscScalar),& xvt2);
82: scat->work2 = xvt2;
83: }
84: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,xvt2,scat->count,disply,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
85: #if defined(PETSC_USE_COMPLEX)
86: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,((PetscObject)ctx)->comm);
87: #else
88: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,((PetscObject)ctx)->comm);
89: #endif
90: if (addv == ADD_VALUES) {
91: for (i=0; i<xx_n; i++) {
92: xvt[i] += xvt2[i];
93: }
94: #if !defined(PETSC_USE_COMPLEX)
95: } else if (addv == MAX_VALUES) {
96: for (i=0; i<xx_n; i++) {
97: xvt[i] = PetscMax(xvt[i],xvt2[i]);
98: }
99: #endif
100: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
101: MPI_Scatterv(xvt,scat->count,disply,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
102: } else {
103: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,0, 0,0,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
104: #if defined(PETSC_USE_COMPLEX)
105: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,((PetscObject)ctx)->comm);
106: #else
107: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,((PetscObject)ctx)->comm);
108: #endif
109: MPI_Scatterv(0,scat->count,disply,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
110: }
111: }
112: } else {
113: PetscScalar *yvt;
114: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
115: PetscInt i;
116: PetscMPIInt *displx = scat->displx;
118: if (addv == INSERT_VALUES) {
119: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,displx,MPIU_SCALAR,((PetscObject)ctx)->comm);
120: } else {
121: if (scat->work1) yvt = scat->work1;
122: else {
123: PetscMalloc(yy_n*sizeof(PetscScalar),&yvt);
124: scat->work1 = yvt;
125: }
126: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,displx,MPIU_SCALAR,((PetscObject)ctx)->comm);
127: if (addv == ADD_VALUES){
128: for (i=0; i<yy_n; i++) {
129: yv[i] += yvt[i];
130: }
131: #if !defined(PETSC_USE_COMPLEX)
132: } else if (addv == MAX_VALUES) {
133: for (i=0; i<yy_n; i++) {
134: yv[i] = PetscMax(yv[i],yvt[i]);
135: }
136: #endif
137: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
138: }
139: }
140: VecRestoreArray(y,&yv);
141: if (x != y) {VecRestoreArray(x,&xv);}
142: return(0);
143: }
145: /*
146: This is special scatter code for when the entire parallel vector is copied to processor 0.
148: */
151: PetscErrorCode VecScatterBegin_MPI_ToOne(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
152: {
154: PetscMPIInt rank;
155: PetscInt yy_n,xx_n;
156: PetscScalar *xv,*yv;
157: MPI_Comm comm;
160: VecGetLocalSize(y,&yy_n);
161: VecGetLocalSize(x,&xx_n);
162: VecGetArray(x,&xv);
163: VecGetArray(y,&yv);
165: PetscObjectGetComm((PetscObject)x,&comm);
166: MPI_Comm_rank(comm,&rank);
168: /* -------- Reverse scatter; spread from processor 0 to other processors */
169: if (mode & SCATTER_REVERSE) {
170: PetscScalar *yvt;
171: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
172: PetscInt i;
173: PetscMPIInt *disply = scat->displx;
175: if (addv == INSERT_VALUES) {
176: MPI_Scatterv(xv,scat->count,disply,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
177: } else {
178: if (scat->work2) yvt = scat->work2;
179: else {
180: PetscMalloc(xx_n*sizeof(PetscScalar),&yvt);
181: scat->work2 = yvt;
182: }
183: MPI_Scatterv(xv,scat->count,disply,MPIU_SCALAR,yvt,yy_n,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
184: if (addv == ADD_VALUES) {
185: for (i=0; i<yy_n; i++) {
186: yv[i] += yvt[i];
187: }
188: #if !defined(PETSC_USE_COMPLEX)
189: } else if (addv == MAX_VALUES) {
190: for (i=0; i<yy_n; i++) {
191: yv[i] = PetscMax(yv[i],yvt[i]);
192: }
193: #endif
194: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
195: }
196: /* --------- Forward scatter; gather all values onto processor 0 */
197: } else {
198: PetscScalar *yvt = 0;
199: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
200: PetscInt i;
201: PetscMPIInt *displx = scat->displx;
203: if (addv == INSERT_VALUES) {
204: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,displx,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
205: } else {
206: if (!rank) {
207: if (scat->work1) yvt = scat->work1;
208: else {
209: PetscMalloc(yy_n*sizeof(PetscScalar),&yvt);
210: scat->work1 = yvt;
211: }
212: }
213: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,displx,MPIU_SCALAR,0,((PetscObject)ctx)->comm);
214: if (!rank) {
215: if (addv == ADD_VALUES) {
216: for (i=0; i<yy_n; i++) {
217: yv[i] += yvt[i];
218: }
219: #if !defined(PETSC_USE_COMPLEX)
220: } else if (addv == MAX_VALUES) {
221: for (i=0; i<yy_n; i++) {
222: yv[i] = PetscMax(yv[i],yvt[i]);
223: }
224: #endif
225: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
226: }
227: }
228: }
229: VecRestoreArray(x,&xv);
230: VecRestoreArray(y,&yv);
231: return(0);
232: }
234: /*
235: The follow to are used for both VecScatterBegin_MPI_ToAll() and VecScatterBegin_MPI_ToOne()
236: */
239: PetscErrorCode VecScatterDestroy_MPI_ToAll(VecScatter ctx)
240: {
241: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
242: PetscErrorCode ierr;
245: PetscFree(scat->work1);
246: PetscFree(scat->work2);
247: PetscFree(scat->displx);
248: PetscFree2(ctx->todata,scat->count);
249: PetscHeaderDestroy(ctx);
250: return(0);
251: }
255: PetscErrorCode VecScatterDestroy_SGtoSG(VecScatter ctx)
256: {
257: PetscErrorCode ierr;
260: PetscFree2(((VecScatter_Seq_General*)ctx->todata)->vslots,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
261: PetscFree2(ctx->todata,ctx->fromdata);
262: PetscHeaderDestroy(ctx);
263: return(0);
264: }
268: PetscErrorCode VecScatterDestroy_SGtoSS(VecScatter ctx)
269: {
270: PetscErrorCode ierr;
273: PetscFree(((VecScatter_Seq_General*)ctx->fromdata)->vslots);
274: PetscFree2(ctx->todata,ctx->fromdata);
275: PetscHeaderDestroy(ctx);
276: return(0);
277: }
281: PetscErrorCode VecScatterDestroy_SStoSG(VecScatter ctx)
282: {
283: PetscErrorCode ierr;
286: PetscFree(((VecScatter_Seq_General*)ctx->todata)->vslots);
287: PetscFree2(ctx->todata,ctx->fromdata);
288: PetscHeaderDestroy(ctx);
289: return(0);
290: }
294: PetscErrorCode VecScatterDestroy_SStoSS(VecScatter ctx)
295: {
299: PetscFree2(ctx->todata,ctx->fromdata);
300: PetscHeaderDestroy(ctx);
301: return(0);
302: }
304: /* -------------------------------------------------------------------------*/
307: PetscErrorCode VecScatterCopy_MPI_ToAll(VecScatter in,VecScatter out)
308: {
309: VecScatter_MPI_ToAll *in_to = (VecScatter_MPI_ToAll*)in->todata,*sto;
310: PetscErrorCode ierr;
311: PetscMPIInt size;
312: PetscInt i;
315: out->begin = in->begin;
316: out->end = in->end;
317: out->copy = in->copy;
318: out->destroy = in->destroy;
319: out->view = in->view;
321: MPI_Comm_size(((PetscObject)out)->comm,&size);
322: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&sto->count);
323: PetscMalloc(size*sizeof(PetscMPIInt),&sto->displx);
324: sto->type = in_to->type;
325: for (i=0; i<size; i++) {
326: sto->count[i] = in_to->count[i];
327: sto->displx[i] = in_to->displx[i];
328: }
329: sto->work1 = 0;
330: sto->work2 = 0;
331: out->todata = (void*)sto;
332: out->fromdata = (void*)0;
333: return(0);
334: }
336: /* --------------------------------------------------------------------------------------*/
337: /*
338: Scatter: sequential general to sequential general
339: */
342: PetscErrorCode VecScatterBegin_SGtoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
343: {
344: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
345: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
346: PetscErrorCode ierr;
347: PetscInt i,n = gen_from->n,*fslots,*tslots;
348: PetscScalar *xv,*yv;
349:
351: VecGetArray(x,&xv);
352: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
353: if (mode & SCATTER_REVERSE){
354: gen_to = (VecScatter_Seq_General*)ctx->fromdata;
355: gen_from = (VecScatter_Seq_General*)ctx->todata;
356: }
357: fslots = gen_from->vslots;
358: tslots = gen_to->vslots;
360: if (addv == INSERT_VALUES) {
361: for (i=0; i<n; i++) {yv[tslots[i]] = xv[fslots[i]];}
362: } else if (addv == ADD_VALUES) {
363: for (i=0; i<n; i++) {yv[tslots[i]] += xv[fslots[i]];}
364: #if !defined(PETSC_USE_COMPLEX)
365: } else if (addv == MAX_VALUES) {
366: for (i=0; i<n; i++) {yv[tslots[i]] = PetscMax(yv[tslots[i]],xv[fslots[i]]);}
367: #endif
368: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
369: VecRestoreArray(x,&xv);
370: if (x != y) {VecRestoreArray(y,&yv);}
371: return(0);
372: }
374: /*
375: Scatter: sequential general to sequential stride 1
376: */
379: PetscErrorCode VecScatterBegin_SGtoSS_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
380: {
381: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
382: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
383: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
384: PetscErrorCode ierr;
385: PetscInt first = gen_to->first;
386: PetscScalar *xv,*yv;
387:
389: VecGetArray(x,&xv);
390: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
391: if (mode & SCATTER_REVERSE){
392: xv += first;
393: if (addv == INSERT_VALUES) {
394: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
395: } else if (addv == ADD_VALUES) {
396: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
397: #if !defined(PETSC_USE_COMPLEX)
398: } else if (addv == MAX_VALUES) {
399: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
400: #endif
401: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
402: } else {
403: yv += first;
404: if (addv == INSERT_VALUES) {
405: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
406: } else if (addv == ADD_VALUES) {
407: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
408: #if !defined(PETSC_USE_COMPLEX)
409: } else if (addv == MAX_VALUES) {
410: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
411: #endif
412: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
413: }
414: VecRestoreArray(x,&xv);
415: if (x != y) {VecRestoreArray(y,&yv);}
416: return(0);
417: }
419: /*
420: Scatter: sequential general to sequential stride
421: */
424: PetscErrorCode VecScatterBegin_SGtoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
425: {
426: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
427: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
428: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
429: PetscErrorCode ierr;
430: PetscInt first = gen_to->first,step = gen_to->step;
431: PetscScalar *xv,*yv;
432:
434: VecGetArray(x,&xv);
435: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
437: if (mode & SCATTER_REVERSE){
438: if (addv == INSERT_VALUES) {
439: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
440: } else if (addv == ADD_VALUES) {
441: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
442: #if !defined(PETSC_USE_COMPLEX)
443: } else if (addv == MAX_VALUES) {
444: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
445: #endif
446: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
447: } else {
448: if (addv == INSERT_VALUES) {
449: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
450: } else if (addv == ADD_VALUES) {
451: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
452: #if !defined(PETSC_USE_COMPLEX)
453: } else if (addv == MAX_VALUES) {
454: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
455: #endif
456: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
457: }
458: VecRestoreArray(x,&xv);
459: if (x != y) {VecRestoreArray(y,&yv);}
460: return(0);
461: }
463: /*
464: Scatter: sequential stride 1 to sequential general
465: */
468: PetscErrorCode VecScatterBegin_SStoSG_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
469: {
470: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
471: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
472: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
473: PetscErrorCode ierr;
474: PetscInt first = gen_from->first;
475: PetscScalar *xv,*yv;
476:
478: VecGetArray(x,&xv);
479: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
481: if (mode & SCATTER_REVERSE){
482: yv += first;
483: if (addv == INSERT_VALUES) {
484: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
485: } else if (addv == ADD_VALUES) {
486: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
487: #if !defined(PETSC_USE_COMPLEX)
488: } else if (addv == MAX_VALUES) {
489: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
490: #endif
491: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
492: } else {
493: xv += first;
494: if (addv == INSERT_VALUES) {
495: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
496: } else if (addv == ADD_VALUES) {
497: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
498: #if !defined(PETSC_USE_COMPLEX)
499: } else if (addv == MAX_VALUES) {
500: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
501: #endif
502: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
503: }
504: VecRestoreArray(x,&xv);
505: if (x != y) {VecRestoreArray(y,&yv);}
506: return(0);
507: }
511: /*
512: Scatter: sequential stride to sequential general
513: */
514: PetscErrorCode VecScatterBegin_SStoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
515: {
516: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
517: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
518: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
519: PetscErrorCode ierr;
520: PetscInt first = gen_from->first,step = gen_from->step;
521: PetscScalar *xv,*yv;
522:
524: VecGetArray(x,&xv);
525: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
527: if (mode & SCATTER_REVERSE){
528: if (addv == INSERT_VALUES) {
529: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
530: } else if (addv == ADD_VALUES) {
531: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
532: #if !defined(PETSC_USE_COMPLEX)
533: } else if (addv == MAX_VALUES) {
534: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
535: #endif
536: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
537: } else {
538: if (addv == INSERT_VALUES) {
539: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
540: } else if (addv == ADD_VALUES) {
541: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
542: #if !defined(PETSC_USE_COMPLEX)
543: } else if (addv == MAX_VALUES) {
544: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
545: #endif
546: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
547: }
548: VecRestoreArray(x,&xv);
549: if (x != y) {VecRestoreArray(y,&yv);}
550: return(0);
551: }
553: /*
554: Scatter: sequential stride to sequential stride
555: */
558: PetscErrorCode VecScatterBegin_SStoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
559: {
560: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
561: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
562: PetscInt i,n = gen_from->n,to_first = gen_to->first,to_step = gen_to->step;
563: PetscErrorCode ierr;
564: PetscInt from_first = gen_from->first,from_step = gen_from->step;
565: PetscScalar *xv,*yv;
566:
568: VecGetArray(x,&xv);
569: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
571: if (mode & SCATTER_REVERSE){
572: from_first = gen_to->first;
573: to_first = gen_from->first;
574: from_step = gen_to->step;
575: to_step = gen_from->step;
576: }
578: if (addv == INSERT_VALUES) {
579: if (to_step == 1 && from_step == 1) {
580: PetscMemcpy(yv+to_first,xv+from_first,n*sizeof(PetscScalar));
581: } else {
582: for (i=0; i<n; i++) {
583: yv[to_first + i*to_step] = xv[from_first+i*from_step];
584: }
585: }
586: } else if (addv == ADD_VALUES) {
587: if (to_step == 1 && from_step == 1) {
588: yv += to_first; xv += from_first;
589: for (i=0; i<n; i++) {
590: yv[i] += xv[i];
591: }
592: } else {
593: for (i=0; i<n; i++) {
594: yv[to_first + i*to_step] += xv[from_first+i*from_step];
595: }
596: }
597: #if !defined(PETSC_USE_COMPLEX)
598: } else if (addv == MAX_VALUES) {
599: if (to_step == 1 && from_step == 1) {
600: yv += to_first; xv += from_first;
601: for (i=0; i<n; i++) {
602: yv[i] = PetscMax(yv[i],xv[i]);
603: }
604: } else {
605: for (i=0; i<n; i++) {
606: yv[to_first + i*to_step] = PetscMax(yv[to_first + i*to_step],xv[from_first+i*from_step]);
607: }
608: }
609: #endif
610: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
611: VecRestoreArray(x,&xv);
612: if (x != y) {VecRestoreArray(y,&yv);}
613: return(0);
614: }
616: /* --------------------------------------------------------------------------*/
621: PetscErrorCode VecScatterCopy_SGToSG(VecScatter in,VecScatter out)
622: {
623: PetscErrorCode ierr;
624: VecScatter_Seq_General *in_to = (VecScatter_Seq_General*)in->todata,*out_to = PETSC_NULL;
625: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
626:
628: out->begin = in->begin;
629: out->end = in->end;
630: out->copy = in->copy;
631: out->destroy = in->destroy;
632: out->view = in->view;
634: PetscMalloc2(1,VecScatter_Seq_General,&out_to,1,VecScatter_Seq_General,&out_from);
635: PetscMalloc2(in_to->n,PetscInt,&out_to->vslots,in_from->n,PetscInt,&out_from->vslots);
636: out_to->n = in_to->n;
637: out_to->type = in_to->type;
638: out_to->nonmatching_computed = PETSC_FALSE;
639: out_to->slots_nonmatching = 0;
640: out_to->is_copy = PETSC_FALSE;
641: PetscMemcpy(out_to->vslots,in_to->vslots,(out_to->n)*sizeof(PetscInt));
643: out_from->n = in_from->n;
644: out_from->type = in_from->type;
645: out_from->nonmatching_computed = PETSC_FALSE;
646: out_from->slots_nonmatching = 0;
647: out_from->is_copy = PETSC_FALSE;
648: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
650: out->todata = (void*)out_to;
651: out->fromdata = (void*)out_from;
652: return(0);
653: }
658: PetscErrorCode VecScatterCopy_SGToSS(VecScatter in,VecScatter out)
659: {
660: PetscErrorCode ierr;
661: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to = PETSC_NULL;
662: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
663:
665: out->begin = in->begin;
666: out->end = in->end;
667: out->copy = in->copy;
668: out->destroy = in->destroy;
669: out->view = in->view;
671: PetscMalloc2(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_General,&out_from);
672: PetscMalloc(in_from->n*sizeof(PetscInt),&out_from->vslots);
673: out_to->n = in_to->n;
674: out_to->type = in_to->type;
675: out_to->first = in_to->first;
676: out_to->step = in_to->step;
677: out_to->type = in_to->type;
679: out_from->n = in_from->n;
680: out_from->type = in_from->type;
681: out_from->nonmatching_computed = PETSC_FALSE;
682: out_from->slots_nonmatching = 0;
683: out_from->is_copy = PETSC_FALSE;
684: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
686: out->todata = (void*)out_to;
687: out->fromdata = (void*)out_from;
688: return(0);
689: }
691: /* --------------------------------------------------------------------------*/
692: /*
693: Scatter: parallel to sequential vector, sequential strides for both.
694: */
697: PetscErrorCode VecScatterCopy_SStoSS(VecScatter in,VecScatter out)
698: {
699: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to = PETSC_NULL;
700: VecScatter_Seq_Stride *in_from = (VecScatter_Seq_Stride*)in->fromdata,*out_from = PETSC_NULL;
701: PetscErrorCode ierr;
704: out->begin = in->begin;
705: out->end = in->end;
706: out->copy = in->copy;
707: out->destroy = in->destroy;
708: out->view = in->view;
710: PetscMalloc2(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_Stride,&out_from);
711: out_to->n = in_to->n;
712: out_to->type = in_to->type;
713: out_to->first = in_to->first;
714: out_to->step = in_to->step;
715: out_to->type = in_to->type;
716: out_from->n = in_from->n;
717: out_from->type = in_from->type;
718: out_from->first = in_from->first;
719: out_from->step = in_from->step;
720: out_from->type = in_from->type;
721: out->todata = (void*)out_to;
722: out->fromdata = (void*)out_from;
723: return(0);
724: }
726: EXTERN PetscErrorCode VecScatterCreate_PtoS(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,PetscInt,VecScatter);
727: EXTERN PetscErrorCode VecScatterCreate_PtoP(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,VecScatter);
728: EXTERN PetscErrorCode VecScatterCreate_StoP(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,PetscInt,VecScatter);
730: /* =======================================================================*/
731: #define VEC_SEQ_ID 0
732: #define VEC_MPI_ID 1
734: /*
735: Blocksizes we have optimized scatters for
736: */
738: #define VecScatterOptimizedBS(mbs) ((2 <= mbs && mbs <= 8) || mbs == 12)
740: PetscErrorCode VecScatterCreateEmpty(MPI_Comm comm,VecScatter *newctx)
741: {
742: VecScatter ctx;
746: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
747: ctx->inuse = PETSC_FALSE;
748: ctx->beginandendtogether = PETSC_FALSE;
749: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether,PETSC_NULL);
750: if (ctx->beginandendtogether) {
751: PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
752: }
753: ctx->packtogether = PETSC_FALSE;
754: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether,PETSC_NULL);
755: if (ctx->packtogether) {
756: PetscInfo(ctx,"Pack all messages before sending\n");
757: }
758: *newctx = ctx;
759: return(0);
760: }
764: /*@C
765: VecScatterCreate - Creates a vector scatter context.
767: Collective on Vec
769: Input Parameters:
770: + xin - a vector that defines the shape (parallel data layout of the vector)
771: of vectors from which we scatter
772: . yin - a vector that defines the shape (parallel data layout of the vector)
773: of vectors to which we scatter
774: . ix - the indices of xin to scatter (if PETSC_NULL scatters all values)
775: - iy - the indices of yin to hold results (if PETSC_NULL fills entire vector yin)
777: Output Parameter:
778: . newctx - location to store the new scatter context
780: Options Database Keys: (uses regular MPI_Sends by default)
781: + -vecscatter_ssend - Uses MPI_Ssend_init() instead of MPI_Send_init()
782: . -vecscatter_rsend - use ready receiver mode for MPI sends
783: . -vecscatter_merge - VecScatterBegin() handles all of the communication, VecScatterEnd() is a nop
784: eliminates the chance for overlap of computation and communication
785: . -vecscatter_sendfirst - Posts sends before receives
786: . -vecscatter_packtogether - Pack all messages before sending, receive all messages before unpacking
787: . -vecscatter_alltoall - Uses MPI all to all communication for scatter
788: - -vecscatter_nopack - Avoid packing to work vector when possible (if used with -vecscatter_alltoall then will use MPI_Alltoallw()
790: $
791: $ --When packing is used--
792: $ MPI Datatypes (no packing) sendfirst merge packtogether persistent* -vecscatter_
793: $
794: $ Message passing Send p X X X always
795: $ Ssend p X X X always _ssend
796: $ Rsend p nonsense X X always _rsend
797: $ AlltoAll v or w X nonsense always X nonsense _alltoall
798: $ MPI_Win p nonsense p p nonsense _window
799: $
800: $ -vecscatter_ _nopack _sendfirst _merge _packtogether
801: $
802: $ Since persistent sends and receives require a constant memory address they can only be used when data is packed into the work vector
803: $ because the in and out array may be different for each call to VecScatterBegin/End().
804: $
805: $ p indicates possible, but not implemented. X indicates implemented
806: $
808: Level: intermediate
810: Notes:
811: In calls to VecScatter() you can use different vectors than the xin and
812: yin you used above; BUT they must have the same parallel data layout, for example,
813: they could be obtained from VecDuplicate().
814: A VecScatter context CANNOT be used in two or more simultaneous scatters;
815: that is you cannot call a second VecScatterBegin() with the same scatter
816: context until the VecScatterEnd() has been called on the first VecScatterBegin().
817: In this case a separate VecScatter is needed for each concurrent scatter.
819: Currently the MPI_Send(), MPI_Ssend() and MPI_Rsend() all use PERSISTENT versions.
820: (this unfortunately requires that the same in and out arrays be used for each use, this
821: is why when not using MPI_alltoallw() we always need to pack the input into the work array before sending
822: and unpack upon receeving instead of using MPI datatypes to avoid the packing/unpacking).
824: Concepts: scatter^between vectors
825: Concepts: gather^between vectors
827: .seealso: VecScatterDestroy(), VecScatterCreateToAll(), VecScatterCreateToZero()
828: @*/
829: PetscErrorCode VecScatterCreate(Vec xin,IS ix,Vec yin,IS iy,VecScatter *newctx)
830: {
831: VecScatter ctx;
833: PetscMPIInt size;
834: PetscInt totalv,xin_type = VEC_SEQ_ID,yin_type = VEC_SEQ_ID,*range;
835: MPI_Comm comm,ycomm;
836: PetscTruth ixblock,iyblock,iystride,islocal,cando,flag;
837: IS tix = 0,tiy = 0;
841: /*
842: Determine if the vectors are "parallel", ie. it shares a comm with other processors, or
843: sequential (it does not share a comm). The difference is that parallel vectors treat the
844: index set as providing indices in the global parallel numbering of the vector, with
845: sequential vectors treat the index set as providing indices in the local sequential
846: numbering
847: */
848: PetscObjectGetComm((PetscObject)xin,&comm);
849: MPI_Comm_size(comm,&size);
850: if (size > 1) {xin_type = VEC_MPI_ID;}
852: PetscObjectGetComm((PetscObject)yin,&ycomm);
853: MPI_Comm_size(ycomm,&size);
854: if (size > 1) {comm = ycomm; yin_type = VEC_MPI_ID;}
855:
856: /* generate the Scatter context */
857: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
858: ctx->inuse = PETSC_FALSE;
860: ctx->beginandendtogether = PETSC_FALSE;
861: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether,PETSC_NULL);
862: if (ctx->beginandendtogether) {
863: PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
864: }
865: ctx->packtogether = PETSC_FALSE;
866: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether,PETSC_NULL);
867: if (ctx->packtogether) {
868: PetscInfo(ctx,"Pack all messages before sending\n");
869: }
871: VecGetLocalSize(xin,&ctx->from_n);
872: VecGetLocalSize(yin,&ctx->to_n);
874: /*
875: if ix or iy is not included; assume just grabbing entire vector
876: */
877: if (!ix && xin_type == VEC_SEQ_ID) {
878: ISCreateStride(comm,ctx->from_n,0,1,&ix);
879: tix = ix;
880: } else if (!ix && xin_type == VEC_MPI_ID) {
881: if (yin_type == VEC_MPI_ID) {
882: PetscInt ntmp, low;
883: VecGetLocalSize(xin,&ntmp);
884: VecGetOwnershipRange(xin,&low,PETSC_NULL);
885: ISCreateStride(comm,ntmp,low,1,&ix);
886: } else{
887: PetscInt Ntmp;
888: VecGetSize(xin,&Ntmp);
889: ISCreateStride(comm,Ntmp,0,1,&ix);
890: }
891: tix = ix;
892: } else if (!ix) {
893: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ix not given, but not Seq or MPI vector");
894: }
896: if (!iy && yin_type == VEC_SEQ_ID) {
897: ISCreateStride(comm,ctx->to_n,0,1,&iy);
898: tiy = iy;
899: } else if (!iy && yin_type == VEC_MPI_ID) {
900: if (xin_type == VEC_MPI_ID) {
901: PetscInt ntmp, low;
902: VecGetLocalSize(yin,&ntmp);
903: VecGetOwnershipRange(yin,&low,PETSC_NULL);
904: ISCreateStride(comm,ntmp,low,1,&iy);
905: } else{
906: PetscInt Ntmp;
907: VecGetSize(yin,&Ntmp);
908: ISCreateStride(comm,Ntmp,0,1,&iy);
909: }
910: tiy = iy;
911: } else if (!iy) {
912: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
913: }
915: /* ===========================================================================================================
916: Check for special cases
917: ==========================================================================================================*/
918: /* ---------------------------------------------------------------------------*/
919: if (xin_type == VEC_SEQ_ID && yin_type == VEC_SEQ_ID) {
920: if (((PetscObject)ix)->type == IS_GENERAL && ((PetscObject)iy)->type == IS_GENERAL){
921: PetscInt nx,ny;
922: const PetscInt *idx,*idy;
923: VecScatter_Seq_General *to = PETSC_NULL,*from = PETSC_NULL;
925: ISGetLocalSize(ix,&nx);
926: ISGetLocalSize(iy,&ny);
927: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
928: ISGetIndices(ix,&idx);
929: ISGetIndices(iy,&idy);
930: PetscMalloc2(1,VecScatter_Seq_General,&to,1,VecScatter_Seq_General,&from);
931: PetscMalloc2(nx,PetscInt,&to->vslots,nx,PetscInt,&from->vslots);
932: to->n = nx;
933: #if defined(PETSC_USE_DEBUG)
934: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
935: #endif
936: PetscMemcpy(to->vslots,idy,nx*sizeof(PetscInt));
937: from->n = nx;
938: #if defined(PETSC_USE_DEBUG)
939: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
940: #endif
941: PetscMemcpy(from->vslots,idx,nx*sizeof(PetscInt));
942: to->type = VEC_SCATTER_SEQ_GENERAL;
943: from->type = VEC_SCATTER_SEQ_GENERAL;
944: ctx->todata = (void*)to;
945: ctx->fromdata = (void*)from;
946: ctx->begin = VecScatterBegin_SGtoSG;
947: ctx->end = 0;
948: ctx->destroy = VecScatterDestroy_SGtoSG;
949: ctx->copy = VecScatterCopy_SGToSG;
950: PetscInfo(xin,"Special case: sequential vector general scatter\n");
951: goto functionend;
952: } else if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
953: PetscInt nx,ny,to_first,to_step,from_first,from_step;
954: VecScatter_Seq_Stride *from8 = PETSC_NULL,*to8 = PETSC_NULL;
956: ISGetLocalSize(ix,&nx);
957: ISGetLocalSize(iy,&ny);
958: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
959: ISStrideGetInfo(iy,&to_first,&to_step);
960: ISStrideGetInfo(ix,&from_first,&from_step);
961: PetscMalloc2(1,VecScatter_Seq_Stride,&to8,1,VecScatter_Seq_Stride,&from8);
962: to8->n = nx;
963: to8->first = to_first;
964: to8->step = to_step;
965: from8->n = nx;
966: from8->first = from_first;
967: from8->step = from_step;
968: to8->type = VEC_SCATTER_SEQ_STRIDE;
969: from8->type = VEC_SCATTER_SEQ_STRIDE;
970: ctx->todata = (void*)to8;
971: ctx->fromdata = (void*)from8;
972: ctx->begin = VecScatterBegin_SStoSS;
973: ctx->end = 0;
974: ctx->destroy = VecScatterDestroy_SStoSS;
975: ctx->copy = VecScatterCopy_SStoSS;
976: PetscInfo(xin,"Special case: sequential vector stride to stride\n");
977: goto functionend;
978: } else if (((PetscObject)ix)->type == IS_GENERAL && ((PetscObject)iy)->type == IS_STRIDE){
979: PetscInt nx,ny,first,step;
980: const PetscInt *idx;
981: VecScatter_Seq_General *from9 = PETSC_NULL;
982: VecScatter_Seq_Stride *to9 = PETSC_NULL;
984: ISGetLocalSize(ix,&nx);
985: ISGetIndices(ix,&idx);
986: ISGetLocalSize(iy,&ny);
987: ISStrideGetInfo(iy,&first,&step);
988: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
989: PetscMalloc2(1,VecScatter_Seq_Stride,&to9,1,VecScatter_Seq_General,&from9);
990: PetscMalloc(nx*sizeof(PetscInt),&from9->vslots);
991: to9->n = nx;
992: to9->first = first;
993: to9->step = step;
994: from9->n = nx;
995: #if defined(PETSC_USE_DEBUG)
996: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
997: #endif
998: PetscMemcpy(from9->vslots,idx,nx*sizeof(PetscInt));
999: ctx->todata = (void*)to9; ctx->fromdata = (void*)from9;
1000: if (step == 1) ctx->begin = VecScatterBegin_SGtoSS_Stride1;
1001: else ctx->begin = VecScatterBegin_SGtoSS;
1002: ctx->destroy = VecScatterDestroy_SGtoSS;
1003: ctx->end = 0;
1004: ctx->copy = VecScatterCopy_SGToSS;
1005: to9->type = VEC_SCATTER_SEQ_STRIDE;
1006: from9->type = VEC_SCATTER_SEQ_GENERAL;
1007: PetscInfo(xin,"Special case: sequential vector general to stride\n");
1008: goto functionend;
1009: } else if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_GENERAL){
1010: PetscInt nx,ny,first,step;
1011: const PetscInt *idy;
1012: VecScatter_Seq_General *to10 = PETSC_NULL;
1013: VecScatter_Seq_Stride *from10 = PETSC_NULL;
1015: ISGetLocalSize(ix,&nx);
1016: ISGetIndices(iy,&idy);
1017: ISGetLocalSize(iy,&ny);
1018: ISStrideGetInfo(ix,&first,&step);
1019: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1020: PetscMalloc2(1,VecScatter_Seq_General,&to10,1,VecScatter_Seq_Stride,&from10);
1021: PetscMalloc(nx*sizeof(PetscInt),&to10->vslots);
1022: from10->n = nx;
1023: from10->first = first;
1024: from10->step = step;
1025: to10->n = nx;
1026: #if defined(PETSC_USE_DEBUG)
1027: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1028: #endif
1029: PetscMemcpy(to10->vslots,idy,nx*sizeof(PetscInt));
1030: ctx->todata = (void*)to10;
1031: ctx->fromdata = (void*)from10;
1032: if (step == 1) ctx->begin = VecScatterBegin_SStoSG_Stride1;
1033: else ctx->begin = VecScatterBegin_SStoSG;
1034: ctx->destroy = VecScatterDestroy_SStoSG;
1035: ctx->end = 0;
1036: ctx->copy = 0;
1037: to10->type = VEC_SCATTER_SEQ_GENERAL;
1038: from10->type = VEC_SCATTER_SEQ_STRIDE;
1039: PetscInfo(xin,"Special case: sequential vector stride to general\n");
1040: goto functionend;
1041: } else {
1042: PetscInt nx,ny;
1043: const PetscInt *idx,*idy;
1044: VecScatter_Seq_General *to11 = PETSC_NULL,*from11 = PETSC_NULL;
1045: PetscTruth idnx,idny;
1047: ISGetLocalSize(ix,&nx);
1048: ISGetLocalSize(iy,&ny);
1049: if (nx != ny) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match, in %D out %D",nx,ny);
1051: ISIdentity(ix,&idnx);
1052: ISIdentity(iy,&idny);
1053: if (idnx && idny) {
1054: VecScatter_Seq_Stride *to13 = PETSC_NULL,*from13 = PETSC_NULL;
1055: PetscMalloc2(1,VecScatter_Seq_Stride,&to13,1,VecScatter_Seq_Stride,&from13);
1056: to13->n = nx;
1057: to13->first = 0;
1058: to13->step = 1;
1059: from13->n = nx;
1060: from13->first = 0;
1061: from13->step = 1;
1062: to13->type = VEC_SCATTER_SEQ_STRIDE;
1063: from13->type = VEC_SCATTER_SEQ_STRIDE;
1064: ctx->todata = (void*)to13;
1065: ctx->fromdata = (void*)from13;
1066: ctx->begin = VecScatterBegin_SStoSS;
1067: ctx->end = 0;
1068: ctx->destroy = VecScatterDestroy_SStoSS;
1069: ctx->copy = VecScatterCopy_SStoSS;
1070: PetscInfo(xin,"Special case: sequential copy\n");
1071: goto functionend;
1072: }
1074: ISGetIndices(iy,&idy);
1075: ISGetIndices(ix,&idx);
1076: PetscMalloc2(1,VecScatter_Seq_General,&to11,1,VecScatter_Seq_General,&from11);
1077: PetscMalloc2(nx,PetscInt,&to11->vslots,nx,PetscInt,&from11->vslots);
1078: to11->n = nx;
1079: #if defined(PETSC_USE_DEBUG)
1080: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1081: #endif
1082: PetscMemcpy(to11->vslots,idy,nx*sizeof(PetscInt));
1083: from11->n = nx;
1084: #if defined(PETSC_USE_DEBUG)
1085: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
1086: #endif
1087: PetscMemcpy(from11->vslots,idx,nx*sizeof(PetscInt));
1088: to11->type = VEC_SCATTER_SEQ_GENERAL;
1089: from11->type = VEC_SCATTER_SEQ_GENERAL;
1090: ctx->todata = (void*)to11;
1091: ctx->fromdata = (void*)from11;
1092: ctx->begin = VecScatterBegin_SGtoSG;
1093: ctx->end = 0;
1094: ctx->destroy = VecScatterDestroy_SGtoSG;
1095: ctx->copy = VecScatterCopy_SGToSG;
1096: ISRestoreIndices(ix,&idx);
1097: ISRestoreIndices(iy,&idy);
1098: PetscInfo(xin,"Sequential vector scatter with block indices\n");
1099: goto functionend;
1100: }
1101: }
1102: /* ---------------------------------------------------------------------------*/
1103: if (xin_type == VEC_MPI_ID && yin_type == VEC_SEQ_ID) {
1105: /* ===========================================================================================================
1106: Check for special cases
1107: ==========================================================================================================*/
1108: islocal = PETSC_FALSE;
1109: /* special case extracting (subset of) local portion */
1110: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1111: PetscInt nx,ny,to_first,to_step,from_first,from_step;
1112: PetscInt start,end;
1113: VecScatter_Seq_Stride *from12 = PETSC_NULL,*to12 = PETSC_NULL;
1115: VecGetOwnershipRange(xin,&start,&end);
1116: ISGetLocalSize(ix,&nx);
1117: ISStrideGetInfo(ix,&from_first,&from_step);
1118: ISGetLocalSize(iy,&ny);
1119: ISStrideGetInfo(iy,&to_first,&to_step);
1120: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1121: if (ix->min >= start && ix->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1122: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1123: if (cando) {
1124: PetscMalloc2(1,VecScatter_Seq_Stride,&to12,1,VecScatter_Seq_Stride,&from12);
1125: to12->n = nx;
1126: to12->first = to_first;
1127: to12->step = to_step;
1128: from12->n = nx;
1129: from12->first = from_first-start;
1130: from12->step = from_step;
1131: to12->type = VEC_SCATTER_SEQ_STRIDE;
1132: from12->type = VEC_SCATTER_SEQ_STRIDE;
1133: ctx->todata = (void*)to12;
1134: ctx->fromdata = (void*)from12;
1135: ctx->begin = VecScatterBegin_SStoSS;
1136: ctx->end = 0;
1137: ctx->destroy = VecScatterDestroy_SStoSS;
1138: ctx->copy = VecScatterCopy_SStoSS;
1139: PetscInfo(xin,"Special case: processors only getting local values\n");
1140: goto functionend;
1141: }
1142: } else {
1143: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1144: }
1146: /* test for special case of all processors getting entire vector */
1147: /* contains check that PetscMPIInt can handle the sizes needed */
1148: totalv = 0;
1149: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1150: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1151: PetscMPIInt *count = PETSC_NULL,*displx;
1152: VecScatter_MPI_ToAll *sto = PETSC_NULL;
1154: ISGetLocalSize(ix,&nx);
1155: ISStrideGetInfo(ix,&from_first,&from_step);
1156: ISGetLocalSize(iy,&ny);
1157: ISStrideGetInfo(iy,&to_first,&to_step);
1158: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1159: VecGetSize(xin,&N);
1160: if (nx != N) {
1161: totalv = 0;
1162: } else if (from_first == 0 && from_step == 1 && from_first == to_first && from_step == to_step){
1163: totalv = 1;
1164: } else totalv = 0;
1165: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1167: #if defined(PETSC_USE_64BIT_INDICES)
1168: if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1169: #else
1170: if (cando) {
1171: #endif
1172: MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1173: PetscMalloc(size*sizeof(PetscMPIInt),&displx);
1174: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1175: range = xin->map->range;
1176: for (i=0; i<size; i++) {
1177: count[i] = PetscMPIIntCast(range[i+1] - range[i]);
1178: displx[i] = PetscMPIIntCast(range[i]);
1179: }
1180: sto->count = count;
1181: sto->displx = displx;
1182: sto->work1 = 0;
1183: sto->work2 = 0;
1184: sto->type = VEC_SCATTER_MPI_TOALL;
1185: ctx->todata = (void*)sto;
1186: ctx->fromdata = 0;
1187: ctx->begin = VecScatterBegin_MPI_ToAll;
1188: ctx->end = 0;
1189: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1190: ctx->copy = VecScatterCopy_MPI_ToAll;
1191: PetscInfo(xin,"Special case: all processors get entire parallel vector\n");
1192: goto functionend;
1193: }
1194: } else {
1195: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1196: }
1198: /* test for special case of processor 0 getting entire vector */
1199: /* contains check that PetscMPIInt can handle the sizes needed */
1200: totalv = 0;
1201: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1202: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1203: PetscMPIInt rank,*count = PETSC_NULL,*displx;
1204: VecScatter_MPI_ToAll *sto = PETSC_NULL;
1206: PetscObjectGetComm((PetscObject)xin,&comm);
1207: MPI_Comm_rank(comm,&rank);
1208: ISGetLocalSize(ix,&nx);
1209: ISStrideGetInfo(ix,&from_first,&from_step);
1210: ISGetLocalSize(iy,&ny);
1211: ISStrideGetInfo(iy,&to_first,&to_step);
1212: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1213: if (!rank) {
1214: VecGetSize(xin,&N);
1215: if (nx != N) {
1216: totalv = 0;
1217: } else if (from_first == 0 && from_step == 1 &&
1218: from_first == to_first && from_step == to_step){
1219: totalv = 1;
1220: } else totalv = 0;
1221: } else {
1222: if (!nx) totalv = 1;
1223: else totalv = 0;
1224: }
1225: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1227: #if defined(PETSC_USE_64BIT_INDICES)
1228: if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1229: #else
1230: if (cando) {
1231: #endif
1232: MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1233: PetscMalloc(size*sizeof(PetscMPIInt),&displx);
1234: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1235: range = xin->map->range;
1236: for (i=0; i<size; i++) {
1237: count[i] = PetscMPIIntCast(range[i+1] - range[i]);
1238: displx[i] = PetscMPIIntCast(range[i]);
1239: }
1240: sto->count = count;
1241: sto->displx = displx;
1242: sto->work1 = 0;
1243: sto->work2 = 0;
1244: sto->type = VEC_SCATTER_MPI_TOONE;
1245: ctx->todata = (void*)sto;
1246: ctx->fromdata = 0;
1247: ctx->begin = VecScatterBegin_MPI_ToOne;
1248: ctx->end = 0;
1249: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1250: ctx->copy = VecScatterCopy_MPI_ToAll;
1251: PetscInfo(xin,"Special case: processor zero gets entire parallel vector, rest get none\n");
1252: goto functionend;
1253: }
1254: } else {
1255: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1256: }
1258: ISBlock(ix,&ixblock);
1259: ISBlock(iy,&iyblock);
1260: ISStride(iy,&iystride);
1261: if (ixblock) {
1262: /* special case block to block */
1263: if (iyblock) {
1264: PetscInt nx,ny,bsx,bsy;
1265: const PetscInt *idx,*idy;
1266: ISBlockGetBlockSize(iy,&bsy);
1267: ISBlockGetBlockSize(ix,&bsx);
1268: if (bsx == bsy && VecScatterOptimizedBS(bsx)) {
1269: ISBlockGetLocalSize(ix,&nx);
1270: ISBlockGetIndices(ix,&idx);
1271: ISBlockGetLocalSize(iy,&ny);
1272: ISBlockGetIndices(iy,&idy);
1273: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1274: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,bsx,ctx);
1275: ISBlockRestoreIndices(ix,&idx);
1276: ISBlockRestoreIndices(iy,&idy);
1277: PetscInfo(xin,"Special case: blocked indices\n");
1278: goto functionend;
1279: }
1280: /* special case block to stride */
1281: } else if (iystride) {
1282: PetscInt ystart,ystride,ysize,bsx;
1283: ISStrideGetInfo(iy,&ystart,&ystride);
1284: ISGetLocalSize(iy,&ysize);
1285: ISBlockGetBlockSize(ix,&bsx);
1286: /* see if stride index set is equivalent to block index set */
1287: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1288: PetscInt nx,il,*idy;
1289: const PetscInt *idx;
1290: ISBlockGetLocalSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1291: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1292: PetscMalloc(nx*sizeof(PetscInt),&idy);
1293: if (nx) {
1294: idy[0] = ystart;
1295: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1296: }
1297: VecScatterCreate_PtoS(nx,idx,nx,idy,xin,yin,bsx,ctx);
1298: PetscFree(idy);
1299: ISBlockRestoreIndices(ix,&idx);
1300: PetscInfo(xin,"Special case: blocked indices to stride\n");
1301: goto functionend;
1302: }
1303: }
1304: }
1305: /* left over general case */
1306: {
1307: PetscInt nx,ny;
1308: const PetscInt *idx,*idy;
1309: ISGetLocalSize(ix,&nx);
1310: ISGetIndices(ix,&idx);
1311: ISGetLocalSize(iy,&ny);
1312: ISGetIndices(iy,&idy);
1313: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1314: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,1,ctx);
1315: ISRestoreIndices(ix,&idx);
1316: ISRestoreIndices(iy,&idy);
1317: PetscInfo(xin,"General case: MPI to Seq\n");
1318: goto functionend;
1319: }
1320: }
1321: /* ---------------------------------------------------------------------------*/
1322: if (xin_type == VEC_SEQ_ID && yin_type == VEC_MPI_ID) {
1323: /* ===========================================================================================================
1324: Check for special cases
1325: ==========================================================================================================*/
1326: /* special case local copy portion */
1327: islocal = PETSC_FALSE;
1328: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1329: PetscInt nx,ny,to_first,to_step,from_step,start,end,from_first;
1330: VecScatter_Seq_Stride *from = PETSC_NULL,*to = PETSC_NULL;
1332: VecGetOwnershipRange(yin,&start,&end);
1333: ISGetLocalSize(ix,&nx);
1334: ISStrideGetInfo(ix,&from_first,&from_step);
1335: ISGetLocalSize(iy,&ny);
1336: ISStrideGetInfo(iy,&to_first,&to_step);
1337: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1338: if (iy->min >= start && iy->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1339: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1340: if (cando) {
1341: PetscMalloc2(1,VecScatter_Seq_Stride,&to,1,VecScatter_Seq_Stride,&from);
1342: to->n = nx;
1343: to->first = to_first-start;
1344: to->step = to_step;
1345: from->n = nx;
1346: from->first = from_first;
1347: from->step = from_step;
1348: to->type = VEC_SCATTER_SEQ_STRIDE;
1349: from->type = VEC_SCATTER_SEQ_STRIDE;
1350: ctx->todata = (void*)to;
1351: ctx->fromdata = (void*)from;
1352: ctx->begin = VecScatterBegin_SStoSS;
1353: ctx->end = 0;
1354: ctx->destroy = VecScatterDestroy_SStoSS;
1355: ctx->copy = VecScatterCopy_SStoSS;
1356: PetscInfo(xin,"Special case: sequential stride to MPI stride\n");
1357: goto functionend;
1358: }
1359: } else {
1360: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1361: }
1362: /* special case block to stride */
1363: if (((PetscObject)ix)->type == IS_BLOCK && ((PetscObject)iy)->type == IS_STRIDE){
1364: PetscInt ystart,ystride,ysize,bsx;
1365: ISStrideGetInfo(iy,&ystart,&ystride);
1366: ISGetLocalSize(iy,&ysize);
1367: ISBlockGetBlockSize(ix,&bsx);
1368: /* see if stride index set is equivalent to block index set */
1369: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1370: PetscInt nx,il,*idy;
1371: const PetscInt *idx;
1372: ISBlockGetLocalSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1373: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1374: PetscMalloc(nx*sizeof(PetscInt),&idy);
1375: if (nx) {
1376: idy[0] = ystart;
1377: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1378: }
1379: VecScatterCreate_StoP(nx,idx,nx,idy,xin,yin,bsx,ctx);
1380: PetscFree(idy);
1381: ISBlockRestoreIndices(ix,&idx);
1382: PetscInfo(xin,"Special case: Blocked indices to stride\n");
1383: goto functionend;
1384: }
1385: }
1387: /* general case */
1388: {
1389: PetscInt nx,ny;
1390: const PetscInt *idx,*idy;
1391: ISGetLocalSize(ix,&nx);
1392: ISGetIndices(ix,&idx);
1393: ISGetLocalSize(iy,&ny);
1394: ISGetIndices(iy,&idy);
1395: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1396: VecScatterCreate_StoP(nx,idx,ny,idy,xin,yin,1,ctx);
1397: ISRestoreIndices(ix,&idx);
1398: ISRestoreIndices(iy,&idy);
1399: PetscInfo(xin,"General case: Seq to MPI\n");
1400: goto functionend;
1401: }
1402: }
1403: /* ---------------------------------------------------------------------------*/
1404: if (xin_type == VEC_MPI_ID && yin_type == VEC_MPI_ID) {
1405: /* no special cases for now */
1406: PetscInt nx,ny;
1407: const PetscInt *idx,*idy;
1408: ISGetLocalSize(ix,&nx);
1409: ISGetIndices(ix,&idx);
1410: ISGetLocalSize(iy,&ny);
1411: ISGetIndices(iy,&idy);
1412: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1413: VecScatterCreate_PtoP(nx,idx,ny,idy,xin,yin,ctx);
1414: ISRestoreIndices(ix,&idx);
1415: ISRestoreIndices(iy,&idy);
1416: PetscInfo(xin,"General case: MPI to MPI\n");
1417: goto functionend;
1418: }
1420: functionend:
1421: *newctx = ctx;
1422: if (tix) {ISDestroy(tix);}
1423: if (tiy) {ISDestroy(tiy);}
1424: flag = PETSC_FALSE;
1425: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_view_info",&flag,PETSC_NULL);
1426: if (flag) {
1427: PetscViewer viewer;
1428: PetscViewerASCIIGetStdout(comm,&viewer);
1429: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1430: VecScatterView(ctx,viewer);
1431: PetscViewerPopFormat(viewer);
1432: }
1433: flag = PETSC_FALSE;
1434: PetscOptionsGetTruth(PETSC_NULL,"-vecscatter_view",&flag,PETSC_NULL);
1435: if (flag) {
1436: PetscViewer viewer;
1437: PetscViewerASCIIGetStdout(comm,&viewer);
1438: VecScatterView(ctx,viewer);
1439: }
1440: return(0);
1441: }
1443: /* ------------------------------------------------------------------*/
1446: /*@
1447: VecScatterGetMerged - Returns true if the scatter is completed in the VecScatterBegin()
1448: and the VecScatterEnd() does nothing
1450: Not Collective
1452: Input Parameter:
1453: . ctx - scatter context created with VecScatterCreate()
1455: Output Parameter:
1456: . flg - PETSC_TRUE if the VecScatterBegin/End() are all done during the VecScatterBegin()
1458: Level: developer
1460: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1461: @*/
1462: PetscErrorCode VecScatterGetMerged(VecScatter ctx,PetscTruth *flg)
1463: {
1466: *flg = ctx->beginandendtogether;
1467: return(0);
1468: }
1472: /*@
1473: VecScatterBegin - Begins a generalized scatter from one vector to
1474: another. Complete the scattering phase with VecScatterEnd().
1476: Collective on VecScatter and Vec
1478: Input Parameters:
1479: + inctx - scatter context generated by VecScatterCreate()
1480: . x - the vector from which we scatter
1481: . y - the vector to which we scatter
1482: . addv - either ADD_VALUES or INSERT_VALUES, with INSERT_VALUES mode any location
1483: not scattered to retains its old value; i.e. the vector is NOT first zeroed.
1484: - mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1485: SCATTER_FORWARD or SCATTER_REVERSE
1488: Level: intermediate
1490: Options Database: See VecScatterCreate()
1492: Notes:
1493: The vectors x and y need not be the same vectors used in the call
1494: to VecScatterCreate(), but x must have the same parallel data layout
1495: as that passed in as the x to VecScatterCreate(), similarly for the y.
1496: Most likely they have been obtained from VecDuplicate().
1498: You cannot change the values in the input vector between the calls to VecScatterBegin()
1499: and VecScatterEnd().
1501: If you use SCATTER_REVERSE the two arguments x and y should be reversed, from
1502: the SCATTER_FORWARD.
1503:
1504: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1506: This scatter is far more general than the conventional
1507: scatter, since it can be a gather or a scatter or a combination,
1508: depending on the indices ix and iy. If x is a parallel vector and y
1509: is sequential, VecScatterBegin() can serve to gather values to a
1510: single processor. Similarly, if y is parallel and x sequential, the
1511: routine can scatter from one processor to many processors.
1513: Concepts: scatter^between vectors
1514: Concepts: gather^between vectors
1516: .seealso: VecScatterCreate(), VecScatterEnd()
1517: @*/
1518: PetscErrorCode VecScatterBegin(VecScatter inctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1519: {
1521: #if defined(PETSC_USE_DEBUG)
1522: PetscInt to_n,from_n;
1523: #endif
1529: if (inctx->inuse) SETERRQ(PETSC_ERR_ARG_WRONGSTATE," Scatter ctx already in use");
1530: CHKMEMQ;
1532: #if defined(PETSC_USE_DEBUG)
1533: /*
1534: Error checking to make sure these vectors match the vectors used
1535: to create the vector scatter context. -1 in the from_n and to_n indicate the
1536: vector lengths are unknown (for example with mapped scatters) and thus
1537: no error checking is performed.
1538: */
1539: if (inctx->from_n >= 0 && inctx->to_n >= 0) {
1540: VecGetLocalSize(x,&from_n);
1541: VecGetLocalSize(y,&to_n);
1542: if (mode & SCATTER_REVERSE) {
1543: if (to_n != inctx->from_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter reverse and vector to != ctx from size)",to_n,inctx->from_n);
1544: if (from_n != inctx->to_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter reverse and vector from != ctx to size)",from_n,inctx->to_n);
1545: } else {
1546: if (to_n != inctx->to_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter forward and vector to != ctx to size)",to_n,inctx->to_n);
1547: if (from_n != inctx->from_n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter forward and vector from != ctx from size)",from_n,inctx->from_n);
1548: }
1549: }
1550: #endif
1552: inctx->inuse = PETSC_TRUE;
1553: PetscLogEventBarrierBegin(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1554: (*inctx->begin)(inctx,x,y,addv,mode);
1555: if (inctx->beginandendtogether && inctx->end) {
1556: inctx->inuse = PETSC_FALSE;
1557: (*inctx->end)(inctx,x,y,addv,mode);
1558: }
1559: PetscLogEventBarrierEnd(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1560: CHKMEMQ;
1561: return(0);
1562: }
1564: /* --------------------------------------------------------------------*/
1567: /*@
1568: VecScatterEnd - Ends a generalized scatter from one vector to another. Call
1569: after first calling VecScatterBegin().
1571: Collective on VecScatter and Vec
1573: Input Parameters:
1574: + ctx - scatter context generated by VecScatterCreate()
1575: . x - the vector from which we scatter
1576: . y - the vector to which we scatter
1577: . addv - either ADD_VALUES or INSERT_VALUES.
1578: - mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1579: SCATTER_FORWARD, SCATTER_REVERSE
1581: Level: intermediate
1583: Notes:
1584: If you use SCATTER_REVERSE the arguments x and y should be reversed, from the SCATTER_FORWARD.
1586: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1588: .seealso: VecScatterBegin(), VecScatterCreate()
1589: @*/
1590: PetscErrorCode VecScatterEnd(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1591: {
1598: ctx->inuse = PETSC_FALSE;
1599: if (!ctx->end) return(0);
1600: CHKMEMQ;
1601: if (!ctx->beginandendtogether) {
1602: PetscLogEventBegin(VEC_ScatterEnd,ctx,x,y,0);
1603: (*(ctx)->end)(ctx,x,y,addv,mode);
1604: PetscLogEventEnd(VEC_ScatterEnd,ctx,x,y,0);
1605: }
1606: CHKMEMQ;
1607: return(0);
1608: }
1612: /*@
1613: VecScatterDestroy - Destroys a scatter context created by
1614: VecScatterCreate().
1616: Collective on VecScatter
1618: Input Parameter:
1619: . ctx - the scatter context
1621: Level: intermediate
1623: .seealso: VecScatterCreate(), VecScatterCopy()
1624: @*/
1625: PetscErrorCode VecScatterDestroy(VecScatter ctx)
1626: {
1631: if (--((PetscObject)ctx)->refct > 0) return(0);
1633: /* if memory was published with AMS then destroy it */
1634: PetscObjectDepublish(ctx);
1636: (*ctx->destroy)(ctx);
1637: return(0);
1638: }
1642: /*@
1643: VecScatterCopy - Makes a copy of a scatter context.
1645: Collective on VecScatter
1647: Input Parameter:
1648: . sctx - the scatter context
1650: Output Parameter:
1651: . ctx - the context copy
1653: Level: advanced
1655: .seealso: VecScatterCreate(), VecScatterDestroy()
1656: @*/
1657: PetscErrorCode VecScatterCopy(VecScatter sctx,VecScatter *ctx)
1658: {
1664: if (!sctx->copy) SETERRQ(PETSC_ERR_SUP,"Cannot copy this type");
1665: PetscHeaderCreate(*ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",((PetscObject)sctx)->comm,VecScatterDestroy,VecScatterView);
1666: (*ctx)->to_n = sctx->to_n;
1667: (*ctx)->from_n = sctx->from_n;
1668: (*sctx->copy)(sctx,*ctx);
1669: return(0);
1670: }
1673: /* ------------------------------------------------------------------*/
1676: /*@
1677: VecScatterView - Views a vector scatter context.
1679: Collective on VecScatter
1681: Input Parameters:
1682: + ctx - the scatter context
1683: - viewer - the viewer for displaying the context
1685: Level: intermediate
1687: @*/
1688: PetscErrorCode VecScatterView(VecScatter ctx,PetscViewer viewer)
1689: {
1694: if (!viewer) {
1695: PetscViewerASCIIGetStdout(((PetscObject)ctx)->comm,&viewer);
1696: }
1698: if (ctx->view) {
1699: (*ctx->view)(ctx,viewer);
1700: }
1701: return(0);
1702: }
1706: /*@C
1707: VecScatterRemap - Remaps the "from" and "to" indices in a
1708: vector scatter context. FOR EXPERTS ONLY!
1710: Collective on VecScatter
1712: Input Parameters:
1713: + scat - vector scatter context
1714: . from - remapping for "from" indices (may be PETSC_NULL)
1715: - to - remapping for "to" indices (may be PETSC_NULL)
1717: Level: developer
1719: Notes: In the parallel case the todata is actually the indices
1720: from which the data is TAKEN! The from stuff is where the
1721: data is finally put. This is VERY VERY confusing!
1723: In the sequential case the todata is the indices where the
1724: data is put and the fromdata is where it is taken from.
1725: This is backwards from the paralllel case! CRY! CRY! CRY!
1727: @*/
1728: PetscErrorCode VecScatterRemap(VecScatter scat,PetscInt *rto,PetscInt *rfrom)
1729: {
1730: VecScatter_Seq_General *to,*from;
1731: VecScatter_MPI_General *mto;
1732: PetscInt i;
1739: from = (VecScatter_Seq_General *)scat->fromdata;
1740: mto = (VecScatter_MPI_General *)scat->todata;
1742: if (mto->type == VEC_SCATTER_MPI_TOALL) SETERRQ(PETSC_ERR_ARG_SIZ,"Not for to all scatter");
1744: if (rto) {
1745: if (mto->type == VEC_SCATTER_MPI_GENERAL) {
1746: /* handle off processor parts */
1747: for (i=0; i<mto->starts[mto->n]; i++) {
1748: mto->indices[i] = rto[mto->indices[i]];
1749: }
1750: /* handle local part */
1751: to = &mto->local;
1752: for (i=0; i<to->n; i++) {
1753: to->vslots[i] = rto[to->vslots[i]];
1754: }
1755: } else if (from->type == VEC_SCATTER_SEQ_GENERAL) {
1756: for (i=0; i<from->n; i++) {
1757: from->vslots[i] = rto[from->vslots[i]];
1758: }
1759: } else if (from->type == VEC_SCATTER_SEQ_STRIDE) {
1760: VecScatter_Seq_Stride *sto = (VecScatter_Seq_Stride*)from;
1761:
1762: /* if the remapping is the identity and stride is identity then skip remap */
1763: if (sto->step == 1 && sto->first == 0) {
1764: for (i=0; i<sto->n; i++) {
1765: if (rto[i] != i) {
1766: SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1767: }
1768: }
1769: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1770: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1771: }
1773: if (rfrom) {
1774: SETERRQ(PETSC_ERR_SUP,"Unable to remap the FROM in scatters yet");
1775: }
1777: /*
1778: Mark then vector lengths as unknown because we do not know the
1779: lengths of the remapped vectors
1780: */
1781: scat->from_n = -1;
1782: scat->to_n = -1;
1784: return(0);
1785: }