Actual source code: vscat.c
1: /*$Id: vscat.c,v 1.173 2001/08/10 03:29:59 bsmith Exp $*/
3: /*
4: Code for creating scatters between vectors. This file
5: includes the code for scattering between sequential vectors and
6: some special cases for parallel scatters.
7: */
9: #include src/vec/is/isimpl.h
10: #include src/vec/vecimpl.h
12: /* Logging support */
13: int VEC_SCATTER_COOKIE;
15: /*
16: Checks if any indices are less than zero and generates an error
17: */
18: static int VecScatterCheckIndices_Private(int nmax,int n,int *idx)
19: {
20: int i;
23: for (i=0; i<n; i++) {
24: if (idx[i] < 0) SETERRQ2(1,"Negative index %d at %d location",idx[i],i);
25: if (idx[i] >= nmax) SETERRQ3(1,"Index %d at %d location greater than max %d",idx[i],i,nmax);
26: }
27: return(0);
28: }
30: /*
31: This is special scatter code for when the entire parallel vector is
32: copied to each processor.
34: This code was written by Cameron Cooper, Occidental College, Fall 1995,
35: will working at ANL as a SERS student.
36: */
37: int VecScatterBegin_MPI_ToAll(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
38: {
39: int ierr,yy_n,xx_n,*range;
40: PetscScalar *xv,*yv;
41: PetscMap map;
44: VecGetLocalSize(y,&yy_n);
45: VecGetLocalSize(x,&xx_n);
46: VecGetArray(y,&yv);
47: if (x != y) {VecGetArray(x,&xv);}
49: if (mode & SCATTER_REVERSE) {
50: PetscScalar *xvt,*xvt2;
51: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
52: int i;
54: if (addv == INSERT_VALUES) {
55: int rstart,rend;
56: /*
57: copy the correct part of the local vector into the local storage of
58: the MPI one. Note: this operation only makes sense if all the local
59: vectors have the same values
60: */
61: VecGetOwnershipRange(y,&rstart,&rend);
62: PetscMemcpy(yv,xv+rstart,yy_n*sizeof(PetscScalar));
63: } else {
64: MPI_Comm comm;
65: int rank;
66: PetscObjectGetComm((PetscObject)y,&comm);
67: MPI_Comm_rank(comm,&rank);
68: if (scat->work1) xvt = scat->work1;
69: else {
70: ierr = PetscMalloc((xx_n+1)*sizeof(PetscScalar),&xvt);
71: scat->work1 = xvt;
72: PetscLogObjectMemory(ctx,xx_n*sizeof(PetscScalar));
73: }
74: if (!rank) { /* I am the zeroth processor, values are accumulated here */
75: if (scat->work2) xvt2 = scat->work2;
76: else {
77: ierr = PetscMalloc((xx_n+1)*sizeof(PetscScalar),& xvt2);
78: scat->work2 = xvt2;
79: PetscLogObjectMemory(ctx,xx_n*sizeof(PetscScalar));
80: }
81: VecGetPetscMap(y,&map);
82: PetscMapGetGlobalRange(map,&range);
83: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,xvt2,scat->count,range,MPIU_SCALAR,0,ctx->comm);
84: #if defined(PETSC_USE_COMPLEX)
85: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,ctx->comm);
86: #else
87: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,ctx->comm);
88: #endif
89: if (addv == ADD_VALUES) {
90: for (i=0; i<xx_n; i++) {
91: xvt[i] += xvt2[i];
92: }
93: #if !defined(PETSC_USE_COMPLEX)
94: } else if (addv == MAX_VALUES) {
95: for (i=0; i<xx_n; i++) {
96: xvt[i] = PetscMax(xvt[i],xvt2[i]);
97: }
98: #endif
99: } else {SETERRQ(1,"Wrong insert option");}
100: MPI_Scatterv(xvt,scat->count,map->range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
101: } else {
102: VecGetPetscMap(y,&map);
103: PetscMapGetGlobalRange(map,&range);
104: MPI_Gatherv(yv,yy_n,MPIU_SCALAR,0, 0,0,MPIU_SCALAR,0,ctx->comm);
105: #if defined(PETSC_USE_COMPLEX)
106: MPI_Reduce(xv,xvt,2*xx_n,MPIU_REAL,MPI_SUM,0,ctx->comm);
107: #else
108: MPI_Reduce(xv,xvt,xx_n,MPIU_SCALAR,MPI_SUM,0,ctx->comm);
109: #endif
110: MPI_Scatterv(0,scat->count,range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
111: }
112: }
113: } else {
114: PetscScalar *yvt;
115: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
116: int i;
118: VecGetPetscMap(x,&map);
119: PetscMapGetGlobalRange(map,&range);
120: if (addv == INSERT_VALUES) {
121: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,range,MPIU_SCALAR,ctx->comm);
122: } else {
123: if (scat->work1) yvt = scat->work1;
124: else {
125: ierr = PetscMalloc((yy_n+1)*sizeof(PetscScalar),&yvt);
126: scat->work1 = yvt;
127: PetscLogObjectMemory(ctx,yy_n*sizeof(PetscScalar));
128: }
129: MPI_Allgatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,map->range,MPIU_SCALAR,ctx->comm);
130: if (addv == ADD_VALUES){
131: for (i=0; i<yy_n; i++) {
132: yv[i] += yvt[i];
133: }
134: #if !defined(PETSC_USE_COMPLEX)
135: } else if (addv == MAX_VALUES) {
136: for (i=0; i<yy_n; i++) {
137: yv[i] = PetscMax(yv[i],yvt[i]);
138: }
139: #endif
140: } else {SETERRQ(1,"Wrong insert option");}
141: }
142: }
143: VecRestoreArray(y,&yv);
144: if (x != y) {VecRestoreArray(x,&xv);}
145: return(0);
146: }
148: /*
149: This is special scatter code for when the entire parallel vector is
150: copied to processor 0.
152: */
153: int VecScatterBegin_MPI_ToOne(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
154: {
155: int rank,ierr,yy_n,xx_n,*range;
156: PetscScalar *xv,*yv;
157: MPI_Comm comm;
158: PetscMap map;
161: VecGetLocalSize(y,&yy_n);
162: VecGetLocalSize(x,&xx_n);
163: VecGetArray(x,&xv);
164: VecGetArray(y,&yv);
166: PetscObjectGetComm((PetscObject)x,&comm);
167: MPI_Comm_rank(comm,&rank);
169: /* -------- Reverse scatter; spread from processor 0 to other processors */
170: if (mode & SCATTER_REVERSE) {
171: PetscScalar *yvt;
172: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
173: int i;
175: VecGetPetscMap(y,&map);
176: PetscMapGetGlobalRange(map,&range);
177: if (addv == INSERT_VALUES) {
178: MPI_Scatterv(xv,scat->count,range,MPIU_SCALAR,yv,yy_n,MPIU_SCALAR,0,ctx->comm);
179: } else {
180: if (scat->work2) yvt = scat->work2;
181: else {
182: ierr = PetscMalloc((xx_n+1)*sizeof(PetscScalar),&yvt);
183: scat->work2 = yvt;
184: PetscLogObjectMemory(ctx,xx_n*sizeof(PetscScalar));
185: }
186: MPI_Scatterv(xv,scat->count,range,MPIU_SCALAR,yvt,yy_n,MPIU_SCALAR,0,ctx->comm);
187: if (addv == ADD_VALUES) {
188: for (i=0; i<yy_n; i++) {
189: yv[i] += yvt[i];
190: }
191: #if !defined(PETSC_USE_COMPLEX)
192: } else if (addv == MAX_VALUES) {
193: for (i=0; i<yy_n; i++) {
194: yv[i] = PetscMax(yv[i],yvt[i]);
195: }
196: #endif
197: } else {SETERRQ(1,"Wrong insert option");}
198: }
199: /* --------- Forward scatter; gather all values onto processor 0 */
200: } else {
201: PetscScalar *yvt = 0;
202: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
203: int i;
205: VecGetPetscMap(x,&map);
206: PetscMapGetGlobalRange(map,&range);
207: if (addv == INSERT_VALUES) {
208: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yv,scat->count,range,MPIU_SCALAR,0,ctx->comm);
209: } else {
210: if (!rank) {
211: if (scat->work1) yvt = scat->work1;
212: else {
213: ierr = PetscMalloc((yy_n+1)*sizeof(PetscScalar),&yvt);
214: scat->work1 = yvt;
215: PetscLogObjectMemory(ctx,yy_n*sizeof(PetscScalar));
216: }
217: }
218: MPI_Gatherv(xv,xx_n,MPIU_SCALAR,yvt,scat->count,range,MPIU_SCALAR,0,ctx->comm);
219: if (!rank) {
220: if (addv == ADD_VALUES) {
221: for (i=0; i<yy_n; i++) {
222: yv[i] += yvt[i];
223: }
224: #if !defined(PETSC_USE_COMPLEX)
225: } else if (addv == MAX_VALUES) {
226: for (i=0; i<yy_n; i++) {
227: yv[i] = PetscMax(yv[i],yvt[i]);
228: }
229: #endif
230: } else {SETERRQ(1,"Wrong insert option");}
231: }
232: }
233: }
234: VecRestoreArray(x,&xv);
235: VecRestoreArray(y,&yv);
236: return(0);
237: }
239: /*
240: The follow to are used for both VecScatterBegin_MPI_ToAll() and VecScatterBegin_MPI_ToOne()
241: */
242: int VecScatterDestroy_MPI_ToAll(VecScatter ctx)
243: {
244: VecScatter_MPI_ToAll *scat = (VecScatter_MPI_ToAll*)ctx->todata;
248: PetscFree(scat->count);
249: if (scat->work1) {PetscFree(scat->work1);}
250: if (scat->work2) {PetscFree(scat->work2);}
251: PetscFree(ctx->todata);
252: PetscLogObjectDestroy(ctx);
253: PetscHeaderDestroy(ctx);
254: return(0);
255: }
257: int VecScatterCopy_MPI_ToAll(VecScatter in,VecScatter out)
258: {
259: VecScatter_MPI_ToAll *in_to = (VecScatter_MPI_ToAll*)in->todata,*sto;
260: int size,i,ierr;
263: out->postrecvs = 0;
264: out->begin = in->begin;
265: out->end = in->end;
266: out->copy = in->copy;
267: out->destroy = in->destroy;
268: out->view = in->view;
270: ierr = PetscNew(VecScatter_MPI_ToAll,&sto);
271: sto->type = in_to->type;
273: MPI_Comm_size(out->comm,&size);
274: PetscMalloc(size*sizeof(int),&sto->count);
275: for (i=0; i<size; i++) {
276: sto->count[i] = in_to->count[i];
277: }
278: sto->work1 = 0;
279: sto->work2 = 0;
280: PetscLogObjectMemory(out,sizeof(VecScatter_MPI_ToAll)+size*sizeof(int));
281: out->todata = (void*)sto;
282: out->fromdata = (void*)0;
283: return(0);
284: }
286: /* --------------------------------------------------------------------------------------*/
287: /*
288: Scatter: sequential general to sequential general
289: */
290: int VecScatterBegin_SGtoSG(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
291: {
292: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
293: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
294: int i,n = gen_from->n,*fslots,*tslots,ierr;
295: PetscScalar *xv,*yv;
296:
298: VecGetArray(x,&xv);
299: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
300: if (mode & SCATTER_REVERSE){
301: gen_to = (VecScatter_Seq_General*)ctx->fromdata;
302: gen_from = (VecScatter_Seq_General*)ctx->todata;
303: }
304: fslots = gen_from->slots;
305: tslots = gen_to->slots;
307: if (addv == INSERT_VALUES) {
308: for (i=0; i<n; i++) {yv[tslots[i]] = xv[fslots[i]];}
309: } else if (addv == ADD_VALUES) {
310: for (i=0; i<n; i++) {yv[tslots[i]] += xv[fslots[i]];}
311: #if !defined(PETSC_USE_COMPLEX)
312: } else if (addv == MAX_VALUES) {
313: for (i=0; i<n; i++) {yv[tslots[i]] = PetscMax(yv[tslots[i]],xv[fslots[i]]);}
314: #endif
315: } else {SETERRQ(1,"Wrong insert option");}
316: VecRestoreArray(x,&xv);
317: if (x != y) {VecRestoreArray(y,&yv);}
318: return(0);
319: }
321: /*
322: Scatter: sequential general to sequential stride 1
323: */
324: int VecScatterBegin_SGtoSS_Stride1(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
325: {
326: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
327: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
328: int i,n = gen_from->n,*fslots = gen_from->slots;
329: int first = gen_to->first,ierr;
330: PetscScalar *xv,*yv;
331:
333: VecGetArray(x,&xv);
334: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
335: if (mode & SCATTER_REVERSE){
336: xv += first;
337: if (addv == INSERT_VALUES) {
338: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
339: } else if (addv == ADD_VALUES) {
340: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
341: #if !defined(PETSC_USE_COMPLEX)
342: } else if (addv == MAX_VALUES) {
343: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
344: #endif
345: } else {SETERRQ(1,"Wrong insert option");}
346: } else {
347: yv += first;
348: if (addv == INSERT_VALUES) {
349: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
350: } else if (addv == ADD_VALUES) {
351: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
352: #if !defined(PETSC_USE_COMPLEX)
353: } else if (addv == MAX_VALUES) {
354: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
355: #endif
356: } else {SETERRQ(1,"Wrong insert option");}
357: }
358: VecRestoreArray(x,&xv);
359: if (x != y) {VecRestoreArray(y,&yv);}
360: return(0);
361: }
363: /*
364: Scatter: sequential general to sequential stride
365: */
366: int VecScatterBegin_SGtoSS(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
367: {
368: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
369: VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
370: int i,n = gen_from->n,*fslots = gen_from->slots;
371: int first = gen_to->first,step = gen_to->step,ierr;
372: PetscScalar *xv,*yv;
373:
375: VecGetArray(x,&xv);
376: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
378: if (mode & SCATTER_REVERSE){
379: if (addv == INSERT_VALUES) {
380: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
381: } else if (addv == ADD_VALUES) {
382: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
383: #if !defined(PETSC_USE_COMPLEX)
384: } else if (addv == MAX_VALUES) {
385: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
386: #endif
387: } else {SETERRQ(1,"Wrong insert option");}
388: } else {
389: if (addv == INSERT_VALUES) {
390: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
391: } else if (addv == ADD_VALUES) {
392: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
393: #if !defined(PETSC_USE_COMPLEX)
394: } else if (addv == MAX_VALUES) {
395: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
396: #endif
397: } else {SETERRQ(1,"Wrong insert option");}
398: }
399: VecRestoreArray(x,&xv);
400: if (x != y) {VecRestoreArray(y,&yv);}
401: return(0);
402: }
404: /*
405: Scatter: sequential stride 1 to sequential general
406: */
407: int VecScatterBegin_SStoSG_Stride1(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
408: {
409: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
410: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
411: int i,n = gen_from->n,*fslots = gen_to->slots;
412: int first = gen_from->first,ierr;
413: PetscScalar *xv,*yv;
414:
416: VecGetArray(x,&xv);
417: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
419: if (mode & SCATTER_REVERSE){
420: yv += first;
421: if (addv == INSERT_VALUES) {
422: for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
423: } else if (addv == ADD_VALUES) {
424: for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
425: #if !defined(PETSC_USE_COMPLEX)
426: } else if (addv == MAX_VALUES) {
427: for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
428: #endif
429: } else {SETERRQ(1,"Wrong insert option");}
430: } else {
431: xv += first;
432: if (addv == INSERT_VALUES) {
433: for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
434: } else if (addv == ADD_VALUES) {
435: for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
436: #if !defined(PETSC_USE_COMPLEX)
437: } else if (addv == MAX_VALUES) {
438: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
439: #endif
440: } else {SETERRQ(1,"Wrong insert option");}
441: }
442: VecRestoreArray(x,&xv);
443: if (x != y) {VecRestoreArray(y,&yv);}
444: return(0);
445: }
447: /*
448: Scatter: sequential stride to sequential general
449: */
450: int VecScatterBegin_SStoSG(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
451: {
452: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
453: VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
454: int i,n = gen_from->n,*fslots = gen_to->slots;
455: int first = gen_from->first,step = gen_from->step,ierr;
456: PetscScalar *xv,*yv;
457:
459: VecGetArray(x,&xv);
460: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
462: if (mode & SCATTER_REVERSE){
463: if (addv == INSERT_VALUES) {
464: for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
465: } else if (addv == ADD_VALUES) {
466: for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
467: #if !defined(PETSC_USE_COMPLEX)
468: } else if (addv == MAX_VALUES) {
469: for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
470: #endif
471: } else {SETERRQ(1,"Wrong insert option");}
472: } else {
473: if (addv == INSERT_VALUES) {
474: for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
475: } else if (addv == ADD_VALUES) {
476: for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
477: #if !defined(PETSC_USE_COMPLEX)
478: } else if (addv == MAX_VALUES) {
479: for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
480: #endif
481: } else {SETERRQ(1,"Wrong insert option");}
482: }
483: VecRestoreArray(x,&xv);
484: if (x != y) {VecRestoreArray(y,&yv);}
485: return(0);
486: }
488: /*
489: Scatter: sequential stride to sequential stride
490: */
491: int VecScatterBegin_SStoSS(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
492: {
493: VecScatter_Seq_Stride *gen_to = (VecScatter_Seq_Stride*)ctx->todata;
494: VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
495: int i,n = gen_from->n,to_first = gen_to->first,to_step = gen_to->step;
496: int from_first = gen_from->first,from_step = gen_from->step,ierr;
497: PetscScalar *xv,*yv;
498:
500: VecGetArray(x,&xv);
501: if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
503: if (mode & SCATTER_REVERSE){
504: from_first = gen_to->first;
505: to_first = gen_from->first;
506: from_step = gen_to->step;
507: to_step = gen_from->step;
508: }
510: if (addv == INSERT_VALUES) {
511: if (to_step == 1 && from_step == 1) {
512: PetscMemcpy(yv+to_first,xv+from_first,n*sizeof(PetscScalar));
513: } else {
514: for (i=0; i<n; i++) {
515: yv[to_first + i*to_step] = xv[from_first+i*from_step];
516: }
517: }
518: } else if (addv == ADD_VALUES) {
519: if (to_step == 1 && from_step == 1) {
520: yv += to_first; xv += from_first;
521: for (i=0; i<n; i++) {
522: yv[i] += xv[i];
523: }
524: } else {
525: for (i=0; i<n; i++) {
526: yv[to_first + i*to_step] += xv[from_first+i*from_step];
527: }
528: }
529: #if !defined(PETSC_USE_COMPLEX)
530: } else if (addv == MAX_VALUES) {
531: if (to_step == 1 && from_step == 1) {
532: yv += to_first; xv += from_first;
533: for (i=0; i<n; i++) {
534: yv[i] = PetscMax(yv[i],xv[i]);
535: }
536: } else {
537: for (i=0; i<n; i++) {
538: yv[to_first + i*to_step] = PetscMax(yv[to_first + i*to_step],xv[from_first+i*from_step]);
539: }
540: }
541: #endif
542: } else {SETERRQ(1,"Wrong insert option");}
543: VecRestoreArray(x,&xv);
544: if (x != y) {VecRestoreArray(y,&yv);}
545: return(0);
546: }
548: /* --------------------------------------------------------------------------*/
551: int VecScatterCopy_SGToSG(VecScatter in,VecScatter out)
552: {
553: int ierr;
554: VecScatter_Seq_General *in_to = (VecScatter_Seq_General*)in->todata,*out_to;
555: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from;
556:
558: out->postrecvs = 0;
559: out->begin = in->begin;
560: out->end = in->end;
561: out->copy = in->copy;
562: out->destroy = in->destroy;
563: out->view = in->view;
565: ierr = PetscMalloc(in_to->n*sizeof(int)+sizeof(VecScatter_Seq_General),&out_to);
566: out_to->n = in_to->n;
567: out_to->type = in_to->type;
568: out_to->nonmatching_computed = PETSC_FALSE;
569: out_to->slots_nonmatching = 0;
570: out_to->is_copy = PETSC_FALSE;
571: out_to->slots = (int*)(out_to + 1);
572: PetscMemcpy(out_to->slots,in_to->slots,(out_to->n)*sizeof(int));
574: ierr = PetscMalloc(in_from->n*sizeof(int)+sizeof(VecScatter_Seq_General),&out_from);
575: out_from->n = in_from->n;
576: out_from->type = in_from->type;
577: out_from->nonmatching_computed = PETSC_FALSE;
578: out_from->slots_nonmatching = 0;
579: out_from->is_copy = PETSC_FALSE;
580: out_from->slots = (int*)(out_from + 1);
581: PetscMemcpy(out_from->slots,in_from->slots,(out_from->n)*sizeof(int));
583: PetscLogObjectMemory(out,2*sizeof(VecScatter_Seq_General)+(out_from->n+out_to->n)*sizeof(int));
584: out->todata = (void*)out_to;
585: out->fromdata = (void*)out_from;
586: return(0);
587: }
589: int VecScatterDestroy_SGtoSG(VecScatter ctx)
590: {
594: PetscFree(ctx->todata);
595: PetscFree(ctx->fromdata);
596: PetscLogObjectDestroy(ctx);
597: PetscHeaderDestroy(ctx);
598: return(0);
599: }
601: int VecScatterCopy_SGToStride(VecScatter in,VecScatter out)
602: {
604: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to;
605: VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from;
606:
608: out->postrecvs = 0;
609: out->begin = in->begin;
610: out->end = in->end;
611: out->copy = in->copy;
612: out->destroy = in->destroy;
613: out->view = in->view;
615: ierr = PetscNew(VecScatter_Seq_Stride,&out_to);
616: out_to->n = in_to->n;
617: out_to->type = in_to->type;
618: out_to->first = in_to->first;
619: out_to->step = in_to->step;
620: out_to->type = in_to->type;
622: ierr = PetscMalloc(in_from->n*sizeof(int)+sizeof(VecScatter_Seq_General),&out_from);
623: out_from->n = in_from->n;
624: out_from->type = in_from->type;
625: out_from->nonmatching_computed = PETSC_FALSE;
626: out_from->slots_nonmatching = 0;
627: out_from->is_copy = PETSC_FALSE;
628: out_from->slots = (int*)(out_from + 1);
629: PetscMemcpy(out_from->slots,in_from->slots,(out_from->n)*sizeof(int));
631: PetscLogObjectMemory(out,sizeof(VecScatter_Seq_General)+sizeof(VecScatter_Seq_Stride)+in_from->n*sizeof(int));
632: out->todata = (void*)out_to;
633: out->fromdata = (void*)out_from;
634: return(0);
635: }
637: /* --------------------------------------------------------------------------*/
638: /*
639: Scatter: parallel to sequential vector, sequential strides for both.
640: */
641: int VecScatterCopy_PStoSS(VecScatter in,VecScatter out)
642: {
643: VecScatter_Seq_Stride *in_to = (VecScatter_Seq_Stride*)in->todata,*out_to;
644: VecScatter_Seq_Stride *in_from = (VecScatter_Seq_Stride*)in->fromdata,*out_from;
645: int ierr;
648: out->postrecvs = 0;
649: out->begin = in->begin;
650: out->end = in->end;
651: out->copy = in->copy;
652: out->destroy = in->destroy;
653: out->view = in->view;
655: ierr = PetscNew(VecScatter_Seq_Stride,&out_to);
656: out_to->n = in_to->n;
657: out_to->type = in_to->type;
658: out_to->first = in_to->first;
659: out_to->step = in_to->step;
660: out_to->type = in_to->type;
661: ierr = PetscNew(VecScatter_Seq_Stride,&out_from);
662: PetscLogObjectMemory(out,2*sizeof(VecScatter_Seq_Stride));
663: out_from->n = in_from->n;
664: out_from->type = in_from->type;
665: out_from->first = in_from->first;
666: out_from->step = in_from->step;
667: out_from->type = in_from->type;
668: out->todata = (void*)out_to;
669: out->fromdata = (void*)out_from;
670: return(0);
671: }
673: EXTERN int VecScatterCreate_PtoS(int,int *,int,int *,Vec,Vec,int,VecScatter);
674: EXTERN int VecScatterCreate_PtoP(int,int *,int,int *,Vec,Vec,VecScatter);
675: EXTERN int VecScatterCreate_StoP(int,int *,int,int *,Vec,VecScatter);
677: /* =======================================================================*/
678: #define VEC_SEQ_ID 0
679: #define VEC_MPI_ID 1
681: /*@C
682: VecScatterCreate - Creates a vector scatter context.
684: Collective on Vec
686: Input Parameters:
687: + xin - a vector that defines the shape (parallel data layout of the vector)
688: of vectors from which we scatter
689: . yin - a vector that defines the shape (parallel data layout of the vector)
690: of vectors to which we scatter
691: . ix - the indices of xin to scatter
692: - iy - the indices of yin to hold results
694: Output Parameter:
695: . newctx - location to store the new scatter context
697: Options Database:
698: + -vecscatter_merge - Merges scatter send and receive (may offer better performance with some MPIs)
699: . -vecscatter_ssend - Uses MPI_Ssend_init() instead of MPI_Send_init() (may offer better performance with some MPIs)
700: . -vecscatter_sendfirst - Posts sends before receives (may offer better performance with some MPIs)
701: . -vecscatter_rr - use ready receiver mode for MPI sends in scatters (rarely used)
702: - -vecscatter_packtogether - Pack all messages before sending, receive all messages before unpacking
703: Level: intermediate
705: Notes:
706: In calls to VecScatter() you can use different vectors than the xin and
707: yin you used above; BUT they must have the same parallel data layout, for example,
708: they could be obtained from VecDuplicate().
709: A VecScatter context CANNOT be used in two or more simultaneous scatters;
710: that is you cannot call a second VecScatterBegin() with the same scatter
711: context until the VecScatterEnd() has been called on the first VecScatterBegin().
712: In this case a separate VecScatter is needed for each concurrent scatter.
714: Concepts: scatter^between vectors
715: Concepts: gather^between vectors
717: .seealso: VecScatterDestroy()
718: @*/
719: int VecScatterCreate(Vec xin,IS ix,Vec yin,IS iy,VecScatter *newctx)
720: {
721: VecScatter ctx;
722: int len,size,cando,totalv,ierr,*range,xin_type = VEC_SEQ_ID,yin_type = VEC_SEQ_ID;
723: PetscTruth flag;
724: MPI_Comm comm,ycomm;
725: PetscTruth ixblock,iyblock,iystride,islocal;
726: IS tix = 0,tiy = 0;
730: /*
731: Determine if the vectors are "parallel", ie. it shares a comm with other processors, or
732: sequential (it does not share a comm). The difference is that parallel vectors treat the
733: index set as providing indices in the global parallel numbering of the vector, with
734: sequential vectors treat the index set as providing indices in the local sequential
735: numbering
736: */
737: PetscObjectGetComm((PetscObject)xin,&comm);
738: MPI_Comm_size(comm,&size);
739: if (size > 1) {xin_type = VEC_MPI_ID;}
741: PetscObjectGetComm((PetscObject)yin,&ycomm);
742: MPI_Comm_size(ycomm,&size);
743: if (size > 1) {comm = ycomm; yin_type = VEC_MPI_ID;}
744:
745: /* generate the Scatter context */
746: PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",comm,VecScatterDestroy,VecScatterView);
747: PetscLogObjectCreate(ctx);
748: PetscLogObjectMemory(ctx,sizeof(struct _p_VecScatter));
749: ctx->inuse = PETSC_FALSE;
751: ctx->beginandendtogether = PETSC_FALSE;
752: PetscOptionsHasName(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether);
753: if (ctx->beginandendtogether) {
754: PetscLogInfo(ctx,"VecScatterCreate:Using combined (merged) vector scatter begin and endn");
755: }
756: PetscOptionsHasName(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether);
757: if (ctx->packtogether) {
758: PetscLogInfo(ctx,"VecScatterCreate:Pack all messages before sendingn");
759: }
761: VecGetLocalSize(xin,&ctx->from_n);
762: VecGetLocalSize(yin,&ctx->to_n);
764: /*
765: if ix or iy is not included; assume just grabbing entire vector
766: */
767: if (!ix && xin_type == VEC_SEQ_ID) {
768: ISCreateStride(comm,ctx->from_n,0,1,&ix);
769: tix = ix;
770: } else if (!ix && xin_type == VEC_MPI_ID) {
771: int bign;
772: VecGetSize(xin,&bign);
773: ISCreateStride(comm,bign,0,1,&ix);
774: tix = ix;
775: } else if (!ix) {
776: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
777: }
779: if (!iy && yin_type == VEC_SEQ_ID) {
780: ISCreateStride(comm,ctx->to_n,0,1,&iy);
781: tiy = iy;
782: } else if (!iy && yin_type == VEC_MPI_ID) {
783: int bign;
784: VecGetSize(yin,&bign);
785: ISCreateStride(comm,bign,0,1,&iy);
786: tiy = iy;
787: } else if (!iy) {
788: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
789: }
791: /* ===========================================================================================================
792: Check for special cases
793: ==========================================================================================================*/
794: /* ---------------------------------------------------------------------------*/
795: if (xin_type == VEC_SEQ_ID && yin_type == VEC_SEQ_ID) {
796: if (ix->type == IS_GENERAL && iy->type == IS_GENERAL){
797: int nx,ny,*idx,*idy;
798: VecScatter_Seq_General *to,*from;
800: ISGetLocalSize(ix,&nx);
801: ISGetLocalSize(iy,&ny);
802: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
803: ISGetIndices(ix,&idx);
804: ISGetIndices(iy,&idy);
805: len = sizeof(VecScatter_Seq_General) + nx*sizeof(int);
806: PetscMalloc(len,&to);
807: PetscLogObjectMemory(ctx,2*len);
808: to->slots = (int*)(to + 1);
809: to->n = nx;
810: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
811: PetscMemcpy(to->slots,idy,nx*sizeof(int));
812: PetscMalloc(len,&from);
813: from->slots = (int*)(from + 1);
814: from->n = nx;
815: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
816: PetscMemcpy(from->slots,idx,nx*sizeof(int));
817: to->type = VEC_SCATTER_SEQ_GENERAL;
818: from->type = VEC_SCATTER_SEQ_GENERAL;
819: ctx->todata = (void*)to;
820: ctx->fromdata = (void*)from;
821: ctx->postrecvs = 0;
822: ctx->begin = VecScatterBegin_SGtoSG;
823: ctx->end = 0;
824: ctx->destroy = VecScatterDestroy_SGtoSG;
825: ctx->copy = VecScatterCopy_SGToSG;
826: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential vector general scattern");
827: goto functionend;
828: } else if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
829: int nx,ny,to_first,to_step,from_first,from_step;
830: VecScatter_Seq_Stride *from8,*to8;
832: ISGetLocalSize(ix,&nx);
833: ISGetLocalSize(iy,&ny);
834: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
835: ISStrideGetInfo(iy,&to_first,&to_step);
836: ISStrideGetInfo(ix,&from_first,&from_step);
837: ierr = PetscNew(VecScatter_Seq_Stride,&to8);
838: to8->n = nx;
839: to8->first = to_first;
840: to8->step = to_step;
841: ierr = PetscNew(VecScatter_Seq_Stride,&from8);
842: PetscLogObjectMemory(ctx,2*sizeof(VecScatter_Seq_Stride));
843: from8->n = nx;
844: from8->first = from_first;
845: from8->step = from_step;
846: to8->type = VEC_SCATTER_SEQ_STRIDE;
847: from8->type = VEC_SCATTER_SEQ_STRIDE;
848: ctx->todata = (void*)to8;
849: ctx->fromdata = (void*)from8;
850: ctx->postrecvs = 0;
851: ctx->begin = VecScatterBegin_SStoSS;
852: ctx->end = 0;
853: ctx->destroy = VecScatterDestroy_SGtoSG;
854: ctx->copy = 0;
855: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential vector stride to striden");
856: goto functionend;
857: } else if (ix->type == IS_GENERAL && iy->type == IS_STRIDE){
858: int nx,ny,*idx,first,step;
859: VecScatter_Seq_General *from9;
860: VecScatter_Seq_Stride *to9;
862: ISGetLocalSize(ix,&nx);
863: ISGetIndices(ix,&idx);
864: ISGetLocalSize(iy,&ny);
865: ISStrideGetInfo(iy,&first,&step);
866: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
867: ierr = PetscNew(VecScatter_Seq_Stride,&to9);
868: to9->n = nx;
869: to9->first = first;
870: to9->step = step;
871: len = sizeof(VecScatter_Seq_General) + nx*sizeof(int);
872: PetscMalloc(len,&from9);
873: PetscLogObjectMemory(ctx,len + sizeof(VecScatter_Seq_Stride));
874: from9->slots = (int*)(from9 + 1);
875: from9->n = nx;
876: ierr = VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
877: ierr = PetscMemcpy(from9->slots,idx,nx*sizeof(int));
878: ctx->todata = (void*)to9; ctx->fromdata = (void*)from9;
879: ctx->postrecvs = 0;
880: if (step == 1) ctx->begin = VecScatterBegin_SGtoSS_Stride1;
881: else ctx->begin = VecScatterBegin_SGtoSS;
882: ctx->destroy = VecScatterDestroy_SGtoSG;
883: ctx->end = 0;
884: ctx->copy = VecScatterCopy_SGToStride;
885: to9->type = VEC_SCATTER_SEQ_STRIDE;
886: from9->type = VEC_SCATTER_SEQ_GENERAL;
887: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential vector general to striden");
888: goto functionend;
889: } else if (ix->type == IS_STRIDE && iy->type == IS_GENERAL){
890: int nx,ny,*idy,first,step;
891: VecScatter_Seq_General *to10;
892: VecScatter_Seq_Stride *from10;
894: ISGetLocalSize(ix,&nx);
895: ISGetIndices(iy,&idy);
896: ISGetLocalSize(iy,&ny);
897: ISStrideGetInfo(ix,&first,&step);
898: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
899: ierr = PetscNew(VecScatter_Seq_Stride,&from10);
900: from10->n = nx;
901: from10->first = first;
902: from10->step = step;
903: len = sizeof(VecScatter_Seq_General) + nx*sizeof(int);
904: PetscMalloc(len,&to10);
905: PetscLogObjectMemory(ctx,len + sizeof(VecScatter_Seq_Stride));
906: to10->slots = (int*)(to10 + 1);
907: to10->n = nx;
908: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
909: PetscMemcpy(to10->slots,idy,nx*sizeof(int));
910: ctx->todata = (void*)to10;
911: ctx->fromdata = (void*)from10;
912: ctx->postrecvs = 0;
913: if (step == 1) ctx->begin = VecScatterBegin_SStoSG_Stride1;
914: else ctx->begin = VecScatterBegin_SStoSG;
915: ctx->destroy = VecScatterDestroy_SGtoSG;
916: ctx->end = 0;
917: ctx->copy = 0;
918: to10->type = VEC_SCATTER_SEQ_GENERAL;
919: from10->type = VEC_SCATTER_SEQ_STRIDE;
920: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential vector stride to generaln");
921: goto functionend;
922: } else {
923: int nx,ny,*idx,*idy;
924: VecScatter_Seq_General *to11,*from11;
925: PetscTruth idnx,idny;
927: ISGetLocalSize(ix,&nx);
928: ISGetLocalSize(iy,&ny);
929: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
931: ISIdentity(ix,&idnx);
932: ISIdentity(iy,&idny);
933: if (idnx && idny) {
934: VecScatter_Seq_Stride *to13,*from13;
935: ierr = PetscNew(VecScatter_Seq_Stride,&to13);
936: to13->n = nx;
937: to13->first = 0;
938: to13->step = 1;
939: ierr = PetscNew(VecScatter_Seq_Stride,&from13);
940: PetscLogObjectMemory(ctx,2*sizeof(VecScatter_Seq_Stride));
941: from13->n = nx;
942: from13->first = 0;
943: from13->step = 1;
944: to13->type = VEC_SCATTER_SEQ_STRIDE;
945: from13->type = VEC_SCATTER_SEQ_STRIDE;
946: ctx->todata = (void*)to13;
947: ctx->fromdata = (void*)from13;
948: ctx->postrecvs = 0;
949: ctx->begin = VecScatterBegin_SStoSS;
950: ctx->end = 0;
951: ctx->destroy = VecScatterDestroy_SGtoSG;
952: ctx->copy = 0;
953: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential copyn");
954: goto functionend;
955: }
957: ISGetIndices(iy,&idy);
958: ISGetIndices(ix,&idx);
959: len = sizeof(VecScatter_Seq_General) + nx*sizeof(int);
960: PetscMalloc(len,&to11);
961: PetscLogObjectMemory(ctx,2*len);
962: to11->slots = (int*)(to11 + 1);
963: to11->n = nx;
964: VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
965: PetscMemcpy(to11->slots,idy,nx*sizeof(int));
966: PetscMalloc(len,&from11);
967: from11->slots = (int*)(from11 + 1);
968: from11->n = nx;
969: VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
970: PetscMemcpy(from11->slots,idx,nx*sizeof(int));
971: to11->type = VEC_SCATTER_SEQ_GENERAL;
972: from11->type = VEC_SCATTER_SEQ_GENERAL;
973: ctx->todata = (void*)to11;
974: ctx->fromdata = (void*)from11;
975: ctx->postrecvs = 0;
976: ctx->begin = VecScatterBegin_SGtoSG;
977: ctx->end = 0;
978: ctx->destroy = VecScatterDestroy_SGtoSG;
979: ctx->copy = VecScatterCopy_SGToSG;
980: ISRestoreIndices(ix,&idx);
981: ISRestoreIndices(iy,&idy);
982: PetscLogInfo(xin,"VecScatterCreate:Sequential vector scatter with block indicesn");
983: goto functionend;
984: }
985: }
986: /* ---------------------------------------------------------------------------*/
987: if (xin_type == VEC_MPI_ID && yin_type == VEC_SEQ_ID) {
989: /* ===========================================================================================================
990: Check for special cases
991: ==========================================================================================================*/
992: islocal = PETSC_FALSE;
993: /* special case extracting (subset of) local portion */
994: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
995: int nx,ny,to_first,to_step,from_first,from_step;
996: int start,end;
997: VecScatter_Seq_Stride *from12,*to12;
999: VecGetOwnershipRange(xin,&start,&end);
1000: ISGetLocalSize(ix,&nx);
1001: ISStrideGetInfo(ix,&from_first,&from_step);
1002: ISGetLocalSize(iy,&ny);
1003: ISStrideGetInfo(iy,&to_first,&to_step);
1004: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1005: if (ix->min >= start && ix->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1006: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1007: if (cando) {
1008: ierr = PetscNew(VecScatter_Seq_Stride,&to12);
1009: to12->n = nx;
1010: to12->first = to_first;
1011: to12->step = to_step;
1012: ierr = PetscNew(VecScatter_Seq_Stride,&from12);
1013: PetscLogObjectMemory(ctx,2*sizeof(VecScatter_Seq_Stride));
1014: from12->n = nx;
1015: from12->first = from_first-start;
1016: from12->step = from_step;
1017: to12->type = VEC_SCATTER_SEQ_STRIDE;
1018: from12->type = VEC_SCATTER_SEQ_STRIDE;
1019: ctx->todata = (void*)to12;
1020: ctx->fromdata = (void*)from12;
1021: ctx->postrecvs = 0;
1022: ctx->begin = VecScatterBegin_SStoSS;
1023: ctx->end = 0;
1024: ctx->destroy = VecScatterDestroy_SGtoSG;
1025: ctx->copy = VecScatterCopy_PStoSS;
1026: PetscLogInfo(xin,"VecScatterCreate:Special case: processors only getting local valuesn");
1027: goto functionend;
1028: }
1029: } else {
1030: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1031: }
1033: /* test for special case of all processors getting entire vector */
1034: totalv = 0;
1035: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1036: int i,nx,ny,to_first,to_step,from_first,from_step,*count,N;
1037: VecScatter_MPI_ToAll *sto;
1039: ISGetLocalSize(ix,&nx);
1040: ISStrideGetInfo(ix,&from_first,&from_step);
1041: ISGetLocalSize(iy,&ny);
1042: ISStrideGetInfo(iy,&to_first,&to_step);
1043: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1044: VecGetSize(xin,&N);
1045: if (nx != N) {
1046: totalv = 0;
1047: } else if (from_first == 0 && from_step == 1 &&
1048: from_first == to_first && from_step == to_step){
1049: totalv = 1;
1050: } else totalv = 0;
1051: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1053: if (cando) {
1054: PetscMap map;
1056: ierr = MPI_Comm_size(ctx->comm,&size);
1057: ierr = PetscNew(VecScatter_MPI_ToAll,&sto);
1058: ierr = PetscMalloc(size*sizeof(int),&count);
1059: ierr = VecGetPetscMap(xin,&map);
1060: ierr = PetscMapGetGlobalRange(map,&range);
1061: for (i=0; i<size; i++) {
1062: count[i] = range[i+1] - range[i];
1063: }
1064: sto->count = count;
1065: sto->work1 = 0;
1066: sto->work2 = 0;
1067: sto->type = VEC_SCATTER_MPI_TOALL;
1068: PetscLogObjectMemory(ctx,sizeof(VecScatter_MPI_ToAll)+size*sizeof(int));
1069: ctx->todata = (void*)sto;
1070: ctx->fromdata = 0;
1071: ctx->postrecvs = 0;
1072: ctx->begin = VecScatterBegin_MPI_ToAll;
1073: ctx->end = 0;
1074: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1075: ctx->copy = VecScatterCopy_MPI_ToAll;
1076: PetscLogInfo(xin,"VecScatterCreate:Special case: all processors get entire parallel vectorn");
1077: goto functionend;
1078: }
1079: } else {
1080: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1081: }
1083: /* test for special case of processor 0 getting entire vector */
1084: totalv = 0;
1085: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1086: int i,nx,ny,to_first,to_step,from_first,from_step,*count,rank,N;
1087: VecScatter_MPI_ToAll *sto;
1089: PetscObjectGetComm((PetscObject)xin,&comm);
1090: MPI_Comm_rank(comm,&rank);
1091: ISGetLocalSize(ix,&nx);
1092: ISStrideGetInfo(ix,&from_first,&from_step);
1093: ISGetLocalSize(iy,&ny);
1094: ISStrideGetInfo(iy,&to_first,&to_step);
1095: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1096: if (!rank) {
1097: VecGetSize(xin,&N);
1098: if (nx != N) {
1099: totalv = 0;
1100: } else if (from_first == 0 && from_step == 1 &&
1101: from_first == to_first && from_step == to_step){
1102: totalv = 1;
1103: } else totalv = 0;
1104: } else {
1105: if (!nx) totalv = 1;
1106: else totalv = 0;
1107: }
1108: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1110: if (cando) {
1111: PetscMap map;
1113: ierr = MPI_Comm_size(ctx->comm,&size);
1114: ierr = PetscNew(VecScatter_MPI_ToAll,&sto);
1115: ierr = PetscMalloc(size*sizeof(int),&count);
1116: ierr = VecGetPetscMap(xin,&map);
1117: ierr = PetscMapGetGlobalRange(map,&range);
1118: for (i=0; i<size; i++) {
1119: count[i] = range[i+1] - range[i];
1120: }
1121: sto->count = count;
1122: sto->work1 = 0;
1123: sto->work2 = 0;
1124: sto->type = VEC_SCATTER_MPI_TOONE;
1125: PetscLogObjectMemory(ctx,sizeof(VecScatter_MPI_ToAll)+size*sizeof(int));
1126: ctx->todata = (void*)sto;
1127: ctx->fromdata = 0;
1128: ctx->postrecvs = 0;
1129: ctx->begin = VecScatterBegin_MPI_ToOne;
1130: ctx->end = 0;
1131: ctx->destroy = VecScatterDestroy_MPI_ToAll;
1132: ctx->copy = VecScatterCopy_MPI_ToAll;
1133: PetscLogInfo(xin,"VecScatterCreate:Special case: processor zero gets entire parallel vector, rest get nonen");
1134: goto functionend;
1135: }
1136: } else {
1137: MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,xin->comm);
1138: }
1140: ISBlock(ix,&ixblock);
1141: ISBlock(iy,&iyblock);
1142: ISStride(iy,&iystride);
1143: if (ixblock) {
1144: /* special case block to block */
1145: if (iyblock) {
1146: int nx,ny,*idx,*idy,bsx,bsy;
1147: ISBlockGetBlockSize(iy,&bsy);
1148: ISBlockGetBlockSize(ix,&bsx);
1149: if (bsx == bsy && (bsx == 12 || bsx == 5 || bsx == 4 || bsx == 3 || bsx == 2)) {
1150: ISBlockGetSize(ix,&nx);
1151: ISBlockGetIndices(ix,&idx);
1152: ISBlockGetSize(iy,&ny);
1153: ISBlockGetIndices(iy,&idy);
1154: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1155: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,bsx,ctx);
1156: ISBlockRestoreIndices(ix,&idx);
1157: ISBlockRestoreIndices(iy,&idy);
1158: PetscLogInfo(xin,"VecScatterCreate:Special case: blocked indicesn");
1159: goto functionend;
1160: }
1161: /* special case block to stride */
1162: } else if (iystride) {
1163: int ystart,ystride,ysize,bsx;
1164: ISStrideGetInfo(iy,&ystart,&ystride);
1165: ISGetLocalSize(iy,&ysize);
1166: ISBlockGetBlockSize(ix,&bsx);
1167: /* see if stride index set is equivalent to block index set */
1168: if (((bsx == 2) || (bsx == 3) || (bsx == 4) || (bsx == 5) || (bsx == 12)) &&
1169: ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1170: int nx,*idx,*idy,il;
1171: ISBlockGetSize(ix,&nx); ISBlockGetIndices(ix,&idx);
1172: if (ysize != bsx*nx) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1173: PetscMalloc((nx+1)*sizeof(int),&idy);
1174: idy[0] = ystart;
1175: for (il=1; il<nx; il++) idy[il] = idy[il-1] + bsx;
1176: VecScatterCreate_PtoS(nx,idx,nx,idy,xin,yin,bsx,ctx);
1177: PetscFree(idy);
1178: ISBlockRestoreIndices(ix,&idx);
1179: PetscLogInfo(xin,"VecScatterCreate:Special case: blocked indices to striden");
1180: goto functionend;
1181: }
1182: }
1183: }
1184: /* left over general case */
1185: {
1186: int nx,ny,*idx,*idy;
1187: ISGetLocalSize(ix,&nx);
1188: ISGetIndices(ix,&idx);
1189: ISGetLocalSize(iy,&ny);
1190: ISGetIndices(iy,&idy);
1191: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1192: VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,1,ctx);
1193: ISRestoreIndices(ix,&idx);
1194: ISRestoreIndices(iy,&idy);
1195: goto functionend;
1196: }
1197: }
1198: /* ---------------------------------------------------------------------------*/
1199: if (xin_type == VEC_SEQ_ID && yin_type == VEC_MPI_ID) {
1200: /* ===========================================================================================================
1201: Check for special cases
1202: ==========================================================================================================*/
1203: /* special case local copy portion */
1204: islocal = PETSC_FALSE;
1205: if (ix->type == IS_STRIDE && iy->type == IS_STRIDE){
1206: int nx,ny,to_first,to_step,from_step,start,end,from_first;
1207: VecScatter_Seq_Stride *from,*to;
1209: VecGetOwnershipRange(yin,&start,&end);
1210: ISGetLocalSize(ix,&nx);
1211: ISStrideGetInfo(ix,&from_first,&from_step);
1212: ISGetLocalSize(iy,&ny);
1213: ISStrideGetInfo(iy,&to_first,&to_step);
1214: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1215: if (iy->min >= start && iy->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1216: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,yin->comm);
1217: if (cando) {
1218: ierr = PetscNew(VecScatter_Seq_Stride,&to);
1219: to->n = nx;
1220: to->first = to_first-start;
1221: to->step = to_step;
1222: ierr = PetscNew(VecScatter_Seq_Stride,&from);
1223: PetscLogObjectMemory(ctx,2*sizeof(VecScatter_Seq_Stride));
1224: from->n = nx;
1225: from->first = from_first;
1226: from->step = from_step;
1227: to->type = VEC_SCATTER_SEQ_STRIDE;
1228: from->type = VEC_SCATTER_SEQ_STRIDE;
1229: ctx->todata = (void*)to;
1230: ctx->fromdata = (void*)from;
1231: ctx->postrecvs = 0;
1232: ctx->begin = VecScatterBegin_SStoSS;
1233: ctx->end = 0;
1234: ctx->destroy = VecScatterDestroy_SGtoSG;
1235: ctx->copy = VecScatterCopy_PStoSS;
1236: PetscLogInfo(xin,"VecScatterCreate:Special case: sequential stride to striden");
1237: goto functionend;
1238: }
1239: } else {
1240: MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,yin->comm);
1241: }
1242: /* general case */
1243: {
1244: int nx,ny,*idx,*idy;
1245: ISGetLocalSize(ix,&nx);
1246: ISGetIndices(ix,&idx);
1247: ISGetLocalSize(iy,&ny);
1248: ISGetIndices(iy,&idy);
1249: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1250: VecScatterCreate_StoP(nx,idx,ny,idy,yin,ctx);
1251: ISRestoreIndices(ix,&idx);
1252: ISRestoreIndices(iy,&idy);
1253: goto functionend;
1254: }
1255: }
1256: /* ---------------------------------------------------------------------------*/
1257: if (xin_type == VEC_MPI_ID && yin_type == VEC_MPI_ID) {
1258: /* no special cases for now */
1259: int nx,ny,*idx,*idy;
1260: ierr = ISGetLocalSize(ix,&nx);
1261: ierr = ISGetIndices(ix,&idx);
1262: ierr = ISGetLocalSize(iy,&ny);
1263: ierr = ISGetIndices(iy,&idy);
1264: if (nx != ny) SETERRQ(PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1265: ierr = VecScatterCreate_PtoP(nx,idx,ny,idy,xin,yin,ctx);
1266: ierr = ISRestoreIndices(ix,&idx);
1267: ierr = ISRestoreIndices(iy,&idy);
1268: goto functionend;
1269: }
1271: functionend:
1272: *newctx = ctx;
1273: if (tix) {ISDestroy(tix);}
1274: if (tiy) {ISDestroy(tiy);}
1275: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view_info",&flag);
1276: if (flag) {
1277: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(comm),PETSC_VIEWER_ASCII_INFO);
1278: VecScatterView(ctx,PETSC_VIEWER_STDOUT_(comm));
1279: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(comm));
1280: }
1281: PetscOptionsHasName(PETSC_NULL,"-vecscatter_view",&flag);
1282: if (flag) {
1283: VecScatterView(ctx,PETSC_VIEWER_STDOUT_(comm));
1284: }
1285: return(0);
1286: }
1288: /* ------------------------------------------------------------------*/
1289: /*@
1290: VecScatterPostRecvs - Posts the receives required for the ready-receiver
1291: version of the VecScatter routines.
1293: Collective on VecScatter
1295: Input Parameters:
1296: + x - the vector from which we scatter (not needed,can be null)
1297: . y - the vector to which we scatter
1298: . addv - either ADD_VALUES or INSERT_VALUES
1299: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1300: SCATTER_FORWARD, SCATTER_REVERSE
1301: - inctx - scatter context generated by VecScatterCreate()
1303: Output Parameter:
1304: . y - the vector to which we scatter
1306: Level: advanced
1308: Notes:
1309: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1310: the SCATTER_FORWARD.
1311: The vectors x and y cannot be the same. y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1313: This scatter is far more general than the conventional
1314: scatter, since it can be a gather or a scatter or a combination,
1315: depending on the indices ix and iy. If x is a parallel vector and y
1316: is sequential, VecScatterBegin() can serve to gather values to a
1317: single processor. Similarly, if y is parallel and x sequential, the
1318: routine can scatter from one processor to many processors.
1320: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1321: @*/
1322: int VecScatterPostRecvs(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter inctx)
1323: {
1330: if (inctx->postrecvs) {
1331: (*inctx->postrecvs)(x,y,addv,mode,inctx);
1332: }
1333: return(0);
1334: }
1336: /* ------------------------------------------------------------------*/
1337: /*@
1338: VecScatterBegin - Begins a generalized scatter from one vector to
1339: another. Complete the scattering phase with VecScatterEnd().
1341: Collective on VecScatter and Vec
1343: Input Parameters:
1344: + x - the vector from which we scatter
1345: . y - the vector to which we scatter
1346: . addv - either ADD_VALUES or INSERT_VALUES
1347: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1348: SCATTER_FORWARD or SCATTER_REVERSE
1349: - inctx - scatter context generated by VecScatterCreate()
1351: Level: intermediate
1353: Notes:
1354: The vectors x and y need not be the same vectors used in the call
1355: to VecScatterCreate(), but x must have the same parallel data layout
1356: as that passed in as the x to VecScatterCreate(), similarly for the y.
1357: Most likely they have been obtained from VecDuplicate().
1359: You cannot change the values in the input vector between the calls to VecScatterBegin()
1360: and VecScatterEnd().
1362: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1363: the SCATTER_FORWARD.
1364:
1365: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1367: This scatter is far more general than the conventional
1368: scatter, since it can be a gather or a scatter or a combination,
1369: depending on the indices ix and iy. If x is a parallel vector and y
1370: is sequential, VecScatterBegin() can serve to gather values to a
1371: single processor. Similarly, if y is parallel and x sequential, the
1372: routine can scatter from one processor to many processors.
1374: Concepts: scatter^between vectors
1375: Concepts: gather^between vectors
1377: .seealso: VecScatterCreate(), VecScatterEnd()
1378: @*/
1379: int VecScatterBegin(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter inctx)
1380: {
1382: #if defined(PETSC_USE_BOPT_g)
1383: int to_n,from_n;
1384: #endif
1389: if (inctx->inuse) SETERRQ(PETSC_ERR_ARG_WRONGSTATE," Scatter ctx already in use");
1390: #if defined(PETSC_USE_BOPT_g)
1391: /*
1392: Error checking to make sure these vectors match the vectors used
1393: to create the vector scatter context. -1 in the from_n and to_n indicate the
1394: vector lengths are unknown (for example with mapped scatters) and thus
1395: no error checking is performed.
1396: */
1397: if (inctx->from_n >= 0 && inctx->to_n >= 0) {
1398: VecGetLocalSize(x,&from_n);
1399: VecGetLocalSize(y,&to_n);
1400: if (mode & SCATTER_REVERSE) {
1401: if (to_n != inctx->from_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter");
1402: if (from_n != inctx->to_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter");
1403: } else {
1404: if (to_n != inctx->to_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter");
1405: if (from_n != inctx->from_n) SETERRQ(PETSC_ERR_ARG_SIZ,"Vector wrong size for scatter");
1406: }
1407: }
1408: #endif
1410: inctx->inuse = PETSC_TRUE;
1411: PetscLogEventBarrierBegin(VEC_ScatterBarrier,0,0,0,0,inctx->comm);
1412: (*inctx->begin)(x,y,addv,mode,inctx);
1413: if (inctx->beginandendtogether && inctx->end) {
1414: inctx->inuse = PETSC_FALSE;
1415: (*inctx->end)(x,y,addv,mode,inctx);
1416: }
1417: PetscLogEventBarrierEnd(VEC_ScatterBarrier,0,0,0,0,inctx->comm);
1418: return(0);
1419: }
1421: /* --------------------------------------------------------------------*/
1422: /*@
1423: VecScatterEnd - Ends a generalized scatter from one vector to another. Call
1424: after first calling VecScatterBegin().
1426: Collective on VecScatter and Vec
1428: Input Parameters:
1429: + x - the vector from which we scatter
1430: . y - the vector to which we scatter
1431: . addv - either ADD_VALUES or INSERT_VALUES.
1432: . mode - the scattering mode, usually SCATTER_FORWARD. The available modes are:
1433: SCATTER_FORWARD, SCATTER_REVERSE
1434: - ctx - scatter context generated by VecScatterCreate()
1436: Level: intermediate
1438: Notes:
1439: If you use SCATTER_REVERSE the first two arguments should be reversed, from
1440: the SCATTER_FORWARD.
1441: y[iy[i]] = x[ix[i]], for i=0,...,ni-1
1443: .seealso: VecScatterBegin(), VecScatterCreate()
1444: @*/
1445: int VecScatterEnd(Vec x,Vec y,InsertMode addv,ScatterMode mode,VecScatter ctx)
1446: {
1452: ctx->inuse = PETSC_FALSE;
1453: if (!ctx->end) return(0);
1454: if (!ctx->beginandendtogether) {
1455: PetscLogEventBegin(VEC_ScatterEnd,ctx,x,y,0);
1456: (*(ctx)->end)(x,y,addv,mode,ctx);
1457: PetscLogEventEnd(VEC_ScatterEnd,ctx,x,y,0);
1458: }
1459: return(0);
1460: }
1462: /*@C
1463: VecScatterDestroy - Destroys a scatter context created by
1464: VecScatterCreate().
1466: Collective on VecScatter
1468: Input Parameter:
1469: . ctx - the scatter context
1471: Level: intermediate
1473: .seealso: VecScatterCreate(), VecScatterCopy()
1474: @*/
1475: int VecScatterDestroy(VecScatter ctx)
1476: {
1481: if (--ctx->refct > 0) return(0);
1483: /* if memory was published with AMS then destroy it */
1484: PetscObjectDepublish(ctx);
1486: (*ctx->destroy)(ctx);
1487: return(0);
1488: }
1490: /*@C
1491: VecScatterCopy - Makes a copy of a scatter context.
1493: Collective on VecScatter
1495: Input Parameter:
1496: . sctx - the scatter context
1498: Output Parameter:
1499: . ctx - the context copy
1501: Level: advanced
1503: .seealso: VecScatterCreate(), VecScatterDestroy()
1504: @*/
1505: int VecScatterCopy(VecScatter sctx,VecScatter *ctx)
1506: {
1512: if (!sctx->copy) SETERRQ(PETSC_ERR_SUP,"Cannot copy this type");
1513: PetscHeaderCreate(*ctx,_p_VecScatter,int,VEC_SCATTER_COOKIE,0,"VecScatter",sctx->comm,VecScatterDestroy,VecScatterView);
1514: PetscLogObjectCreate(*ctx);
1515: PetscLogObjectMemory(*ctx,sizeof(struct _p_VecScatter));
1516: (*ctx)->to_n = sctx->to_n;
1517: (*ctx)->from_n = sctx->from_n;
1518: (*sctx->copy)(sctx,*ctx);
1519: return(0);
1520: }
1523: /* ------------------------------------------------------------------*/
1524: /*@
1525: VecScatterView - Views a vector scatter context.
1527: Collective on VecScatter
1529: Input Parameters:
1530: + ctx - the scatter context
1531: - viewer - the viewer for displaying the context
1533: Level: intermediate
1535: @*/
1536: int VecScatterView(VecScatter ctx,PetscViewer viewer)
1537: {
1542: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ctx->comm);
1544: if (!ctx->view) SETERRQ(PETSC_ERR_SUP,"Cannot view this type of scatter context yet");
1546: (*ctx->view)(ctx,viewer);
1547: return(0);
1548: }
1550: /*@
1551: VecScatterRemap - Remaps the "from" and "to" indices in a
1552: vector scatter context. FOR EXPERTS ONLY!
1554: Collective on VecScatter
1556: Input Parameters:
1557: + scat - vector scatter context
1558: . from - remapping for "from" indices (may be PETSC_NULL)
1559: - to - remapping for "to" indices (may be PETSC_NULL)
1561: Level: developer
1563: Notes: In the parallel case the todata is actually the indices
1564: from which the data is TAKEN! The from stuff is where the
1565: data is finally put. This is VERY VERY confusing!
1567: In the sequential case the todata is the indices where the
1568: data is put and the fromdata is where it is taken from.
1569: This is backwards from the paralllel case! CRY! CRY! CRY!
1571: @*/
1572: int VecScatterRemap(VecScatter scat,int *rto,int *rfrom)
1573: {
1574: VecScatter_Seq_General *to,*from;
1575: VecScatter_MPI_General *mto;
1576: int i;
1583: from = (VecScatter_Seq_General *)scat->fromdata;
1584: mto = (VecScatter_MPI_General *)scat->todata;
1586: if (mto->type == VEC_SCATTER_MPI_TOALL) SETERRQ(PETSC_ERR_ARG_SIZ,"Not for to all scatter");
1588: if (rto) {
1589: if (mto->type == VEC_SCATTER_MPI_GENERAL) {
1590: /* handle off processor parts */
1591: for (i=0; i<mto->starts[mto->n]; i++) {
1592: mto->indices[i] = rto[mto->indices[i]];
1593: }
1594: /* handle local part */
1595: to = &mto->local;
1596: for (i=0; i<to->n; i++) {
1597: to->slots[i] = rto[to->slots[i]];
1598: }
1599: } else if (from->type == VEC_SCATTER_SEQ_GENERAL) {
1600: for (i=0; i<from->n; i++) {
1601: from->slots[i] = rto[from->slots[i]];
1602: }
1603: } else if (from->type == VEC_SCATTER_SEQ_STRIDE) {
1604: VecScatter_Seq_Stride *sto = (VecScatter_Seq_Stride*)from;
1605:
1606: /* if the remapping is the identity and stride is identity then skip remap */
1607: if (sto->step == 1 && sto->first == 0) {
1608: for (i=0; i<sto->n; i++) {
1609: if (rto[i] != i) {
1610: SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1611: }
1612: }
1613: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1614: } else SETERRQ(PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1615: }
1617: if (rfrom) {
1618: SETERRQ(PETSC_ERR_SUP,"Unable to remap the FROM in scatters yet");
1619: }
1621: /*
1622: Mark then vector lengths as unknown because we do not know the
1623: lengths of the remapped vectors
1624: */
1625: scat->from_n = -1;
1626: scat->to_n = -1;
1628: return(0);
1629: }