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: PetscFree4(ctx->todata,((VecScatter_Seq_General*)ctx->todata)->vslots,ctx->fromdata,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
261: PetscHeaderDestroy(ctx);
262: return(0);
263: }
267: PetscErrorCode VecScatterDestroy_SGtoSS(VecScatter ctx)
268: {
269: PetscErrorCode ierr;
272: PetscFree3(ctx->todata,ctx->fromdata,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
273: PetscHeaderDestroy(ctx);
274: return(0);
275: }
279: PetscErrorCode VecScatterDestroy_SStoSG(VecScatter ctx)
280: {
281: PetscErrorCode ierr;
284: PetscFree3(ctx->todata,((VecScatter_Seq_General*)ctx->todata)->vslots,ctx->fromdata);
285: PetscHeaderDestroy(ctx);
286: return(0);
287: }
291: PetscErrorCode VecScatterDestroy_SStoSS(VecScatter ctx)
292: {
296: PetscFree2(ctx->todata,ctx->fromdata);
297: PetscHeaderDestroy(ctx);
298: return(0);
299: }
301: /* -------------------------------------------------------------------------*/
304: PetscErrorCode VecScatterCopy_MPI_ToAll(VecScatter in,VecScatter out)
305: {
306: VecScatter_MPI_ToAll *in_to = (VecScatter_MPI_ToAll*)in->todata,*sto;
307: PetscErrorCode ierr;
308: PetscMPIInt size;
309: PetscInt i;
312: out->begin = in->begin;
313: out->end = in->end;
314: out->copy = in->copy;
315: out->destroy = in->destroy;
316: out->view = in->view;
318: MPI_Comm_size(((PetscObject)out)->comm,&size);
319: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&sto->count);
320: PetscMalloc(size*sizeof(PetscMPIInt),&sto->displx);
321: sto->type = in_to->type;
322: for (i=0; i<size; i++) {
323: sto->count[i] = in_to->count[i];
324: sto->displx[i] = in_to->displx[i];
325: }
326: sto->work1 = 0;
327: sto->work2 = 0;
328: out->todata = (void*)sto;
329: out->fromdata = (void*)0;
330: return(0);
331: }
333: /* --------------------------------------------------------------------------------------*/
334: /*
335: Scatter: sequential general to sequential general
336: */
339: PetscErrorCode VecScatterBegin_SGtoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
340: {
341: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
342: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
343: PetscErrorCode ierr;
344: PetscInt i,n = gen_from->n,*fslots,*tslots;
345: PetscScalar *xv,*yv;
346:
348: VecGetArray(x,&xv);
349: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
350: if (mode & SCATTER_REVERSE){
351: gen_to = (VecScatter_Seq_General*)ctx->fromdata;
352: gen_from = (VecScatter_Seq_General*)ctx->todata;
353: }
354: fslots = gen_from->vslots;
355: tslots = gen_to->vslots;
357: if (addv == INSERT_VALUES) {
358: for (i=0; i<n; i++) {yv[tslots[i]] = xv[fslots[i]];}
359: } else if (addv == ADD_VALUES) {
360: for (i=0; i<n; i++) {yv[tslots[i]] += xv[fslots[i]];}
361: #if !defined(PETSC_USE_COMPLEX)
362: } else if (addv == MAX_VALUES) {
363: for (i=0; i<n; i++) {yv[tslots[i]] = PetscMax(yv[tslots[i]],xv[fslots[i]]);}
364: #endif
365: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
366: VecRestoreArray(x,&xv);
367: if (x != y) {VecRestoreArray(y,&yv);}
368: return(0);
369: }
371: /*
372: Scatter: sequential general to sequential stride 1
373: */
376: PetscErrorCode VecScatterBegin_SGtoSS_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
377: {
378: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
379: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
380: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
381: PetscErrorCode ierr;
382: PetscInt first = gen_to->first;
383: PetscScalar *xv,*yv;
384:
386: VecGetArray(x,&xv);
387: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
388: if (mode & SCATTER_REVERSE){
389: xv += first;
390: if (addv == INSERT_VALUES) {
391: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
392: } else if (addv == ADD_VALUES) {
393: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
394: #if !defined(PETSC_USE_COMPLEX)
395: } else if (addv == MAX_VALUES) {
396: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
397: #endif
398: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
399: } else {
400: yv += first;
401: if (addv == INSERT_VALUES) {
402: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
403: } else if (addv == ADD_VALUES) {
404: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
405: #if !defined(PETSC_USE_COMPLEX)
406: } else if (addv == MAX_VALUES) {
407: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
408: #endif
409: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
410: }
411: VecRestoreArray(x,&xv);
412: if (x != y) {VecRestoreArray(y,&yv);}
413: return(0);
414: }
416: /*
417: Scatter: sequential general to sequential stride
418: */
421: PetscErrorCode VecScatterBegin_SGtoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
422: {
423: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
424: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
425: PetscInt i,n = gen_from->n,*fslots = gen_from->vslots;
426: PetscErrorCode ierr;
427: PetscInt first = gen_to->first,step = gen_to->step;
428: PetscScalar *xv,*yv;
429:
431: VecGetArray(x,&xv);
432: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
434: if (mode & SCATTER_REVERSE){
435: if (addv == INSERT_VALUES) {
436: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
437: } else if (addv == ADD_VALUES) {
438: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
439: #if !defined(PETSC_USE_COMPLEX)
440: } else if (addv == MAX_VALUES) {
441: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
442: #endif
443: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
444: } else {
445: if (addv == INSERT_VALUES) {
446: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
447: } else if (addv == ADD_VALUES) {
448: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
449: #if !defined(PETSC_USE_COMPLEX)
450: } else if (addv == MAX_VALUES) {
451: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
452: #endif
453: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
454: }
455: VecRestoreArray(x,&xv);
456: if (x != y) {VecRestoreArray(y,&yv);}
457: return(0);
458: }
460: /*
461: Scatter: sequential stride 1 to sequential general
462: */
465: PetscErrorCode VecScatterBegin_SStoSG_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
466: {
467: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
468: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
469: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
470: PetscErrorCode ierr;
471: PetscInt first = gen_from->first;
472: PetscScalar *xv,*yv;
473:
475: VecGetArray(x,&xv);
476: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
478: if (mode & SCATTER_REVERSE){
479: yv += first;
480: if (addv == INSERT_VALUES) {
481: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
482: } else if (addv == ADD_VALUES) {
483: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
484: #if !defined(PETSC_USE_COMPLEX)
485: } else if (addv == MAX_VALUES) {
486: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
487: #endif
488: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
489: } else {
490: xv += first;
491: if (addv == INSERT_VALUES) {
492: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
493: } else if (addv == ADD_VALUES) {
494: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
495: #if !defined(PETSC_USE_COMPLEX)
496: } else if (addv == MAX_VALUES) {
497: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
498: #endif
499: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
500: }
501: VecRestoreArray(x,&xv);
502: if (x != y) {VecRestoreArray(y,&yv);}
503: return(0);
504: }
508: /*
509: Scatter: sequential stride to sequential general
510: */
511: PetscErrorCode VecScatterBegin_SStoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
512: {
513: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
514: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
515: PetscInt i,n = gen_from->n,*fslots = gen_to->vslots;
516: PetscErrorCode ierr;
517: PetscInt first = gen_from->first,step = gen_from->step;
518: PetscScalar *xv,*yv;
519:
521: VecGetArray(x,&xv);
522: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
524: if (mode & SCATTER_REVERSE){
525: if (addv == INSERT_VALUES) {
526: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
527: } else if (addv == ADD_VALUES) {
528: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
529: #if !defined(PETSC_USE_COMPLEX)
530: } else if (addv == MAX_VALUES) {
531: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
532: #endif
533: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
534: } else {
535: if (addv == INSERT_VALUES) {
536: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
537: } else if (addv == ADD_VALUES) {
538: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
539: #if !defined(PETSC_USE_COMPLEX)
540: } else if (addv == MAX_VALUES) {
541: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
542: #endif
543: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
544: }
545: VecRestoreArray(x,&xv);
546: if (x != y) {VecRestoreArray(y,&yv);}
547: return(0);
548: }
550: /*
551: Scatter: sequential stride to sequential stride
552: */
555: PetscErrorCode VecScatterBegin_SStoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
556: {
557: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
558: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
559: PetscInt i,n = gen_from->n,to_first = gen_to->first,to_step = gen_to->step;
560: PetscErrorCode ierr;
561: PetscInt from_first = gen_from->first,from_step = gen_from->step;
562: PetscScalar *xv,*yv;
563:
565: VecGetArray(x,&xv);
566: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
568: if (mode & SCATTER_REVERSE){
569: from_first = gen_to->first;
570: to_first = gen_from->first;
571: from_step = gen_to->step;
572: to_step = gen_from->step;
573: }
575: if (addv == INSERT_VALUES) {
576: if (to_step == 1 && from_step == 1) {
577: PetscMemcpy(yv+to_first,xv+from_first,n*sizeof(PetscScalar));
578: } else {
579: for (i=0; i<n; i++) {
580: yv[to_first + i*to_step] = xv[from_first+i*from_step];
581: }
582: }
583: } else if (addv == ADD_VALUES) {
584: if (to_step == 1 && from_step == 1) {
585: yv += to_first; xv += from_first;
586: for (i=0; i<n; i++) {
587: yv[i] += xv[i];
588: }
589: } else {
590: for (i=0; i<n; i++) {
591: yv[to_first + i*to_step] += xv[from_first+i*from_step];
592: }
593: }
594: #if !defined(PETSC_USE_COMPLEX)
595: } else if (addv == MAX_VALUES) {
596: if (to_step == 1 && from_step == 1) {
597: yv += to_first; xv += from_first;
598: for (i=0; i<n; i++) {
599: yv[i] = PetscMax(yv[i],xv[i]);
600: }
601: } else {
602: for (i=0; i<n; i++) {
603: yv[to_first + i*to_step] = PetscMax(yv[to_first + i*to_step],xv[from_first+i*from_step]);
604: }
605: }
606: #endif
607: } else {SETERRQ(PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");}
608: VecRestoreArray(x,&xv);
609: if (x != y) {VecRestoreArray(y,&yv);}
610: return(0);
611: }
613: /* --------------------------------------------------------------------------*/
618: PetscErrorCode VecScatterCopy_SGToSG(VecScatter in,VecScatter out)
619: {
620: PetscErrorCode ierr;
621: VecScatter_Seq_General *in_to = (VecScatter_Seq_General*)in->todata,*out_to = PETSC_NULL;
622: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
623:
625: out->begin = in->begin;
626: out->end = in->end;
627: out->copy = in->copy;
628: out->destroy = in->destroy;
629: out->view = in->view;
631: 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);
632: out_to->n = in_to->n;
633: out_to->type = in_to->type;
634: out_to->nonmatching_computed = PETSC_FALSE;
635: out_to->slots_nonmatching = 0;
636: out_to->is_copy = PETSC_FALSE;
637: PetscMemcpy(out_to->vslots,in_to->vslots,(out_to->n)*sizeof(PetscInt));
639: out_from->n = in_from->n;
640: out_from->type = in_from->type;
641: out_from->nonmatching_computed = PETSC_FALSE;
642: out_from->slots_nonmatching = 0;
643: out_from->is_copy = PETSC_FALSE;
644: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
646: out->todata = (void*)out_to;
647: out->fromdata = (void*)out_from;
648: return(0);
649: }
654: PetscErrorCode VecScatterCopy_SGToSS(VecScatter in,VecScatter out)
655: {
656: PetscErrorCode ierr;
657: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to = PETSC_NULL;
658: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
659:
661: out->begin = in->begin;
662: out->end = in->end;
663: out->copy = in->copy;
664: out->destroy = in->destroy;
665: out->view = in->view;
667: PetscMalloc3(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_General,&out_from,in_from->n,PetscInt,&out_from->vslots);
668: out_to->n = in_to->n;
669: out_to->type = in_to->type;
670: out_to->first = in_to->first;
671: out_to->step = in_to->step;
672: out_to->type = in_to->type;
674: out_from->n = in_from->n;
675: out_from->type = in_from->type;
676: out_from->nonmatching_computed = PETSC_FALSE;
677: out_from->slots_nonmatching = 0;
678: out_from->is_copy = PETSC_FALSE;
679: PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));
681: out->todata = (void*)out_to;
682: out->fromdata = (void*)out_from;
683: return(0);
684: }
686: /* --------------------------------------------------------------------------*/
687: /*
688: Scatter: parallel to sequential vector, sequential strides for both.
689: */
692: PetscErrorCode VecScatterCopy_SStoSS(VecScatter in,VecScatter out)
693: {
694: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to = PETSC_NULL;
695: VecScatter_Seq_Stride *in_from = (VecScatter_Seq_Stride*)in->fromdata,*out_from = PETSC_NULL;
696: PetscErrorCode ierr;
699: out->begin = in->begin;
700: out->end = in->end;
701: out->copy = in->copy;
702: out->destroy = in->destroy;
703: out->view = in->view;
705: PetscMalloc2(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_Stride,&out_from);
706: out_to->n = in_to->n;
707: out_to->type = in_to->type;
708: out_to->first = in_to->first;
709: out_to->step = in_to->step;
710: out_to->type = in_to->type;
711: out_from->n = in_from->n;
712: out_from->type = in_from->type;
713: out_from->first = in_from->first;
714: out_from->step = in_from->step;
715: out_from->type = in_from->type;
716: out->todata = (void*)out_to;
717: out->fromdata = (void*)out_from;
718: return(0);
719: }
721: EXTERN PetscErrorCode VecScatterCreate_PtoS(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,PetscInt,VecScatter);
722: EXTERN PetscErrorCode VecScatterCreate_PtoP(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,VecScatter);
723: EXTERN PetscErrorCode VecScatterCreate_StoP(PetscInt,const PetscInt *,PetscInt,const PetscInt *,Vec,Vec,PetscInt,VecScatter);
725: /* =======================================================================*/
726: #define VEC_SEQ_ID 0
727: #define VEC_MPI_ID 1
729: /*
730: Blocksizes we have optimized scatters for
731: */
733: #define VecScatterOptimizedBS(mbs) ((2 <= mbs && mbs <= 8) || mbs == 12)
735: PetscErrorCode VecScatterCreateEmpty(MPI_Comm comm,VecScatter *newctx)
736: {
737: VecScatter ctx;
741: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
742: ctx->inuse = PETSC_FALSE;
743: ctx->beginandendtogether = PETSC_FALSE;
744: PetscOptionsHasName(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether);
745: if (ctx->beginandendtogether) {
746: PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
747: }
748: PetscOptionsHasName(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether);
749: if (ctx->packtogether) {
750: PetscInfo(ctx,"Pack all messages before sending\n");
751: }
752: *newctx = ctx;
753: return(0);
754: }
758: /*@C
759: VecScatterCreate - Creates a vector scatter context.
761: Collective on Vec
763: Input Parameters:
764: + xin - a vector that defines the shape (parallel data layout of the vector)
765: of vectors from which we scatter
766: . yin - a vector that defines the shape (parallel data layout of the vector)
767: of vectors to which we scatter
768: . ix - the indices of xin to scatter (if PETSC_NULL scatters all values)
769: - iy - the indices of yin to hold results (if PETSC_NULL fills entire vector yin)
771: Output Parameter:
772: . newctx - location to store the new scatter context
774: Options Database Keys: (uses regular MPI_Sends by default)
775: + -vecscatter_ssend - Uses MPI_Ssend_init() instead of MPI_Send_init()
776: . -vecscatter_rsend - use ready receiver mode for MPI sends
777: . -vecscatter_merge - VecScatterBegin() handles all of the communication, VecScatterEnd() is a nop
778: eliminates the chance for overlap of computation and communication
779: . -vecscatter_sendfirst - Posts sends before receives
780: . -vecscatter_packtogether - Pack all messages before sending, receive all messages before unpacking
781: . -vecscatter_alltoall - Uses MPI all to all communication for scatter
782: - -vecscatter_nopack - Avoid packing to work vector when possible (if used with -vecscatter_alltoall then will use MPI_Alltoallw()
784: $
785: $ --When packing is used--
786: $ MPI Datatypes (no packing) sendfirst merge packtogether persistent* -vecscatter_
787: $
788: $ Message passing Send p X X X always
789: $ Ssend p X X X always _ssend
790: $ Rsend p nonsense X X always _rsend
791: $ AlltoAll v or w X nonsense always X nonsense _alltoall
792: $ MPI_Win p nonsense p p nonsense _window
793: $
794: $ -vecscatter_ _nopack _sendfirst _merge _packtogether
795: $
796: $ Since persistent sends and receives require a constant memory address they can only be used when data is packed into the work vector
797: $ because the in and out array may be different for each call to VecScatterBegin/End().
798: $
799: $ p indicates possible, but not implemented. X indicates implemented
800: $
802: Level: intermediate
804: Notes:
805: In calls to VecScatter() you can use different vectors than the xin and
806: yin you used above; BUT they must have the same parallel data layout, for example,
807: they could be obtained from VecDuplicate().
808: A VecScatter context CANNOT be used in two or more simultaneous scatters;
809: that is you cannot call a second VecScatterBegin() with the same scatter
810: context until the VecScatterEnd() has been called on the first VecScatterBegin().
811: In this case a separate VecScatter is needed for each concurrent scatter.
813: Currently the MPI_Send(), MPI_Ssend() and MPI_Rsend() all use PERSISTENT versions.
814: (this unfortunately requires that the same in and out arrays be used for each use, this
815: is why when not using MPI_alltoallw() we always need to pack the input into the work array before sending
816: and unpack upon receeving instead of using MPI datatypes to avoid the packing/unpacking).
818: Concepts: scatter^between vectors
819: Concepts: gather^between vectors
821: .seealso: VecScatterDestroy(), VecScatterCreateToAll(), VecScatterCreateToZero()
822: @*/
823: PetscErrorCode VecScatterCreate(Vec xin,IS ix,Vec yin,IS iy,VecScatter *newctx)
824: {
825: VecScatter ctx;
827: PetscMPIInt size;
828: PetscInt totalv,xin_type = VEC_SEQ_ID,yin_type = VEC_SEQ_ID,*range;
829: MPI_Comm comm,ycomm;
830: PetscTruth ixblock,iyblock,iystride,islocal,cando,flag;
831: IS tix = 0,tiy = 0;
835: /*
836: Determine if the vectors are "parallel", ie. it shares a comm with other processors, or
837: sequential (it does not share a comm). The difference is that parallel vectors treat the
838: index set as providing indices in the global parallel numbering of the vector, with
839: sequential vectors treat the index set as providing indices in the local sequential
840: numbering
841: */
842: PetscObjectGetComm((PetscObject)xin,&comm);
843: MPI_Comm_size(comm,&size);
844: if (size > 1) {xin_type = VEC_MPI_ID;}
846: PetscObjectGetComm((PetscObject)yin,&ycomm);
847: MPI_Comm_size(ycomm,&size);
848: if (size > 1) {comm = ycomm; yin_type = VEC_MPI_ID;}
849:
850: /* generate the Scatter context */
851: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
852: ctx->inuse = PETSC_FALSE;
854: ctx->beginandendtogether = PETSC_FALSE;
855: PetscOptionsHasName(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether);
856: if (ctx->beginandendtogether) {
857: PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
858: }
859: PetscOptionsHasName(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether);
860: if (ctx->packtogether) {
861: PetscInfo(ctx,"Pack all messages before sending\n");
862: }
864: VecGetLocalSize(xin,&ctx->from_n);
865: VecGetLocalSize(yin,&ctx->to_n);
867: /*
868: if ix or iy is not included; assume just grabbing entire vector
869: */
870: if (!ix && xin_type == VEC_SEQ_ID) {
871: ISCreateStride(comm,ctx->from_n,0,1,&ix);
872: tix = ix;
873: } else if (!ix && xin_type == VEC_MPI_ID) {
874: if (yin_type == VEC_MPI_ID) {
875: PetscInt ntmp, low;
876: VecGetLocalSize(xin,&ntmp);
877: VecGetOwnershipRange(xin,&low,PETSC_NULL);
878: ISCreateStride(comm,ntmp,low,1,&ix);
879: } else{
880: PetscInt Ntmp;
881: VecGetSize(xin,&Ntmp);
882: ISCreateStride(comm,Ntmp,0,1,&ix);
883: }
884: tix = ix;
885: } else if (!ix) {
886: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ix not given, but not Seq or MPI vector");
887: }
889: if (!iy && yin_type == VEC_SEQ_ID) {
890: ISCreateStride(comm,ctx->to_n,0,1,&iy);
891: tiy = iy;
892: } else if (!iy && yin_type == VEC_MPI_ID) {
893: if (xin_type == VEC_MPI_ID) {
894: PetscInt ntmp, low;
895: VecGetLocalSize(yin,&ntmp);
896: VecGetOwnershipRange(yin,&low,PETSC_NULL);
897: ISCreateStride(comm,ntmp,low,1,&iy);
898: } else{
899: PetscInt Ntmp;
900: VecGetSize(yin,&Ntmp);
901: ISCreateStride(comm,Ntmp,0,1,&iy);
902: }
903: tiy = iy;
904: } else if (!iy) {
905: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
906: }
908: /* ===========================================================================================================
909: Check for special cases
910: ==========================================================================================================*/
911: /* ---------------------------------------------------------------------------*/
912: if (xin_type == VEC_SEQ_ID && yin_type == VEC_SEQ_ID) {
913: if (((PetscObject)ix)->type == IS_GENERAL && ((PetscObject)iy)->type == IS_GENERAL){
914: PetscInt nx,ny;
915: const PetscInt *idx,*idy;
916: VecScatter_Seq_General *to = PETSC_NULL,*from = PETSC_NULL;
918: ISGetLocalSize(ix,&nx);
919: ISGetLocalSize(iy,&ny);
920: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
921: ISGetIndices(ix,&idx);
922: ISGetIndices(iy,&idy);
923: PetscMalloc4(1,VecScatter_Seq_General,&to,nx,PetscInt,&to->vslots,1,VecScatter_Seq_General,&from,nx,PetscInt,&from->vslots);
924: to->n = nx;
925: #if defined(PETSC_USE_DEBUG)
926: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
927: #endif
928: PetscMemcpy(to->vslots,idy,nx*sizeof(PetscInt));
929: from->n = nx;
930: #if defined(PETSC_USE_DEBUG)
931: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
932: #endif
933: PetscMemcpy(from->vslots,idx,nx*sizeof(PetscInt));
934: to->type = VEC_SCATTER_SEQ_GENERAL;
935: from->type = VEC_SCATTER_SEQ_GENERAL;
936: ctx->todata = (void*)to;
937: ctx->fromdata = (void*)from;
938: ctx->begin = VecScatterBegin_SGtoSG;
939: ctx->end = 0;
940: ctx->destroy = VecScatterDestroy_SGtoSG;
941: ctx->copy = VecScatterCopy_SGToSG;
942: PetscInfo(xin,"Special case: sequential vector general scatter\n");
943: goto functionend;
944: } else if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
945: PetscInt nx,ny,to_first,to_step,from_first,from_step;
946: VecScatter_Seq_Stride *from8 = PETSC_NULL,*to8 = PETSC_NULL;
948: ISGetLocalSize(ix,&nx);
949: ISGetLocalSize(iy,&ny);
950: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
951: ISStrideGetInfo(iy,&to_first,&to_step);
952: ISStrideGetInfo(ix,&from_first,&from_step);
953: PetscMalloc2(1,VecScatter_Seq_Stride,&to8,1,VecScatter_Seq_Stride,&from8);
954: to8->n = nx;
955: to8->first = to_first;
956: to8->step = to_step;
957: from8->n = nx;
958: from8->first = from_first;
959: from8->step = from_step;
960: to8->type = VEC_SCATTER_SEQ_STRIDE;
961: from8->type = VEC_SCATTER_SEQ_STRIDE;
962: ctx->todata = (void*)to8;
963: ctx->fromdata = (void*)from8;
964: ctx->begin = VecScatterBegin_SStoSS;
965: ctx->end = 0;
966: ctx->destroy = VecScatterDestroy_SStoSS;
967: ctx->copy = VecScatterCopy_SStoSS;
968: PetscInfo(xin,"Special case: sequential vector stride to stride\n");
969: goto functionend;
970: } else if (((PetscObject)ix)->type == IS_GENERAL && ((PetscObject)iy)->type == IS_STRIDE){
971: PetscInt nx,ny,first,step;
972: const PetscInt *idx;
973: VecScatter_Seq_General *from9 = PETSC_NULL;
974: VecScatter_Seq_Stride *to9 = PETSC_NULL;
976: ISGetLocalSize(ix,&nx);
977: ISGetIndices(ix,&idx);
978: ISGetLocalSize(iy,&ny);
979: ISStrideGetInfo(iy,&first,&step);
980: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
981: PetscMalloc3(1,VecScatter_Seq_Stride,&to9,1,VecScatter_Seq_General,&from9,nx,PetscInt,&from9->vslots);
982: to9->n = nx;
983: to9->first = first;
984: to9->step = step;
985: from9->n = nx;
986: #if defined(PETSC_USE_DEBUG)
987: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
988: #endif
989: PetscMemcpy(from9->vslots,idx,nx*sizeof(PetscInt));
990: ctx->todata = (void*)to9; ctx->fromdata = (void*)from9;
991: if (step == 1) ctx->begin = VecScatterBegin_SGtoSS_Stride1;
992: else ctx->begin = VecScatterBegin_SGtoSS;
993: ctx->destroy = VecScatterDestroy_SGtoSS;
994: ctx->end = 0;
995: ctx->copy = VecScatterCopy_SGToSS;
996: to9->type = VEC_SCATTER_SEQ_STRIDE;
997: from9->type = VEC_SCATTER_SEQ_GENERAL;
998: PetscInfo(xin,"Special case: sequential vector general to stride\n");
999: goto functionend;
1000: } else if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_GENERAL){
1001: PetscInt nx,ny,first,step;
1002: const PetscInt *idy;
1003: VecScatter_Seq_General *to10 = PETSC_NULL;
1004: VecScatter_Seq_Stride *from10 = PETSC_NULL;
1006: ISGetLocalSize(ix,&nx);
1007: ISGetIndices(iy,&idy);
1008: ISGetLocalSize(iy,&ny);
1009: ISStrideGetInfo(ix,&first,&step);
1010: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1011: PetscMalloc3(1,VecScatter_Seq_General,&to10,nx,PetscInt,&to10->vslots,1,VecScatter_Seq_Stride,&from10);
1012: from10->n = nx;
1013: from10->first = first;
1014: from10->step = step;
1015: to10->n = nx;
1016: #if defined(PETSC_USE_DEBUG)
1017: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1018: #endif
1019: PetscMemcpy(to10->vslots,idy,nx*sizeof(PetscInt));
1020: ctx->todata = (void*)to10;
1021: ctx->fromdata = (void*)from10;
1022: if (step == 1) ctx->begin = VecScatterBegin_SStoSG_Stride1;
1023: else ctx->begin = VecScatterBegin_SStoSG;
1024: ctx->destroy = VecScatterDestroy_SStoSG;
1025: ctx->end = 0;
1026: ctx->copy = 0;
1027: to10->type = VEC_SCATTER_SEQ_GENERAL;
1028: from10->type = VEC_SCATTER_SEQ_STRIDE;
1029: PetscInfo(xin,"Special case: sequential vector stride to general\n");
1030: goto functionend;
1031: } else {
1032: PetscInt nx,ny;
1033: const PetscInt *idx,*idy;
1034: VecScatter_Seq_General *to11 = PETSC_NULL,*from11 = PETSC_NULL;
1035: PetscTruth idnx,idny;
1037: ISGetLocalSize(ix,&nx);
1038: ISGetLocalSize(iy,&ny);
1039: if (nx != ny) SETERRQ2(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match, in %D out %D",nx,ny);
1041: ISIdentity(ix,&idnx);
1042: ISIdentity(iy,&idny);
1043: if (idnx && idny) {
1044: VecScatter_Seq_Stride *to13 = PETSC_NULL,*from13 = PETSC_NULL;
1045: PetscMalloc2(1,VecScatter_Seq_Stride,&to13,1,VecScatter_Seq_Stride,&from13);
1046: to13->n = nx;
1047: to13->first = 0;
1048: to13->step = 1;
1049: from13->n = nx;
1050: from13->first = 0;
1051: from13->step = 1;
1052: to13->type = VEC_SCATTER_SEQ_STRIDE;
1053: from13->type = VEC_SCATTER_SEQ_STRIDE;
1054: ctx->todata = (void*)to13;
1055: ctx->fromdata = (void*)from13;
1056: ctx->begin = VecScatterBegin_SStoSS;
1057: ctx->end = 0;
1058: ctx->destroy = VecScatterDestroy_SStoSS;
1059: ctx->copy = VecScatterCopy_SStoSS;
1060: PetscInfo(xin,"Special case: sequential copy\n");
1061: goto functionend;
1062: }
1064: ISGetIndices(iy,&idy);
1065: ISGetIndices(ix,&idx);
1066: PetscMalloc4(1,VecScatter_Seq_General,&to11,nx,PetscInt,&to11->vslots,1,VecScatter_Seq_General,&from11,nx,PetscInt,&from11->vslots);
1067: to11->n = nx;
1068: #if defined(PETSC_USE_DEBUG)
1069: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1070: #endif
1071: PetscMemcpy(to11->vslots,idy,nx*sizeof(PetscInt));
1072: from11->n = nx;
1073: #if defined(PETSC_USE_DEBUG)
1074: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
1075: #endif
1076: PetscMemcpy(from11->vslots,idx,nx*sizeof(PetscInt));
1077: to11->type = VEC_SCATTER_SEQ_GENERAL;
1078: from11->type = VEC_SCATTER_SEQ_GENERAL;
1079: ctx->todata = (void*)to11;
1080: ctx->fromdata = (void*)from11;
1081: ctx->begin = VecScatterBegin_SGtoSG;
1082: ctx->end = 0;
1083: ctx->destroy = VecScatterDestroy_SGtoSG;
1084: ctx->copy = VecScatterCopy_SGToSG;
1085: ISRestoreIndices(ix,&idx);
1086: ISRestoreIndices(iy,&idy);
1087: PetscInfo(xin,"Sequential vector scatter with block indices\n");
1088: goto functionend;
1089: }
1090: }
1091: /* ---------------------------------------------------------------------------*/
1092: if (xin_type == VEC_MPI_ID && yin_type == VEC_SEQ_ID) {
1094: /* ===========================================================================================================
1095: Check for special cases
1096: ==========================================================================================================*/
1097: islocal = PETSC_FALSE;
1098: /* special case extracting (subset of) local portion */
1099: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1100: PetscInt nx,ny,to_first,to_step,from_first,from_step;
1101: PetscInt start,end;
1102: VecScatter_Seq_Stride *from12 = PETSC_NULL,*to12 = PETSC_NULL;
1104: VecGetOwnershipRange(xin,&start,&end);
1105: ISGetLocalSize(ix,&nx);
1106: ISStrideGetInfo(ix,&from_first,&from_step);
1107: ISGetLocalSize(iy,&ny);
1108: ISStrideGetInfo(iy,&to_first,&to_step);
1109: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1110: if (ix->min >= start && ix->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1111: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1112: if (cando) {
1113: PetscMalloc2(1,VecScatter_Seq_Stride,&to12,1,VecScatter_Seq_Stride,&from12);
1114: to12->n = nx;
1115: to12->first = to_first;
1116: to12->step = to_step;
1117: from12->n = nx;
1118: from12->first = from_first-start;
1119: from12->step = from_step;
1120: to12->type = VEC_SCATTER_SEQ_STRIDE;
1121: from12->type = VEC_SCATTER_SEQ_STRIDE;
1122: ctx->todata = (void*)to12;
1123: ctx->fromdata = (void*)from12;
1124: ctx->begin = VecScatterBegin_SStoSS;
1125: ctx->end = 0;
1126: ctx->destroy = VecScatterDestroy_SStoSS;
1127: ctx->copy = VecScatterCopy_SStoSS;
1128: PetscInfo(xin,"Special case: processors only getting local values\n");
1129: goto functionend;
1130: }
1131: } else {
1132: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1133: }
1135: /* test for special case of all processors getting entire vector */
1136: /* contains check that PetscMPIInt can handle the sizes needed */
1137: totalv = 0;
1138: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1139: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1140: PetscMPIInt *count = PETSC_NULL,*displx;
1141: VecScatter_MPI_ToAll *sto = PETSC_NULL;
1143: ISGetLocalSize(ix,&nx);
1144: ISStrideGetInfo(ix,&from_first,&from_step);
1145: ISGetLocalSize(iy,&ny);
1146: ISStrideGetInfo(iy,&to_first,&to_step);
1147: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1148: VecGetSize(xin,&N);
1149: if (nx != N) {
1150: totalv = 0;
1151: } else if (from_first == 0 && from_step == 1 && from_first == to_first && from_step == to_step){
1152: totalv = 1;
1153: } else totalv = 0;
1154: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1156: #if defined(PETSC_USE_64BIT_INDICES)
1157: if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1158: #else
1159: if (cando) {
1160: #endif
1161: MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1162: PetscMalloc(size*sizeof(PetscMPIInt),&displx);
1163: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1164: range = xin->map->range;
1165: for (i=0; i<size; i++) {
1166: count[i] = PetscMPIIntCast(range[i+1] - range[i]);
1167: displx[i] = PetscMPIIntCast(range[i]);
1168: }
1169: sto->count = count;
1170: sto->displx = displx;
1171: sto->work1 = 0;
1172: sto->work2 = 0;
1173: sto->type = VEC_SCATTER_MPI_TOALL;
1174: ctx->todata = (void*)sto;
1175: ctx->fromdata = 0;
1176: ctx->begin = VecScatterBegin_MPI_ToAll;
1177: ctx->end = 0;
1178: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1179: ctx->copy = VecScatterCopy_MPI_ToAll;
1180: PetscInfo(xin,"Special case: all processors get entire parallel vector\n");
1181: goto functionend;
1182: }
1183: } else {
1184: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1185: }
1187: /* test for special case of processor 0 getting entire vector */
1188: /* contains check that PetscMPIInt can handle the sizes needed */
1189: totalv = 0;
1190: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1191: PetscInt i,nx,ny,to_first,to_step,from_first,from_step,N;
1192: PetscMPIInt rank,*count = PETSC_NULL,*displx;
1193: VecScatter_MPI_ToAll *sto = PETSC_NULL;
1195: PetscObjectGetComm((PetscObject)xin,&comm);
1196: MPI_Comm_rank(comm,&rank);
1197: ISGetLocalSize(ix,&nx);
1198: ISStrideGetInfo(ix,&from_first,&from_step);
1199: ISGetLocalSize(iy,&ny);
1200: ISStrideGetInfo(iy,&to_first,&to_step);
1201: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1202: if (!rank) {
1203: VecGetSize(xin,&N);
1204: if (nx != N) {
1205: totalv = 0;
1206: } else if (from_first == 0 && from_step == 1 &&
1207: from_first == to_first && from_step == to_step){
1208: totalv = 1;
1209: } else totalv = 0;
1210: } else {
1211: if (!nx) totalv = 1;
1212: else totalv = 0;
1213: }
1214: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1216: #if defined(PETSC_USE_64BIT_INDICES)
1217: if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1218: #else
1219: if (cando) {
1220: #endif
1221: MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1222: PetscMalloc(size*sizeof(PetscMPIInt),&displx);
1223: PetscMalloc2(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count);
1224: range = xin->map->range;
1225: for (i=0; i<size; i++) {
1226: count[i] = PetscMPIIntCast(range[i+1] - range[i]);
1227: displx[i] = PetscMPIIntCast(range[i]);
1228: }
1229: sto->count = count;
1230: sto->displx = displx;
1231: sto->work1 = 0;
1232: sto->work2 = 0;
1233: sto->type = VEC_SCATTER_MPI_TOONE;
1234: ctx->todata = (void*)sto;
1235: ctx->fromdata = 0;
1236: ctx->begin = VecScatterBegin_MPI_ToOne;
1237: ctx->end = 0;
1238: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1239: ctx->copy = VecScatterCopy_MPI_ToAll;
1240: PetscInfo(xin,"Special case: processor zero gets entire parallel vector, rest get none\n");
1241: goto functionend;
1242: }
1243: } else {
1244: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1245: }
1247: ISBlock(ix,&ixblock);
1248: ISBlock(iy,&iyblock);
1249: ISStride(iy,&iystride);
1250: if (ixblock) {
1251: /* special case block to block */
1252: if (iyblock) {
1253: PetscInt nx,ny,bsx,bsy;
1254: const PetscInt *idx,*idy;
1255: ISBlockGetBlockSize(iy,&bsy);
1256: ISBlockGetBlockSize(ix,&bsx);
1257: if (bsx == bsy && VecScatterOptimizedBS(bsx)) {
1258: ISBlockGetLocalSize(ix,&nx);
1259: ISBlockGetIndices(ix,&idx);
1260: ISBlockGetLocalSize(iy,&ny);
1261: ISBlockGetIndices(iy,&idy);
1262: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1263: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,bsx,ctx);
1264: ISBlockRestoreIndices(ix,&idx);
1265: ISBlockRestoreIndices(iy,&idy);
1266: PetscInfo(xin,"Special case: blocked indices\n");
1267: goto functionend;
1268: }
1269: /* special case block to stride */
1270: } else if (iystride) {
1271: PetscInt ystart,ystride,ysize,bsx;
1272: ISStrideGetInfo(iy,&ystart,&ystride);
1273: ISGetLocalSize(iy,&ysize);
1274: ISBlockGetBlockSize(ix,&bsx);
1275: /* see if stride index set is equivalent to block index set */
1276: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1277: PetscInt nx,il,*idy;
1278: const PetscInt *idx;
1279: ISBlockGetLocalSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1280: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1281: PetscMalloc(nx*sizeof(PetscInt),&idy);
1282: if (nx) {
1283: idy[0] = ystart;
1284: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1285: }
1286: VecScatterCreate_PtoS(nx,idx,nx,idy,xin,yin,bsx,ctx);
1287: PetscFree(idy);
1288: ISBlockRestoreIndices(ix,&idx);
1289: PetscInfo(xin,"Special case: blocked indices to stride\n");
1290: goto functionend;
1291: }
1292: }
1293: }
1294: /* left over general case */
1295: {
1296: PetscInt nx,ny;
1297: const PetscInt *idx,*idy;
1298: ISGetLocalSize(ix,&nx);
1299: ISGetIndices(ix,&idx);
1300: ISGetLocalSize(iy,&ny);
1301: ISGetIndices(iy,&idy);
1302: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1303: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,1,ctx);
1304: ISRestoreIndices(ix,&idx);
1305: ISRestoreIndices(iy,&idy);
1306: PetscInfo(xin,"General case: MPI to Seq\n");
1307: goto functionend;
1308: }
1309: }
1310: /* ---------------------------------------------------------------------------*/
1311: if (xin_type == VEC_SEQ_ID && yin_type == VEC_MPI_ID) {
1312: /* ===========================================================================================================
1313: Check for special cases
1314: ==========================================================================================================*/
1315: /* special case local copy portion */
1316: islocal = PETSC_FALSE;
1317: if (((PetscObject)ix)->type == IS_STRIDE && ((PetscObject)iy)->type == IS_STRIDE){
1318: PetscInt nx,ny,to_first,to_step,from_step,start,end,from_first;
1319: VecScatter_Seq_Stride *from = PETSC_NULL,*to = PETSC_NULL;
1321: VecGetOwnershipRange(yin,&start,&end);
1322: ISGetLocalSize(ix,&nx);
1323: ISStrideGetInfo(ix,&from_first,&from_step);
1324: ISGetLocalSize(iy,&ny);
1325: ISStrideGetInfo(iy,&to_first,&to_step);
1326: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1327: if (iy->min >= start && iy->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1328: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1329: if (cando) {
1330: PetscMalloc2(1,VecScatter_Seq_Stride,&to,1,VecScatter_Seq_Stride,&from);
1331: to->n = nx;
1332: to->first = to_first-start;
1333: to->step = to_step;
1334: from->n = nx;
1335: from->first = from_first;
1336: from->step = from_step;
1337: to->type = VEC_SCATTER_SEQ_STRIDE;
1338: from->type = VEC_SCATTER_SEQ_STRIDE;
1339: ctx->todata = (void*)to;
1340: ctx->fromdata = (void*)from;
1341: ctx->begin = VecScatterBegin_SStoSS;
1342: ctx->end = 0;
1343: ctx->destroy = VecScatterDestroy_SStoSS;
1344: ctx->copy = VecScatterCopy_SStoSS;
1345: PetscInfo(xin,"Special case: sequential stride to MPI stride\n");
1346: goto functionend;
1347: }
1348: } else {
1349: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1350: }
1351: /* special case block to stride */
1352: if (((PetscObject)ix)->type == IS_BLOCK && ((PetscObject)iy)->type == IS_STRIDE){
1353: PetscInt ystart,ystride,ysize,bsx;
1354: ISStrideGetInfo(iy,&ystart,&ystride);
1355: ISGetLocalSize(iy,&ysize);
1356: ISBlockGetBlockSize(ix,&bsx);
1357: /* see if stride index set is equivalent to block index set */
1358: if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1359: PetscInt nx,il,*idy;
1360: const PetscInt *idx;
1361: ISBlockGetLocalSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1362: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1363: PetscMalloc(nx*sizeof(PetscInt),&idy);
1364: if (nx) {
1365: idy[0] = ystart;
1366: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1367: }
1368: VecScatterCreate_StoP(nx,idx,nx,idy,xin,yin,bsx,ctx);
1369: PetscFree(idy);
1370: ISBlockRestoreIndices(ix,&idx);
1371: PetscInfo(xin,"Special case: Blocked indices to stride\n");
1372: goto functionend;
1373: }
1374: }
1376: /* general case */
1377: {
1378: PetscInt nx,ny;
1379: const PetscInt *idx,*idy;
1380: ISGetLocalSize(ix,&nx);
1381: ISGetIndices(ix,&idx);
1382: ISGetLocalSize(iy,&ny);
1383: ISGetIndices(iy,&idy);
1384: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1385: VecScatterCreate_StoP(nx,idx,ny,idy,xin,yin,1,ctx);
1386: ISRestoreIndices(ix,&idx);
1387: ISRestoreIndices(iy,&idy);
1388: PetscInfo(xin,"General case: Seq to MPI\n");
1389: goto functionend;
1390: }
1391: }
1392: /* ---------------------------------------------------------------------------*/
1393: if (xin_type == VEC_MPI_ID && yin_type == VEC_MPI_ID) {
1394: /* no special cases for now */
1395: PetscInt nx,ny;
1396: const PetscInt *idx,*idy;
1397: ISGetLocalSize(ix,&nx);
1398: ISGetIndices(ix,&idx);
1399: ISGetLocalSize(iy,&ny);
1400: ISGetIndices(iy,&idy);
1401: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1402: VecScatterCreate_PtoP(nx,idx,ny,idy,xin,yin,ctx);
1403: ISRestoreIndices(ix,&idx);
1404: ISRestoreIndices(iy,&idy);
1405: PetscInfo(xin,"General case: MPI to MPI\n");
1406: goto functionend;
1407: }
1409: functionend:
1410: *newctx = ctx;
1411: if (tix) {ISDestroy(tix);}
1412: if (tiy) {ISDestroy(tiy);}
1413: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view_info",&flag);
1414: if (flag) {
1415: PetscViewer viewer;
1416: PetscViewerASCIIGetStdout(comm,&viewer);
1417: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1418: VecScatterView(ctx,viewer);
1419: PetscViewerPopFormat(viewer);
1420: }
1421: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view",&flag);
1422: if (flag) {
1423: PetscViewer viewer;
1424: PetscViewerASCIIGetStdout(comm,&viewer);
1425: VecScatterView(ctx,viewer);
1426: }
1427: return(0);
1428: }
1430: /* ------------------------------------------------------------------*/
1433: /*@
1434: VecScatterGetMerged - Returns true if the scatter is completed in the VecScatterBegin()
1435: and the VecScatterEnd() does nothing
1437: Not Collective
1439: Input Parameter:
1440: . ctx - scatter context created with VecScatterCreate()
1442: Output Parameter:
1443: . flg - PETSC_TRUE if the VecScatterBegin/End() are all done during the VecScatterBegin()
1445: Level: developer
1447: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1448: @*/
1449: PetscErrorCode VecScatterGetMerged(VecScatter ctx,PetscTruth *flg)
1450: {
1453: *flg = ctx->beginandendtogether;
1454: return(0);
1455: }
1459: /*@
1460: VecScatterBegin - Begins a generalized scatter from one vector to
1461: another. Complete the scattering phase with VecScatterEnd().
1463: Collective on VecScatter and Vec
1465: Input Parameters:
1466: + inctx - scatter context generated by VecScatterCreate()
1467: . x - the vector from which we scatter
1468: . y - the vector to which we scatter
1469: . addv - either ADD_VALUES or INSERT_VALUES, with INSERT_VALUES mode any location
1470: not scattered to retains its old value; i.e. the vector is NOT first zeroed.
1471: - mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1472: SCATTER_FORWARD or SCATTER_REVERSE
1475: Level: intermediate
1477: Options Database: See VecScatterCreate()
1479: Notes:
1480: The vectors x and y need not be the same vectors used in the call
1481: to VecScatterCreate(), but x must have the same parallel data layout
1482: as that passed in as the x to VecScatterCreate(), similarly for the y.
1483: Most likely they have been obtained from VecDuplicate().
1485: You cannot change the values in the input vector between the calls to VecScatterBegin()
1486: and VecScatterEnd().
1488: If you use SCATTER_REVERSE the two arguments x and y should be reversed, from
1489: the SCATTER_FORWARD.
1490:
1491: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1493: This scatter is far more general than the conventional
1494: scatter, since it can be a gather or a scatter or a combination,
1495: depending on the indices ix and iy. If x is a parallel vector and y
1496: is sequential, VecScatterBegin() can serve to gather values to a
1497: single processor. Similarly, if y is parallel and x sequential, the
1498: routine can scatter from one processor to many processors.
1500: Concepts: scatter^between vectors
1501: Concepts: gather^between vectors
1503: .seealso: VecScatterCreate(), VecScatterEnd()
1504: @*/
1505: PetscErrorCode VecScatterBegin(VecScatter inctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1506: {
1508: #if defined(PETSC_USE_DEBUG)
1509: PetscInt to_n,from_n;
1510: #endif
1516: if (inctx->inuse) SETERRQ(PETSC_ERR_ARG_WRONGSTATE," Scatter ctx already in use");
1517: CHKMEMQ;
1519: #if defined(PETSC_USE_DEBUG)
1520: /*
1521: Error checking to make sure these vectors match the vectors used
1522: to create the vector scatter context. -1 in the from_n and to_n indicate the
1523: vector lengths are unknown (for example with mapped scatters) and thus
1524: no error checking is performed.
1525: */
1526: if (inctx->from_n >= 0 && inctx->to_n >= 0) {
1527: VecGetLocalSize(x,&from_n);
1528: VecGetLocalSize(y,&to_n);
1529: if (mode & SCATTER_REVERSE) {
1530: 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);
1531: 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);
1532: } else {
1533: 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);
1534: 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);
1535: }
1536: }
1537: #endif
1539: inctx->inuse = PETSC_TRUE;
1540: PetscLogEventBarrierBegin(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1541: (*inctx->begin)(inctx,x,y,addv,mode);
1542: if (inctx->beginandendtogether && inctx->end) {
1543: inctx->inuse = PETSC_FALSE;
1544: (*inctx->end)(inctx,x,y,addv,mode);
1545: }
1546: PetscLogEventBarrierEnd(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1547: CHKMEMQ;
1548: return(0);
1549: }
1551: /* --------------------------------------------------------------------*/
1554: /*@
1555: VecScatterEnd - Ends a generalized scatter from one vector to another. Call
1556: after first calling VecScatterBegin().
1558: Collective on VecScatter and Vec
1560: Input Parameters:
1561: + ctx - scatter context generated by VecScatterCreate()
1562: . x - the vector from which we scatter
1563: . y - the vector to which we scatter
1564: . addv - either ADD_VALUES or INSERT_VALUES.
1565: - mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1566: SCATTER_FORWARD, SCATTER_REVERSE
1568: Level: intermediate
1570: Notes:
1571: If you use SCATTER_REVERSE the arguments x and y should be reversed, from the SCATTER_FORWARD.
1573: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1575: .seealso: VecScatterBegin(), VecScatterCreate()
1576: @*/
1577: PetscErrorCode VecScatterEnd(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1578: {
1585: ctx->inuse = PETSC_FALSE;
1586: if (!ctx->end) return(0);
1587: CHKMEMQ;
1588: if (!ctx->beginandendtogether) {
1589: PetscLogEventBegin(VEC_ScatterEnd,ctx,x,y,0);
1590: (*(ctx)->end)(ctx,x,y,addv,mode);
1591: PetscLogEventEnd(VEC_ScatterEnd,ctx,x,y,0);
1592: }
1593: CHKMEMQ;
1594: return(0);
1595: }
1599: /*@
1600: VecScatterDestroy - Destroys a scatter context created by
1601: VecScatterCreate().
1603: Collective on VecScatter
1605: Input Parameter:
1606: . ctx - the scatter context
1608: Level: intermediate
1610: .seealso: VecScatterCreate(), VecScatterCopy()
1611: @*/
1612: PetscErrorCode VecScatterDestroy(VecScatter ctx)
1613: {
1618: if (--((PetscObject)ctx)->refct > 0) return(0);
1620: /* if memory was published with AMS then destroy it */
1621: PetscObjectDepublish(ctx);
1623: (*ctx->destroy)(ctx);
1624: return(0);
1625: }
1629: /*@
1630: VecScatterCopy - Makes a copy of a scatter context.
1632: Collective on VecScatter
1634: Input Parameter:
1635: . sctx - the scatter context
1637: Output Parameter:
1638: . ctx - the context copy
1640: Level: advanced
1642: .seealso: VecScatterCreate(), VecScatterDestroy()
1643: @*/
1644: PetscErrorCode VecScatterCopy(VecScatter sctx,VecScatter *ctx)
1645: {
1651: if (!sctx->copy) SETERRQ(PETSC_ERR_SUP,"Cannot copy this type");
1652: PetscHeaderCreate(*ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",((PetscObject)sctx)->comm,VecScatterDestroy,VecScatterView);
1653: (*ctx)->to_n = sctx->to_n;
1654: (*ctx)->from_n = sctx->from_n;
1655: (*sctx->copy)(sctx,*ctx);
1656: return(0);
1657: }
1660: /* ------------------------------------------------------------------*/
1663: /*@
1664: VecScatterView - Views a vector scatter context.
1666: Collective on VecScatter
1668: Input Parameters:
1669: + ctx - the scatter context
1670: - viewer - the viewer for displaying the context
1672: Level: intermediate
1674: @*/
1675: PetscErrorCode VecScatterView(VecScatter ctx,PetscViewer viewer)
1676: {
1681: if (!viewer) {
1682: PetscViewerASCIIGetStdout(((PetscObject)ctx)->comm,&viewer);
1683: }
1685: if (ctx->view) {
1686: (*ctx->view)(ctx,viewer);
1687: }
1688: return(0);
1689: }
1693: /*@
1694: VecScatterRemap - Remaps the "from" and "to" indices in a
1695: vector scatter context. FOR EXPERTS ONLY!
1697: Collective on VecScatter
1699: Input Parameters:
1700: + scat - vector scatter context
1701: . from - remapping for "from" indices (may be PETSC_NULL)
1702: - to - remapping for "to" indices (may be PETSC_NULL)
1704: Level: developer
1706: Notes: In the parallel case the todata is actually the indices
1707: from which the data is TAKEN! The from stuff is where the
1708: data is finally put. This is VERY VERY confusing!
1710: In the sequential case the todata is the indices where the
1711: data is put and the fromdata is where it is taken from.
1712: This is backwards from the paralllel case! CRY! CRY! CRY!
1714: @*/
1715: PetscErrorCode VecScatterRemap(VecScatter scat,PetscInt *rto,PetscInt *rfrom)
1716: {
1717: VecScatter_Seq_General *to,*from;
1718: VecScatter_MPI_General *mto;
1719: PetscInt i;
1726: from = (VecScatter_Seq_General *)scat->fromdata;
1727: mto = (VecScatter_MPI_General *)scat->todata;
1729: if (mto->type == VEC_SCATTER_MPI_TOALL) SETERRQ(PETSC_ERR_ARG_SIZ,"Not for to all scatter");
1731: if (rto) {
1732: if (mto->type == VEC_SCATTER_MPI_GENERAL) {
1733: /* handle off processor parts */
1734: for (i=0; i<mto->starts[mto->n]; i++) {
1735: mto->indices[i] = rto[mto->indices[i]];
1736: }
1737: /* handle local part */
1738: to = &mto->local;
1739: for (i=0; i<to->n; i++) {
1740: to->vslots[i] = rto[to->vslots[i]];
1741: }
1742: } else if (from->type == VEC_SCATTER_SEQ_GENERAL) {
1743: for (i=0; i<from->n; i++) {
1744: from->vslots[i] = rto[from->vslots[i]];
1745: }
1746: } else if (from->type == VEC_SCATTER_SEQ_STRIDE) {
1747: VecScatter_Seq_Stride *sto = (VecScatter_Seq_Stride*)from;
1748:
1749: /* if the remapping is the identity and stride is identity then skip remap */
1750: if (sto->step == 1 && sto->first == 0) {
1751: for (i=0; i<sto->n; i++) {
1752: if (rto[i] != i) {
1753: SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1754: }
1755: }
1756: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1757: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1758: }
1760: if (rfrom) {
1761: SETERRQ(PETSC_ERR_SUP,"Unable to remap the FROM in scatters yet");
1762: }
1764: /*
1765: Mark then vector lengths as unknown because we do not know the
1766: lengths of the remapped vectors
1767: */
1768: scat->from_n = -1;
1769: scat->to_n = -1;
1771: return(0);
1772: }