Actual source code: vscat.c

  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>              /*I "petscis.h" I*/
  9: #include <private/vecimpl.h>             /*I "petscvec.h" I*/

 11: /* Logging support */
 12: PetscClassId  VEC_SCATTER_CLASSID;

 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_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative index %D at %D location",idx[i],i);
 27:     if (idx[i] >= nmax) SETERRQ3(PETSC_COMM_SELF,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,MPIU_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_COMM_SELF,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,MPIU_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_COMM_SELF,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_COMM_SELF,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_COMM_SELF,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:   PetscFree3(ctx->todata,scat->count,scat->displx);
248:   return(0);
249: }

253: PetscErrorCode VecScatterDestroy_SGtoSG(VecScatter ctx)
254: {
255:   PetscErrorCode         ierr;

258:   PetscFree2(((VecScatter_Seq_General*)ctx->todata)->vslots,((VecScatter_Seq_General*)ctx->fromdata)->vslots);
259:   PetscFree2(ctx->todata,ctx->fromdata);
260:   return(0);
261: }

265: PetscErrorCode VecScatterDestroy_SGtoSS(VecScatter ctx)
266: {
267:   PetscErrorCode         ierr;

270:   PetscFree(((VecScatter_Seq_General*)ctx->fromdata)->vslots);
271:   PetscFree2(ctx->todata,ctx->fromdata);
272:   return(0);
273: }

277: PetscErrorCode VecScatterDestroy_SStoSG(VecScatter ctx)
278: {
279:   PetscErrorCode         ierr;

282:   PetscFree(((VecScatter_Seq_General*)ctx->todata)->vslots);
283:   PetscFree2(ctx->todata,ctx->fromdata);
284:   return(0);
285: }

289: PetscErrorCode VecScatterDestroy_SStoSS(VecScatter ctx)
290: {

294:   PetscFree2(ctx->todata,ctx->fromdata);
295:   return(0);
296: }

298: /* -------------------------------------------------------------------------*/
301: PetscErrorCode VecScatterCopy_MPI_ToAll(VecScatter in,VecScatter out)
302: {
303:   VecScatter_MPI_ToAll *in_to = (VecScatter_MPI_ToAll*)in->todata,*sto;
304:   PetscErrorCode       ierr;
305:   PetscMPIInt          size,*count,*displx;
306:   PetscInt             i;

309:   out->begin          = in->begin;
310:   out->end            = in->end;
311:   out->copy           = in->copy;
312:   out->destroy        = in->destroy;
313:   out->view           = in->view;

315:   MPI_Comm_size(((PetscObject)out)->comm,&size);
316:   PetscMalloc3(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count,size,PetscMPIInt,&displx);
317:   sto->type           = in_to->type;
318:   sto->count          = count;
319:   sto->displx         = displx;
320:   for (i=0; i<size; i++) {
321:     sto->count[i]  = in_to->count[i];
322:     sto->displx[i] = in_to->displx[i];
323:   }
324:   sto->work1         = 0;
325:   sto->work2         = 0;
326:   out->todata        = (void*)sto;
327:   out->fromdata      = (void*)0;
328:   return(0);
329: }

331: /* --------------------------------------------------------------------------------------*/
332: /* 
333:         Scatter: sequential general to sequential general 
334: */
337: PetscErrorCode VecScatterBegin_SGtoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
338: {
339:   VecScatter_Seq_General *gen_to = (VecScatter_Seq_General*)ctx->todata;
340:   VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
341:   PetscErrorCode         ierr;
342:   PetscInt               i,n = gen_from->n,*fslots,*tslots;
343:   PetscScalar            *xv,*yv;
344: 
346:   VecGetArray(x,&xv);
347:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
348:   if (mode & SCATTER_REVERSE){
349:     gen_to   = (VecScatter_Seq_General*)ctx->fromdata;
350:     gen_from = (VecScatter_Seq_General*)ctx->todata;
351:   }
352:   fslots = gen_from->vslots;
353:   tslots = gen_to->vslots;

355:   if (addv == INSERT_VALUES) {
356:     for (i=0; i<n; i++) {yv[tslots[i]] = xv[fslots[i]];}
357:   } else if (addv == ADD_VALUES) {
358:     for (i=0; i<n; i++) {yv[tslots[i]] += xv[fslots[i]];}
359: #if !defined(PETSC_USE_COMPLEX)
360:   } else  if (addv == MAX_VALUES) {
361:     for (i=0; i<n; i++) {yv[tslots[i]] = PetscMax(yv[tslots[i]],xv[fslots[i]]);}
362: #endif
363:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
364:   VecRestoreArray(x,&xv);
365:   if (x != y) {VecRestoreArray(y,&yv);}
366:   return(0);
367: }

369: /* 
370:     Scatter: sequential general to sequential stride 1 
371: */
374: PetscErrorCode VecScatterBegin_SGtoSS_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
375: {
376:   VecScatter_Seq_Stride  *gen_to   = (VecScatter_Seq_Stride*)ctx->todata;
377:   VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
378:   PetscInt               i,n = gen_from->n,*fslots = gen_from->vslots;
379:   PetscErrorCode         ierr;
380:   PetscInt               first = gen_to->first;
381:   PetscScalar            *xv,*yv;
382: 
384:   VecGetArray(x,&xv);
385:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}
386:   if (mode & SCATTER_REVERSE){
387:     xv += first;
388:     if (addv == INSERT_VALUES) {
389:       for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
390:     } else  if (addv == ADD_VALUES) {
391:       for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
392: #if !defined(PETSC_USE_COMPLEX)
393:     } else  if (addv == MAX_VALUES) {
394:       for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
395: #endif
396:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
397:   } else {
398:     yv += first;
399:     if (addv == INSERT_VALUES) {
400:       for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
401:     } else  if (addv == ADD_VALUES) {
402:       for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
403: #if !defined(PETSC_USE_COMPLEX)
404:     } else if (addv == MAX_VALUES) {
405:       for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
406: #endif
407:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
408:   }
409:   VecRestoreArray(x,&xv);
410:   if (x != y) {VecRestoreArray(y,&yv);}
411:   return(0);
412: }

414: /* 
415:    Scatter: sequential general to sequential stride 
416: */
419: PetscErrorCode VecScatterBegin_SGtoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
420: {
421:   VecScatter_Seq_Stride  *gen_to   = (VecScatter_Seq_Stride*)ctx->todata;
422:   VecScatter_Seq_General *gen_from = (VecScatter_Seq_General*)ctx->fromdata;
423:   PetscInt               i,n = gen_from->n,*fslots = gen_from->vslots;
424:   PetscErrorCode         ierr;
425:   PetscInt               first = gen_to->first,step = gen_to->step;
426:   PetscScalar            *xv,*yv;
427: 
429:   VecGetArray(x,&xv);
430:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}

432:   if (mode & SCATTER_REVERSE){
433:     if (addv == INSERT_VALUES) {
434:       for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
435:     } else if (addv == ADD_VALUES) {
436:       for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
437: #if !defined(PETSC_USE_COMPLEX)
438:     } else if (addv == MAX_VALUES) {
439:       for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
440: #endif
441:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
442:   } else {
443:     if (addv == INSERT_VALUES) {
444:       for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
445:     } else if (addv == ADD_VALUES) {
446:       for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
447: #if !defined(PETSC_USE_COMPLEX)
448:     } else if (addv == MAX_VALUES) {
449:       for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
450: #endif
451:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
452:   }
453:   VecRestoreArray(x,&xv);
454:   if (x != y) {VecRestoreArray(y,&yv);}
455:   return(0);
456: }

458: /* 
459:     Scatter: sequential stride 1 to sequential general 
460: */
463: PetscErrorCode VecScatterBegin_SStoSG_Stride1(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
464: {
465:   VecScatter_Seq_Stride  *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
466:   VecScatter_Seq_General *gen_to   = (VecScatter_Seq_General*)ctx->todata;
467:   PetscInt               i,n = gen_from->n,*fslots = gen_to->vslots;
468:   PetscErrorCode         ierr;
469:   PetscInt               first = gen_from->first;
470:   PetscScalar            *xv,*yv;
471: 
473:   VecGetArray(x,&xv);
474:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}

476:   if (mode & SCATTER_REVERSE){
477:     yv += first;
478:     if (addv == INSERT_VALUES) {
479:       for (i=0; i<n; i++) {yv[i] = xv[fslots[i]];}
480:     } else  if (addv == ADD_VALUES) {
481:       for (i=0; i<n; i++) {yv[i] += xv[fslots[i]];}
482: #if !defined(PETSC_USE_COMPLEX)
483:     } else  if (addv == MAX_VALUES) {
484:       for (i=0; i<n; i++) {yv[i] = PetscMax(yv[i],xv[fslots[i]]);}
485: #endif
486:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
487:   } else {
488:     xv += first;
489:     if (addv == INSERT_VALUES) {
490:       for (i=0; i<n; i++) {yv[fslots[i]] = xv[i];}
491:     } else  if (addv == ADD_VALUES) {
492:       for (i=0; i<n; i++) {yv[fslots[i]] += xv[i];}
493: #if !defined(PETSC_USE_COMPLEX)
494:     } else  if (addv == MAX_VALUES) {
495:       for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[i]);}
496: #endif
497:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
498:   }
499:   VecRestoreArray(x,&xv);
500:   if (x != y) {VecRestoreArray(y,&yv);}
501:   return(0);
502: }

506: /* 
507:    Scatter: sequential stride to sequential general 
508: */
509: PetscErrorCode VecScatterBegin_SStoSG(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
510: {
511:   VecScatter_Seq_Stride  *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
512:   VecScatter_Seq_General *gen_to   = (VecScatter_Seq_General*)ctx->todata;
513:   PetscInt               i,n = gen_from->n,*fslots = gen_to->vslots;
514:   PetscErrorCode         ierr;
515:   PetscInt               first = gen_from->first,step = gen_from->step;
516:   PetscScalar            *xv,*yv;
517: 
519:   VecGetArray(x,&xv);
520:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}

522:   if (mode & SCATTER_REVERSE){
523:     if (addv == INSERT_VALUES) {
524:       for (i=0; i<n; i++) {yv[first + i*step] = xv[fslots[i]];}
525:     } else  if (addv == ADD_VALUES) {
526:       for (i=0; i<n; i++) {yv[first + i*step] += xv[fslots[i]];}
527: #if !defined(PETSC_USE_COMPLEX)
528:     } else  if (addv == MAX_VALUES) {
529:       for (i=0; i<n; i++) {yv[first + i*step] = PetscMax(yv[first + i*step],xv[fslots[i]]);}
530: #endif
531:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
532:   } else {
533:     if (addv == INSERT_VALUES) {
534:       for (i=0; i<n; i++) {yv[fslots[i]] = xv[first + i*step];}
535:     } else  if (addv == ADD_VALUES) {
536:       for (i=0; i<n; i++) {yv[fslots[i]] += xv[first + i*step];}
537: #if !defined(PETSC_USE_COMPLEX)
538:     } else  if (addv == MAX_VALUES) {
539:       for (i=0; i<n; i++) {yv[fslots[i]] = PetscMax(yv[fslots[i]],xv[first + i*step]);}
540: #endif
541:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
542:   }
543:   VecRestoreArray(x,&xv);
544:   if (x != y) {VecRestoreArray(y,&yv);}
545:   return(0);
546: }

548: /* 
549:      Scatter: sequential stride to sequential stride 
550: */
553: PetscErrorCode VecScatterBegin_SStoSS(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
554: {
555:   VecScatter_Seq_Stride *gen_to   = (VecScatter_Seq_Stride*)ctx->todata;
556:   VecScatter_Seq_Stride *gen_from = (VecScatter_Seq_Stride*)ctx->fromdata;
557:   PetscInt              i,n = gen_from->n,to_first = gen_to->first,to_step = gen_to->step;
558:   PetscErrorCode        ierr;
559:   PetscInt              from_first = gen_from->first,from_step = gen_from->step;
560:   PetscScalar           *xv,*yv;
561: 
563:   VecGetArrayRead(x,(const PetscScalar **)&xv);
564:   if (x != y) {VecGetArray(y,&yv);} else {yv = xv;}

566:   if (mode & SCATTER_REVERSE){
567:     from_first = gen_to->first;
568:     to_first   = gen_from->first;
569:     from_step  = gen_to->step;
570:     to_step    = gen_from->step;
571:   }

573:   if (addv == INSERT_VALUES) {
574:     if (to_step == 1 && from_step == 1) {
575:       PetscMemcpy(yv+to_first,xv+from_first,n*sizeof(PetscScalar));
576:     } else  {
577:       for (i=0; i<n; i++) {
578:         yv[to_first + i*to_step] = xv[from_first+i*from_step];
579:       }
580:     }
581:   } else if (addv == ADD_VALUES) {
582:     if (to_step == 1 && from_step == 1) {
583:       yv += to_first; xv += from_first;
584:       for (i=0; i<n; i++) {
585:         yv[i] += xv[i];
586:       }
587:     } else {
588:       for (i=0; i<n; i++) {
589:         yv[to_first + i*to_step] += xv[from_first+i*from_step];
590:       }
591:     }
592: #if !defined(PETSC_USE_COMPLEX)
593:   } else if (addv == MAX_VALUES) {
594:     if (to_step == 1 && from_step == 1) {
595:       yv += to_first; xv += from_first;
596:       for (i=0; i<n; i++) {
597:         yv[i] = PetscMax(yv[i],xv[i]);
598:       }
599:     } else {
600:       for (i=0; i<n; i++) {
601:         yv[to_first + i*to_step] = PetscMax(yv[to_first + i*to_step],xv[from_first+i*from_step]);
602:       }
603:     }
604: #endif
605:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
606:   VecRestoreArrayRead(x,(const PetscScalar **)&xv);
607:   if (x != y) {VecRestoreArray(y,&yv);}
608:   return(0);
609: }

611: /* --------------------------------------------------------------------------*/


616: PetscErrorCode VecScatterCopy_SGToSG(VecScatter in,VecScatter out)
617: {
618:   PetscErrorCode         ierr;
619:   VecScatter_Seq_General *in_to   = (VecScatter_Seq_General*)in->todata,*out_to = PETSC_NULL;
620:   VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
621: 
623:   out->begin         = in->begin;
624:   out->end           = in->end;
625:   out->copy          = in->copy;
626:   out->destroy       = in->destroy;
627:   out->view          = in->view;

629:   PetscMalloc2(1,VecScatter_Seq_General,&out_to,1,VecScatter_Seq_General,&out_from);
630:   PetscMalloc2(in_to->n,PetscInt,&out_to->vslots,in_from->n,PetscInt,&out_from->vslots);
631:   out_to->n                      = in_to->n;
632:   out_to->type                   = in_to->type;
633:   out_to->nonmatching_computed   = PETSC_FALSE;
634:   out_to->slots_nonmatching      = 0;
635:   out_to->is_copy                = PETSC_FALSE;
636:   PetscMemcpy(out_to->vslots,in_to->vslots,(out_to->n)*sizeof(PetscInt));

638:   out_from->n                    = in_from->n;
639:   out_from->type                 = in_from->type;
640:   out_from->nonmatching_computed = PETSC_FALSE;
641:   out_from->slots_nonmatching    = 0;
642:   out_from->is_copy              = PETSC_FALSE;
643:   PetscMemcpy(out_from->vslots,in_from->vslots,(out_from->n)*sizeof(PetscInt));

645:   out->todata     = (void*)out_to;
646:   out->fromdata   = (void*)out_from;
647:   return(0);
648: }


653: PetscErrorCode VecScatterCopy_SGToSS(VecScatter in,VecScatter out)
654: {
655:   PetscErrorCode         ierr;
656:   VecScatter_Seq_Stride  *in_to   = (VecScatter_Seq_Stride*)in->todata,*out_to = PETSC_NULL;
657:   VecScatter_Seq_General *in_from = (VecScatter_Seq_General*)in->fromdata,*out_from = PETSC_NULL;
658: 
660:   out->begin         = in->begin;
661:   out->end           = in->end;
662:   out->copy          = in->copy;
663:   out->destroy       = in->destroy;
664:   out->view          = in->view;

666:   PetscMalloc2(1,VecScatter_Seq_Stride,&out_to,1,VecScatter_Seq_General,&out_from);
667:   PetscMalloc(in_from->n*sizeof(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: }


725: /* =======================================================================*/
726: #define VEC_SEQ_ID 0
727: #define VEC_MPI_ID 1
728: #define IS_GENERAL_ID 0
729: #define IS_STRIDE_ID  1
730: #define IS_BLOCK_ID   2

732: /*
733:    Blocksizes we have optimized scatters for 
734: */

736: #define VecScatterOptimizedBS(mbs) ((2 <= mbs && mbs <= 8) || mbs == 12)

738: PetscErrorCode  VecScatterCreateEmpty(MPI_Comm comm,VecScatter *newctx)
739: {
740:   VecScatter     ctx;

744:   PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_CLASSID,0,"VecScatter","VecScatter","Vec",comm,VecScatterDestroy,VecScatterView);
745:   ctx->inuse               = PETSC_FALSE;
746:   ctx->beginandendtogether = PETSC_FALSE;
747:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether,PETSC_NULL);
748:   if (ctx->beginandendtogether) {
749:     PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
750:   }
751:   ctx->packtogether = PETSC_FALSE;
752:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether,PETSC_NULL);
753:   if (ctx->packtogether) {
754:     PetscInfo(ctx,"Pack all messages before sending\n");
755:   }
756:   *newctx = ctx;
757:   return(0);
758: }

762: /*@C
763:    VecScatterCreate - Creates a vector scatter context.

765:    Collective on Vec

767:    Input Parameters:
768: +  xin - a vector that defines the shape (parallel data layout of the vector)
769:          of vectors from which we scatter
770: .  yin - a vector that defines the shape (parallel data layout of the vector)
771:          of vectors to which we scatter
772: .  ix - the indices of xin to scatter (if PETSC_NULL scatters all values)
773: -  iy - the indices of yin to hold results (if PETSC_NULL fills entire vector yin)

775:    Output Parameter:
776: .  newctx - location to store the new scatter context

778:    Options Database Keys: (uses regular MPI_Sends by default)
779: +  -vecscatter_view         - Prints detail of communications
780: .  -vecscatter_view_info    - Print less details about communication
781: .  -vecscatter_ssend        - Uses MPI_Ssend_init() instead of MPI_Send_init() 
782: .  -vecscatter_rsend           - use ready receiver mode for MPI sends 
783: .  -vecscatter_merge        - VecScatterBegin() handles all of the communication, VecScatterEnd() is a nop 
784:                               eliminates the chance for overlap of computation and communication 
785: .  -vecscatter_sendfirst    - Posts sends before receives 
786: .  -vecscatter_packtogether - Pack all messages before sending, receive all messages before unpacking
787: .  -vecscatter_alltoall     - Uses MPI all to all communication for scatter
788: .  -vecscatter_window       - Use MPI 2 window operations to move data
789: -  -vecscatter_nopack       - Avoid packing to work vector when possible (if used with -vecscatter_alltoall then will use MPI_Alltoallw()

791: $
792: $                                                                                    --When packing is used--
793: $                               MPI Datatypes (no packing)  sendfirst   merge        packtogether  persistent*    
794: $                                _nopack                   _sendfirst    _merge      _packtogether                -vecscatter_
795: $ ----------------------------------------------------------------------------------------------------------------------------
796: $    Message passing    Send       p                           X            X           X         always
797: $                      Ssend       p                           X            X           X         always          _ssend
798: $                      Rsend       p                        nonsense        X           X         always          _rsend
799: $    AlltoAll  v or w              X                        nonsense     always         X         nonsense        _alltoall
800: $    MPI_Win                       p                        nonsense        p           p         nonsense        _window
801: $                              
802: $   Since persistent sends and receives require a constant memory address they can only be used when data is packed into the work vector
803: $   because the in and out array may be different for each call to VecScatterBegin/End().
804: $
805: $    p indicates possible, but not implemented. X indicates implemented
806: $

808:     Level: intermediate

810:   Notes:
811:    In calls to VecScatter() you can use different vectors than the xin and 
812:    yin you used above; BUT they must have the same parallel data layout, for example,
813:    they could be obtained from VecDuplicate().
814:    A VecScatter context CANNOT be used in two or more simultaneous scatters;
815:    that is you cannot call a second VecScatterBegin() with the same scatter
816:    context until the VecScatterEnd() has been called on the first VecScatterBegin().
817:    In this case a separate VecScatter is needed for each concurrent scatter.

819:    Currently the MPI_Send(), MPI_Ssend() and MPI_Rsend() all use PERSISTENT versions.
820:    (this unfortunately requires that the same in and out arrays be used for each use, this
821:     is why when not using MPI_alltoallw() we always need to pack the input into the work array before sending
822:     and unpack upon receeving instead of using MPI datatypes to avoid the packing/unpacking).

824:    Concepts: scatter^between vectors
825:    Concepts: gather^between vectors

827: .seealso: VecScatterDestroy(), VecScatterCreateToAll(), VecScatterCreateToZero()
828: @*/
829: PetscErrorCode  VecScatterCreate(Vec xin,IS ix,Vec yin,IS iy,VecScatter *newctx)
830: {
831:   VecScatter     ctx;
833:   PetscMPIInt    size;
834:   PetscInt       totalv,xin_type = VEC_SEQ_ID,yin_type = VEC_SEQ_ID,*range;
835:   PetscInt       ix_type = IS_GENERAL_ID,iy_type = IS_GENERAL_ID;
836:   MPI_Comm       comm,ycomm;
837:   PetscBool      ixblock,iyblock,iystride,islocal,cando,flag;
838:   IS             tix = 0,tiy = 0;

841:   if (!ix && !iy) SETERRQ(((PetscObject)xin)->comm,PETSC_ERR_SUP,"Cannot pass default in for both input and output indices");

843:   /*
844:       Determine if the vectors are "parallel", ie. it shares a comm with other processors, or
845:       sequential (it does not share a comm). The difference is that parallel vectors treat the 
846:       index set as providing indices in the global parallel numbering of the vector, with 
847:       sequential vectors treat the index set as providing indices in the local sequential
848:       numbering
849:   */
850:   PetscObjectGetComm((PetscObject)xin,&comm);
851:   MPI_Comm_size(comm,&size);
852:   if (size > 1) {xin_type = VEC_MPI_ID;}

854:   PetscObjectGetComm((PetscObject)yin,&ycomm);
855:   MPI_Comm_size(ycomm,&size);
856:   if (size > 1) {comm = ycomm; yin_type = VEC_MPI_ID;}


859: 
860:   /* generate the Scatter context */
861:   PetscHeaderCreate(ctx,_p_VecScatter,int,VEC_SCATTER_CLASSID,0,"VecScatter","VecScatter","Vec",comm,VecScatterDestroy,VecScatterView);
862:   ctx->inuse               = PETSC_FALSE;

864:   ctx->beginandendtogether = PETSC_FALSE;
865:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_merge",&ctx->beginandendtogether,PETSC_NULL);
866:   if (ctx->beginandendtogether) {
867:     PetscInfo(ctx,"Using combined (merged) vector scatter begin and end\n");
868:   }
869:   ctx->packtogether = PETSC_FALSE;
870:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_packtogether",&ctx->packtogether,PETSC_NULL);
871:   if (ctx->packtogether) {
872:     PetscInfo(ctx,"Pack all messages before sending\n");
873:   }

875:   VecGetLocalSize(xin,&ctx->from_n);
876:   VecGetLocalSize(yin,&ctx->to_n);

878:   /*
879:       if ix or iy is not included; assume just grabbing entire vector
880:   */
881:   if (!ix && xin_type == VEC_SEQ_ID) {
882:     ISCreateStride(comm,ctx->from_n,0,1,&ix);
883:     tix  = ix;
884:   } else if (!ix && xin_type == VEC_MPI_ID) {
885:     if (yin_type == VEC_MPI_ID) {
886:       PetscInt ntmp, low;
887:       VecGetLocalSize(xin,&ntmp);
888:       VecGetOwnershipRange(xin,&low,PETSC_NULL);
889:       ISCreateStride(comm,ntmp,low,1,&ix);
890:     } else{
891:       PetscInt Ntmp;
892:       VecGetSize(xin,&Ntmp);
893:       ISCreateStride(comm,Ntmp,0,1,&ix);
894:     }
895:     tix  = ix;
896:   } else if (!ix) {
897:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ix not given, but not Seq or MPI vector");
898:   }

900:   if (!iy && yin_type == VEC_SEQ_ID) {
901:     ISCreateStride(comm,ctx->to_n,0,1,&iy);
902:     tiy  = iy;
903:   } else if (!iy && yin_type == VEC_MPI_ID) {
904:     if (xin_type == VEC_MPI_ID) {
905:       PetscInt ntmp, low;
906:       VecGetLocalSize(yin,&ntmp);
907:       VecGetOwnershipRange(yin,&low,PETSC_NULL);
908:       ISCreateStride(comm,ntmp,low,1,&iy);
909:     } else{
910:       PetscInt Ntmp;
911:       VecGetSize(yin,&Ntmp);
912:       ISCreateStride(comm,Ntmp,0,1,&iy);
913:     }
914:     tiy  = iy;
915:   } else if (!iy) {
916:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"iy not given, but not Seq or MPI vector");
917:   }

919:   /*
920:      Determine types of index sets
921:   */
922:   PetscTypeCompare((PetscObject)ix,ISBLOCK,&flag);
923:   if (flag) ix_type = IS_BLOCK_ID;
924:   PetscTypeCompare((PetscObject)iy,ISBLOCK,&flag);
925:   if (flag) iy_type = IS_BLOCK_ID;
926:   PetscTypeCompare((PetscObject)ix,ISSTRIDE,&flag);
927:   if (flag) ix_type = IS_STRIDE_ID;
928:   PetscTypeCompare((PetscObject)iy,ISSTRIDE,&flag);
929:   if (flag) iy_type = IS_STRIDE_ID;

931:   /* ===========================================================================================================
932:         Check for special cases
933:      ==========================================================================================================*/
934:   /* ---------------------------------------------------------------------------*/
935:   if (xin_type == VEC_SEQ_ID && yin_type == VEC_SEQ_ID) {
936:     if (ix_type == IS_GENERAL_ID && iy_type == IS_GENERAL_ID){
937:       PetscInt               nx,ny;
938:       const PetscInt         *idx,*idy;
939:       VecScatter_Seq_General *to = PETSC_NULL,*from = PETSC_NULL;

941:       ISGetLocalSize(ix,&nx);
942:       ISGetLocalSize(iy,&ny);
943:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
944:       ISGetIndices(ix,&idx);
945:       ISGetIndices(iy,&idy);
946:       PetscMalloc2(1,VecScatter_Seq_General,&to,1,VecScatter_Seq_General,&from);
947:       PetscMalloc2(nx,PetscInt,&to->vslots,nx,PetscInt,&from->vslots);
948:       to->n             = nx;
949: #if defined(PETSC_USE_DEBUG)
950:       VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
951: #endif
952:       PetscMemcpy(to->vslots,idy,nx*sizeof(PetscInt));
953:       from->n           = nx;
954: #if defined(PETSC_USE_DEBUG)
955:       VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
956: #endif
957:        PetscMemcpy(from->vslots,idx,nx*sizeof(PetscInt));
958:       to->type          = VEC_SCATTER_SEQ_GENERAL;
959:       from->type        = VEC_SCATTER_SEQ_GENERAL;
960:       ctx->todata       = (void*)to;
961:       ctx->fromdata     = (void*)from;
962:       ctx->begin        = VecScatterBegin_SGtoSG;
963:       ctx->end          = 0;
964:       ctx->destroy      = VecScatterDestroy_SGtoSG;
965:       ctx->copy         = VecScatterCopy_SGToSG;
966:       PetscInfo(xin,"Special case: sequential vector general scatter\n");
967:       goto functionend;
968:     } else if (ix_type == IS_STRIDE_ID &&  iy_type == IS_STRIDE_ID){
969:       PetscInt               nx,ny,to_first,to_step,from_first,from_step;
970:       VecScatter_Seq_Stride  *from8 = PETSC_NULL,*to8 = PETSC_NULL;

972:       ISGetLocalSize(ix,&nx);
973:       ISGetLocalSize(iy,&ny);
974:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
975:       ISStrideGetInfo(iy,&to_first,&to_step);
976:       ISStrideGetInfo(ix,&from_first,&from_step);
977:       PetscMalloc2(1,VecScatter_Seq_Stride,&to8,1,VecScatter_Seq_Stride,&from8);
978:       to8->n             = nx;
979:       to8->first         = to_first;
980:       to8->step          = to_step;
981:       from8->n           = nx;
982:       from8->first       = from_first;
983:       from8->step        = from_step;
984:       to8->type          = VEC_SCATTER_SEQ_STRIDE;
985:       from8->type        = VEC_SCATTER_SEQ_STRIDE;
986:       ctx->todata       = (void*)to8;
987:       ctx->fromdata     = (void*)from8;
988:       ctx->begin        = VecScatterBegin_SStoSS;
989:       ctx->end          = 0;
990:       ctx->destroy      = VecScatterDestroy_SStoSS;
991:       ctx->copy         = VecScatterCopy_SStoSS;
992:       PetscInfo(xin,"Special case: sequential vector stride to stride\n");
993:       goto functionend;
994:     } else if (ix_type == IS_GENERAL_ID && iy_type == IS_STRIDE_ID){
995:       PetscInt               nx,ny,first,step;
996:       const PetscInt         *idx;
997:       VecScatter_Seq_General *from9 = PETSC_NULL;
998:       VecScatter_Seq_Stride  *to9 = PETSC_NULL;

1000:       ISGetLocalSize(ix,&nx);
1001:       ISGetIndices(ix,&idx);
1002:       ISGetLocalSize(iy,&ny);
1003:       ISStrideGetInfo(iy,&first,&step);
1004:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1005:       PetscMalloc2(1,VecScatter_Seq_Stride,&to9,1,VecScatter_Seq_General,&from9);
1006:       PetscMalloc(nx*sizeof(PetscInt),&from9->vslots);
1007:       to9->n         = nx;
1008:       to9->first     = first;
1009:       to9->step      = step;
1010:       from9->n       = nx;
1011: #if defined(PETSC_USE_DEBUG)
1012:       VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
1013: #endif
1014:       PetscMemcpy(from9->vslots,idx,nx*sizeof(PetscInt));
1015:       ctx->todata    = (void*)to9; ctx->fromdata = (void*)from9;
1016:       if (step == 1)  ctx->begin = VecScatterBegin_SGtoSS_Stride1;
1017:       else            ctx->begin = VecScatterBegin_SGtoSS;
1018:       ctx->destroy   = VecScatterDestroy_SGtoSS;
1019:       ctx->end       = 0;
1020:       ctx->copy      = VecScatterCopy_SGToSS;
1021:       to9->type      = VEC_SCATTER_SEQ_STRIDE;
1022:       from9->type    = VEC_SCATTER_SEQ_GENERAL;
1023:       PetscInfo(xin,"Special case: sequential vector general to stride\n");
1024:       goto functionend;
1025:     } else if (ix_type == IS_STRIDE_ID && iy_type == IS_GENERAL_ID){
1026:       PetscInt               nx,ny,first,step;
1027:       const PetscInt         *idy;
1028:       VecScatter_Seq_General *to10 = PETSC_NULL;
1029:       VecScatter_Seq_Stride  *from10 = PETSC_NULL;

1031:       ISGetLocalSize(ix,&nx);
1032:       ISGetIndices(iy,&idy);
1033:       ISGetLocalSize(iy,&ny);
1034:       ISStrideGetInfo(ix,&first,&step);
1035:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1036:       PetscMalloc2(1,VecScatter_Seq_General,&to10,1,VecScatter_Seq_Stride,&from10);
1037:       PetscMalloc(nx*sizeof(PetscInt),&to10->vslots);
1038:       from10->n         = nx;
1039:       from10->first     = first;
1040:       from10->step      = step;
1041:       to10->n           = nx;
1042: #if defined(PETSC_USE_DEBUG)
1043:       VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1044: #endif
1045:       PetscMemcpy(to10->vslots,idy,nx*sizeof(PetscInt));
1046:       ctx->todata     = (void*)to10;
1047:       ctx->fromdata   = (void*)from10;
1048:       if (step == 1) ctx->begin = VecScatterBegin_SStoSG_Stride1;
1049:       else           ctx->begin = VecScatterBegin_SStoSG;
1050:       ctx->destroy    = VecScatterDestroy_SStoSG;
1051:       ctx->end        = 0;
1052:       ctx->copy       = 0;
1053:       to10->type      = VEC_SCATTER_SEQ_GENERAL;
1054:       from10->type    = VEC_SCATTER_SEQ_STRIDE;
1055:       PetscInfo(xin,"Special case: sequential vector stride to general\n");
1056:       goto functionend;
1057:     } else {
1058:       PetscInt               nx,ny;
1059:       const PetscInt         *idx,*idy;
1060:       VecScatter_Seq_General *to11 = PETSC_NULL,*from11 = PETSC_NULL;
1061:       PetscBool              idnx,idny;

1063:       ISGetLocalSize(ix,&nx);
1064:       ISGetLocalSize(iy,&ny);
1065:       if (nx != ny) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match, in %D out %D",nx,ny);

1067:       ISIdentity(ix,&idnx);
1068:       ISIdentity(iy,&idny);
1069:       if (idnx && idny) {
1070:         VecScatter_Seq_Stride *to13 = PETSC_NULL,*from13 = PETSC_NULL;
1071:         PetscMalloc2(1,VecScatter_Seq_Stride,&to13,1,VecScatter_Seq_Stride,&from13);
1072:         to13->n           = nx;
1073:         to13->first       = 0;
1074:         to13->step        = 1;
1075:         from13->n         = nx;
1076:         from13->first     = 0;
1077:         from13->step      = 1;
1078:         to13->type        = VEC_SCATTER_SEQ_STRIDE;
1079:         from13->type      = VEC_SCATTER_SEQ_STRIDE;
1080:         ctx->todata       = (void*)to13;
1081:         ctx->fromdata     = (void*)from13;
1082:         ctx->begin        = VecScatterBegin_SStoSS;
1083:         ctx->end          = 0;
1084:         ctx->destroy      = VecScatterDestroy_SStoSS;
1085:         ctx->copy         = VecScatterCopy_SStoSS;
1086:         PetscInfo(xin,"Special case: sequential copy\n");
1087:         goto functionend;
1088:       }

1090:       ISGetIndices(iy,&idy);
1091:       ISGetIndices(ix,&idx);
1092:       PetscMalloc2(1,VecScatter_Seq_General,&to11,1,VecScatter_Seq_General,&from11);
1093:       PetscMalloc2(nx,PetscInt,&to11->vslots,nx,PetscInt,&from11->vslots);
1094:       to11->n           = nx;
1095: #if defined(PETSC_USE_DEBUG)
1096:       VecScatterCheckIndices_Private(ctx->to_n,ny,idy);
1097: #endif
1098:       PetscMemcpy(to11->vslots,idy,nx*sizeof(PetscInt));
1099:       from11->n         = nx;
1100: #if defined(PETSC_USE_DEBUG)
1101:       VecScatterCheckIndices_Private(ctx->from_n,nx,idx);
1102: #endif
1103:       PetscMemcpy(from11->vslots,idx,nx*sizeof(PetscInt));
1104:       to11->type        = VEC_SCATTER_SEQ_GENERAL;
1105:       from11->type      = VEC_SCATTER_SEQ_GENERAL;
1106:       ctx->todata       = (void*)to11;
1107:       ctx->fromdata     = (void*)from11;
1108:       ctx->begin        = VecScatterBegin_SGtoSG;
1109:       ctx->end          = 0;
1110:       ctx->destroy      = VecScatterDestroy_SGtoSG;
1111:       ctx->copy         = VecScatterCopy_SGToSG;
1112:       ISRestoreIndices(ix,&idx);
1113:       ISRestoreIndices(iy,&idy);
1114:       PetscInfo(xin,"Sequential vector scatter with block indices\n");
1115:       goto functionend;
1116:     }
1117:   }
1118:   /* ---------------------------------------------------------------------------*/
1119:   if (xin_type == VEC_MPI_ID && yin_type == VEC_SEQ_ID) {

1121:   /* ===========================================================================================================
1122:         Check for special cases
1123:      ==========================================================================================================*/
1124:     islocal = PETSC_FALSE;
1125:     /* special case extracting (subset of) local portion */
1126:     if (ix_type == IS_STRIDE_ID && iy_type == IS_STRIDE_ID){
1127:       PetscInt              nx,ny,to_first,to_step,from_first,from_step;
1128:       PetscInt              start,end;
1129:       VecScatter_Seq_Stride *from12 = PETSC_NULL,*to12 = PETSC_NULL;

1131:       VecGetOwnershipRange(xin,&start,&end);
1132:       ISGetLocalSize(ix,&nx);
1133:       ISStrideGetInfo(ix,&from_first,&from_step);
1134:       ISGetLocalSize(iy,&ny);
1135:       ISStrideGetInfo(iy,&to_first,&to_step);
1136:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1137:       if (ix->min >= start && ix->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1138:       MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1139:       if (cando) {
1140:         PetscMalloc2(1,VecScatter_Seq_Stride,&to12,1,VecScatter_Seq_Stride,&from12);
1141:         to12->n             = nx;
1142:         to12->first         = to_first;
1143:         to12->step          = to_step;
1144:         from12->n           = nx;
1145:         from12->first       = from_first-start;
1146:         from12->step        = from_step;
1147:         to12->type          = VEC_SCATTER_SEQ_STRIDE;
1148:         from12->type        = VEC_SCATTER_SEQ_STRIDE;
1149:         ctx->todata         = (void*)to12;
1150:         ctx->fromdata       = (void*)from12;
1151:         ctx->begin          = VecScatterBegin_SStoSS;
1152:         ctx->end            = 0;
1153:         ctx->destroy        = VecScatterDestroy_SStoSS;
1154:         ctx->copy           = VecScatterCopy_SStoSS;
1155:         PetscInfo(xin,"Special case: processors only getting local values\n");
1156:         goto functionend;
1157:       }
1158:     } else {
1159:       MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1160:     }

1162:     /* test for special case of all processors getting entire vector */
1163:     /* contains check that PetscMPIInt can handle the sizes needed */
1164:     totalv = 0;
1165:     if (ix_type == IS_STRIDE_ID && iy_type == IS_STRIDE_ID){
1166:       PetscInt             i,nx,ny,to_first,to_step,from_first,from_step,N;
1167:       PetscMPIInt          *count = PETSC_NULL,*displx;
1168:       VecScatter_MPI_ToAll *sto = PETSC_NULL;

1170:       ISGetLocalSize(ix,&nx);
1171:       ISStrideGetInfo(ix,&from_first,&from_step);
1172:       ISGetLocalSize(iy,&ny);
1173:       ISStrideGetInfo(iy,&to_first,&to_step);
1174:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1175:       VecGetSize(xin,&N);
1176:       if (nx != N) {
1177:         totalv = 0;
1178:       } else if (from_first == 0 && from_step == 1 && from_first == to_first && from_step == to_step){
1179:         totalv = 1;
1180:       } else totalv = 0;
1181:       MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);

1183: #if defined(PETSC_USE_64BIT_INDICES)
1184:       if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1185: #else
1186:       if (cando) {
1187: #endif
1188:         MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1189:         PetscMalloc3(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count,size,PetscMPIInt,&displx);
1190:         range = xin->map->range;
1191:         for (i=0; i<size; i++) {
1192:           count[i]  = PetscMPIIntCast(range[i+1] - range[i]);
1193:           displx[i] = PetscMPIIntCast(range[i]);
1194:         }
1195:         sto->count        = count;
1196:         sto->displx       = displx;
1197:         sto->work1        = 0;
1198:         sto->work2        = 0;
1199:         sto->type         = VEC_SCATTER_MPI_TOALL;
1200:         ctx->todata       = (void*)sto;
1201:         ctx->fromdata     = 0;
1202:         ctx->begin        = VecScatterBegin_MPI_ToAll;
1203:         ctx->end          = 0;
1204:         ctx->destroy      = VecScatterDestroy_MPI_ToAll;
1205:         ctx->copy         = VecScatterCopy_MPI_ToAll;
1206:         PetscInfo(xin,"Special case: all processors get entire parallel vector\n");
1207:         goto functionend;
1208:       }
1209:     } else {
1210:       MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1211:     }

1213:     /* test for special case of processor 0 getting entire vector */
1214:     /* contains check that PetscMPIInt can handle the sizes needed */
1215:     totalv = 0;
1216:     if (ix_type == IS_STRIDE_ID && iy_type == IS_STRIDE_ID){
1217:       PetscInt             i,nx,ny,to_first,to_step,from_first,from_step,N;
1218:       PetscMPIInt          rank,*count = PETSC_NULL,*displx;
1219:       VecScatter_MPI_ToAll *sto = PETSC_NULL;

1221:       PetscObjectGetComm((PetscObject)xin,&comm);
1222:       MPI_Comm_rank(comm,&rank);
1223:       ISGetLocalSize(ix,&nx);
1224:       ISStrideGetInfo(ix,&from_first,&from_step);
1225:       ISGetLocalSize(iy,&ny);
1226:       ISStrideGetInfo(iy,&to_first,&to_step);
1227:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1228:       if (!rank) {
1229:         VecGetSize(xin,&N);
1230:         if (nx != N) {
1231:           totalv = 0;
1232:         } else if (from_first == 0        && from_step == 1 &&
1233:                    from_first == to_first && from_step == to_step){
1234:           totalv = 1;
1235:         } else totalv = 0;
1236:       } else {
1237:         if (!nx) totalv = 1;
1238:         else     totalv = 0;
1239:       }
1240:       MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);

1242: #if defined(PETSC_USE_64BIT_INDICES)
1243:       if (cando && (yin->map->N < PETSC_MPI_INT_MAX)) {
1244: #else
1245:       if (cando) {
1246: #endif
1247:         MPI_Comm_size(((PetscObject)ctx)->comm,&size);
1248:         PetscMalloc3(1,VecScatter_MPI_ToAll,&sto,size,PetscMPIInt,&count,size,PetscMPIInt,&displx);
1249:         range = xin->map->range;
1250:         for (i=0; i<size; i++) {
1251:           count[i] = PetscMPIIntCast(range[i+1] - range[i]);
1252:           displx[i] = PetscMPIIntCast(range[i]);
1253:         }
1254:         sto->count        = count;
1255:         sto->displx       = displx;
1256:         sto->work1        = 0;
1257:         sto->work2        = 0;
1258:         sto->type         = VEC_SCATTER_MPI_TOONE;
1259:         ctx->todata       = (void*)sto;
1260:         ctx->fromdata     = 0;
1261:         ctx->begin        = VecScatterBegin_MPI_ToOne;
1262:         ctx->end          = 0;
1263:         ctx->destroy      = VecScatterDestroy_MPI_ToAll;
1264:         ctx->copy         = VecScatterCopy_MPI_ToAll;
1265:         PetscInfo(xin,"Special case: processor zero gets entire parallel vector, rest get none\n");
1266:         goto functionend;
1267:       }
1268:     } else {
1269:       MPI_Allreduce(&totalv,&cando,1,MPI_INT,MPI_LAND,((PetscObject)xin)->comm);
1270:     }

1272:     PetscTypeCompare((PetscObject)ix,ISBLOCK,&ixblock);
1273:     PetscTypeCompare((PetscObject)iy,ISBLOCK,&iyblock);
1274:     PetscTypeCompare((PetscObject)iy,ISSTRIDE,&iystride);
1275:     if (ixblock) {
1276:       /* special case block to block */
1277:       if (iyblock) {
1278:         PetscInt       nx,ny,bsx,bsy;
1279:         const PetscInt *idx,*idy;
1280:         ISGetBlockSize(iy,&bsy);
1281:         ISGetBlockSize(ix,&bsx);
1282:         if (bsx == bsy && VecScatterOptimizedBS(bsx)) {
1283:           ISBlockGetLocalSize(ix,&nx);
1284:           ISBlockGetIndices(ix,&idx);
1285:           ISBlockGetLocalSize(iy,&ny);
1286:           ISBlockGetIndices(iy,&idy);
1287:           if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1288:           VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,bsx,ctx);
1289:           ISBlockRestoreIndices(ix,&idx);
1290:           ISBlockRestoreIndices(iy,&idy);
1291:           PetscInfo(xin,"Special case: blocked indices\n");
1292:           goto functionend;
1293:         }
1294:       /* special case block to stride */
1295:       } else if (iystride) {
1296:         PetscInt ystart,ystride,ysize,bsx;
1297:         ISStrideGetInfo(iy,&ystart,&ystride);
1298:         ISGetLocalSize(iy,&ysize);
1299:         ISGetBlockSize(ix,&bsx);
1300:         /* see if stride index set is equivalent to block index set */
1301:         if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1302:           PetscInt       nx,il,*idy;
1303:           const PetscInt *idx;
1304:           ISBlockGetLocalSize(ix,&nx);
1305:           ISBlockGetIndices(ix,&idx);
1306:           if (ysize != bsx*nx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1307:           PetscMalloc(nx*sizeof(PetscInt),&idy);
1308:           if (nx) {
1309:             idy[0] = ystart/bsx;
1310:             for (il=1; il<nx; il++) idy[il] = idy[il-1] + 1;
1311:           }
1312:           VecScatterCreate_PtoS(nx,idx,nx,idy,xin,yin,bsx,ctx);
1313:           PetscFree(idy);
1314:           ISBlockRestoreIndices(ix,&idx);
1315:           PetscInfo(xin,"Special case: blocked indices to stride\n");
1316:           goto functionend;
1317:         }
1318:       }
1319:     }
1320:     /* left over general case */
1321:     {
1322:       PetscInt       nx,ny;
1323:       const PetscInt *idx,*idy;
1324:       ISGetLocalSize(ix,&nx);
1325:       ISGetIndices(ix,&idx);
1326:       ISGetLocalSize(iy,&ny);
1327:       ISGetIndices(iy,&idy);
1328:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1329:       VecScatterCreate_PtoS(nx,idx,ny,idy,xin,yin,1,ctx);
1330:       ISRestoreIndices(ix,&idx);
1331:       ISRestoreIndices(iy,&idy);
1332:       PetscInfo(xin,"General case: MPI to Seq\n");
1333:       goto functionend;
1334:     }
1335:   }
1336:   /* ---------------------------------------------------------------------------*/
1337:   if (xin_type == VEC_SEQ_ID && yin_type == VEC_MPI_ID) {
1338:   /* ===========================================================================================================
1339:         Check for special cases
1340:      ==========================================================================================================*/
1341:     /* special case local copy portion */
1342:     islocal = PETSC_FALSE;
1343:     if (ix_type == IS_STRIDE_ID && iy_type == IS_STRIDE_ID){
1344:       PetscInt              nx,ny,to_first,to_step,from_step,start,end,from_first;
1345:       VecScatter_Seq_Stride *from = PETSC_NULL,*to = PETSC_NULL;

1347:       VecGetOwnershipRange(yin,&start,&end);
1348:       ISGetLocalSize(ix,&nx);
1349:       ISStrideGetInfo(ix,&from_first,&from_step);
1350:       ISGetLocalSize(iy,&ny);
1351:       ISStrideGetInfo(iy,&to_first,&to_step);
1352:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1353:       if (iy->min >= start && iy->max < end) islocal = PETSC_TRUE; else islocal = PETSC_FALSE;
1354:       MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1355:       if (cando) {
1356:         PetscMalloc2(1,VecScatter_Seq_Stride,&to,1,VecScatter_Seq_Stride,&from);
1357:         to->n             = nx;
1358:         to->first         = to_first-start;
1359:         to->step          = to_step;
1360:         from->n           = nx;
1361:         from->first       = from_first;
1362:         from->step        = from_step;
1363:         to->type          = VEC_SCATTER_SEQ_STRIDE;
1364:         from->type        = VEC_SCATTER_SEQ_STRIDE;
1365:         ctx->todata       = (void*)to;
1366:         ctx->fromdata     = (void*)from;
1367:         ctx->begin        = VecScatterBegin_SStoSS;
1368:         ctx->end          = 0;
1369:         ctx->destroy      = VecScatterDestroy_SStoSS;
1370:         ctx->copy         = VecScatterCopy_SStoSS;
1371:         PetscInfo(xin,"Special case: sequential stride to MPI stride\n");
1372:         goto functionend;
1373:       }
1374:     } else {
1375:       MPI_Allreduce(&islocal,&cando,1,MPI_INT,MPI_LAND,((PetscObject)yin)->comm);
1376:     }
1377:       /* special case block to stride */
1378:     if (ix_type == IS_BLOCK_ID && iy_type == IS_STRIDE_ID){
1379:       PetscInt ystart,ystride,ysize,bsx;
1380:       ISStrideGetInfo(iy,&ystart,&ystride);
1381:       ISGetLocalSize(iy,&ysize);
1382:       ISGetBlockSize(ix,&bsx);
1383:       /* see if stride index set is equivalent to block index set */
1384:       if (VecScatterOptimizedBS(bsx) && ((ystart % bsx) == 0) && (ystride == 1) && ((ysize % bsx) == 0)) {
1385:         PetscInt       nx,il,*idy;
1386:         const PetscInt *idx;
1387:         ISBlockGetLocalSize(ix,&nx);
1388:         ISBlockGetIndices(ix,&idx);
1389:         if (ysize != bsx*nx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1390:         PetscMalloc(nx*sizeof(PetscInt),&idy);
1391:         if (nx) {
1392:           idy[0] = ystart/bsx;
1393:           for (il=1; il<nx; il++) idy[il] = idy[il-1] + 1;
1394:         }
1395:         VecScatterCreate_StoP(nx,idx,nx,idy,xin,yin,bsx,ctx);
1396:         PetscFree(idy);
1397:         ISBlockRestoreIndices(ix,&idx);
1398:         PetscInfo(xin,"Special case: Blocked indices to stride\n");
1399:         goto functionend;
1400:       }
1401:     }

1403:     /* general case */
1404:     {
1405:       PetscInt       nx,ny;
1406:       const PetscInt *idx,*idy;
1407:       ISGetLocalSize(ix,&nx);
1408:       ISGetIndices(ix,&idx);
1409:       ISGetLocalSize(iy,&ny);
1410:       ISGetIndices(iy,&idy);
1411:       if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1412:       VecScatterCreate_StoP(nx,idx,ny,idy,xin,yin,1,ctx);
1413:       ISRestoreIndices(ix,&idx);
1414:       ISRestoreIndices(iy,&idy);
1415:       PetscInfo(xin,"General case: Seq to MPI\n");
1416:       goto functionend;
1417:     }
1418:   }
1419:   /* ---------------------------------------------------------------------------*/
1420:   if (xin_type == VEC_MPI_ID && yin_type == VEC_MPI_ID) {
1421:     /* no special cases for now */
1422:     PetscInt       nx,ny;
1423:     const PetscInt *idx,*idy;
1424:     ISGetLocalSize(ix,&nx);
1425:     ISGetIndices(ix,&idx);
1426:     ISGetLocalSize(iy,&ny);
1427:     ISGetIndices(iy,&idy);
1428:     if (nx != ny) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Local scatter sizes don't match");
1429:     VecScatterCreate_PtoP(nx,idx,ny,idy,xin,yin,1,ctx);
1430:     ISRestoreIndices(ix,&idx);
1431:     ISRestoreIndices(iy,&idy);
1432:     PetscInfo(xin,"General case: MPI to MPI\n");
1433:     goto functionend;
1434:   }

1436:   functionend:
1437:   *newctx = ctx;
1438:   ISDestroy(&tix);
1439:   ISDestroy(&tiy);
1440:   flag = PETSC_FALSE;
1441:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_view_info",&flag,PETSC_NULL);
1442:   if (flag) {
1443:     PetscViewer viewer;
1444:     PetscViewerASCIIGetStdout(comm,&viewer);
1445:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
1446:     VecScatterView(ctx,viewer);
1447:     PetscViewerPopFormat(viewer);
1448:   }
1449:   flag = PETSC_FALSE;
1450:   PetscOptionsGetBool(PETSC_NULL,"-vecscatter_view",&flag,PETSC_NULL);
1451:   if (flag) {
1452:     PetscViewer viewer;
1453:     PetscViewerASCIIGetStdout(comm,&viewer);
1454:     VecScatterView(ctx,viewer);
1455:   }
1456:   return(0);
1457: }

1459: /* ------------------------------------------------------------------*/
1462: /*@
1463:    VecScatterGetMerged - Returns true if the scatter is completed in the VecScatterBegin()
1464:       and the VecScatterEnd() does nothing

1466:    Not Collective

1468:    Input Parameter:
1469: .   ctx - scatter context created with VecScatterCreate()

1471:    Output Parameter:
1472: .   flg - PETSC_TRUE if the VecScatterBegin/End() are all done during the VecScatterBegin()

1474:    Level: developer

1476: .seealso: VecScatterCreate(), VecScatterEnd(), VecScatterBegin()
1477: @*/
1478: PetscErrorCode  VecScatterGetMerged(VecScatter ctx,PetscBool  *flg)
1479: {
1482:   *flg = ctx->beginandendtogether;
1483:   return(0);
1484: }

1488: /*@
1489:    VecScatterBegin - Begins a generalized scatter from one vector to
1490:    another. Complete the scattering phase with VecScatterEnd().

1492:    Neighbor-wise Collective on VecScatter and Vec

1494:    Input Parameters:
1495: +  inctx - scatter context generated by VecScatterCreate()
1496: .  x - the vector from which we scatter
1497: .  y - the vector to which we scatter
1498: .  addv - either ADD_VALUES or INSERT_VALUES, with INSERT_VALUES mode any location 
1499:           not scattered to retains its old value; i.e. the vector is NOT first zeroed.
1500: -  mode - the scattering mode, usually SCATTER_FORWARD.  The available modes are:
1501:     SCATTER_FORWARD or SCATTER_REVERSE


1504:    Level: intermediate

1506:    Options Database: See VecScatterCreate()

1508:    Notes:
1509:    The vectors x and y need not be the same vectors used in the call 
1510:    to VecScatterCreate(), but x must have the same parallel data layout
1511:    as that passed in as the x to VecScatterCreate(), similarly for the y.
1512:    Most likely they have been obtained from VecDuplicate().

1514:    You cannot change the values in the input vector between the calls to VecScatterBegin()
1515:    and VecScatterEnd().

1517:    If you use SCATTER_REVERSE the two arguments x and y should be reversed, from 
1518:    the SCATTER_FORWARD.
1519:    
1520:    y[iy[i]] = x[ix[i]], for i=0,...,ni-1

1522:    This scatter is far more general than the conventional
1523:    scatter, since it can be a gather or a scatter or a combination,
1524:    depending on the indices ix and iy.  If x is a parallel vector and y
1525:    is sequential, VecScatterBegin() can serve to gather values to a
1526:    single processor.  Similarly, if y is parallel and x sequential, the
1527:    routine can scatter from one processor to many processors.

1529:    Concepts: scatter^between vectors
1530:    Concepts: gather^between vectors

1532: .seealso: VecScatterCreate(), VecScatterEnd()
1533: @*/
1534: PetscErrorCode  VecScatterBegin(VecScatter inctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1535: {
1537: #if defined(PETSC_USE_DEBUG)
1538:   PetscInt      to_n,from_n;
1539: #endif

1545:   if (inctx->inuse) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE," Scatter ctx already in use");
1546:   CHKMEMQ;

1548: #if defined(PETSC_USE_DEBUG)
1549:   /*
1550:      Error checking to make sure these vectors match the vectors used
1551:    to create the vector scatter context. -1 in the from_n and to_n indicate the
1552:    vector lengths are unknown (for example with mapped scatters) and thus 
1553:    no error checking is performed.
1554:   */
1555:   if (inctx->from_n >= 0 && inctx->to_n >= 0) {
1556:     VecGetLocalSize(x,&from_n);
1557:     VecGetLocalSize(y,&to_n);
1558:     if (mode & SCATTER_REVERSE) {
1559:       if (to_n != inctx->from_n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter reverse and vector to != ctx from size)",to_n,inctx->from_n);
1560:       if (from_n != inctx->to_n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter reverse and vector from != ctx to size)",from_n,inctx->to_n);
1561:     } else {
1562:       if (to_n != inctx->to_n)     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter forward and vector to != ctx to size)",to_n,inctx->to_n);
1563:       if (from_n != inctx->from_n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Vector wrong size %D for scatter %D (scatter forward and vector from != ctx from size)",from_n,inctx->from_n);
1564:     }
1565:   }
1566: #endif

1568:   inctx->inuse = PETSC_TRUE;
1569:   PetscLogEventBarrierBegin(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1570:   (*inctx->begin)(inctx,x,y,addv,mode);
1571:   if (inctx->beginandendtogether && inctx->end) {
1572:     inctx->inuse = PETSC_FALSE;
1573:     (*inctx->end)(inctx,x,y,addv,mode);
1574:   }
1575:   PetscLogEventBarrierEnd(VEC_ScatterBarrier,0,0,0,0,((PetscObject)inctx)->comm);
1576:   CHKMEMQ;
1577:   return(0);
1578: }

1580: /* --------------------------------------------------------------------*/
1583: /*@
1584:    VecScatterEnd - Ends a generalized scatter from one vector to another.  Call
1585:    after first calling VecScatterBegin().

1587:    Neighbor-wise Collective on VecScatter and Vec

1589:    Input Parameters:
1590: +  ctx - scatter context generated by VecScatterCreate()
1591: .  x - the vector from which we scatter
1592: .  y - the vector to which we scatter
1593: .  addv - either ADD_VALUES or INSERT_VALUES.
1594: -  mode - the scattering mode, usually SCATTER_FORWARD.  The available modes are:
1595:      SCATTER_FORWARD, SCATTER_REVERSE

1597:    Level: intermediate

1599:    Notes:
1600:    If you use SCATTER_REVERSE the arguments x and y should be reversed, from the SCATTER_FORWARD.

1602:    y[iy[i]] = x[ix[i]], for i=0,...,ni-1

1604: .seealso: VecScatterBegin(), VecScatterCreate()
1605: @*/
1606: PetscErrorCode  VecScatterEnd(VecScatter ctx,Vec x,Vec y,InsertMode addv,ScatterMode mode)
1607: {

1614:   ctx->inuse = PETSC_FALSE;
1615:   if (!ctx->end) return(0);
1616:   CHKMEMQ;
1617:   if (!ctx->beginandendtogether) {
1618:     PetscLogEventBegin(VEC_ScatterEnd,ctx,x,y,0);
1619:     (*(ctx)->end)(ctx,x,y,addv,mode);
1620:     PetscLogEventEnd(VEC_ScatterEnd,ctx,x,y,0);
1621:   }
1622:   CHKMEMQ;
1623:   return(0);
1624: }

1628: /*@C
1629:    VecScatterDestroy - Destroys a scatter context created by 
1630:    VecScatterCreate().

1632:    Collective on VecScatter

1634:    Input Parameter:
1635: .  ctx - the scatter context

1637:    Level: intermediate

1639: .seealso: VecScatterCreate(), VecScatterCopy()
1640: @*/
1641: PetscErrorCode  VecScatterDestroy(VecScatter *ctx)
1642: {

1646:   if (!*ctx) return(0);
1648:   if ((*ctx)->inuse) SETERRQ(((PetscObject)(*ctx))->comm,PETSC_ERR_ARG_WRONGSTATE,"Scatter context is in use");
1649:   if (--((PetscObject)(*ctx))->refct > 0) {*ctx = 0; return(0);}

1651:   /* if memory was published with AMS then destroy it */
1652:   PetscObjectDepublish((*ctx));

1654:   (*(*ctx)->destroy)(*ctx);
1655:   PetscHeaderDestroy(ctx);
1656:   return(0);
1657: }

1661: /*@
1662:    VecScatterCopy - Makes a copy of a scatter context.

1664:    Collective on VecScatter

1666:    Input Parameter:
1667: .  sctx - the scatter context

1669:    Output Parameter:
1670: .  ctx - the context copy

1672:    Level: advanced

1674: .seealso: VecScatterCreate(), VecScatterDestroy()
1675: @*/
1676: PetscErrorCode  VecScatterCopy(VecScatter sctx,VecScatter *ctx)
1677: {

1683:   if (!sctx->copy) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot copy this type");
1684:   PetscHeaderCreate(*ctx,_p_VecScatter,int,VEC_SCATTER_CLASSID,0,"VecScatter","VecScatter","Vec",((PetscObject)sctx)->comm,VecScatterDestroy,VecScatterView);
1685:   (*ctx)->to_n   = sctx->to_n;
1686:   (*ctx)->from_n = sctx->from_n;
1687:   (*sctx->copy)(sctx,*ctx);
1688:   return(0);
1689: }


1692: /* ------------------------------------------------------------------*/
1695: /*@
1696:    VecScatterView - Views a vector scatter context.

1698:    Collective on VecScatter

1700:    Input Parameters:
1701: +  ctx - the scatter context
1702: -  viewer - the viewer for displaying the context

1704:    Level: intermediate

1706: @*/
1707: PetscErrorCode  VecScatterView(VecScatter ctx,PetscViewer viewer)
1708: {

1713:   if (!viewer) {
1714:     PetscViewerASCIIGetStdout(((PetscObject)ctx)->comm,&viewer);
1715:   }
1717:   if (ctx->view) {
1718:     (*ctx->view)(ctx,viewer);
1719:   }
1720:   return(0);
1721: }

1725: /*@C
1726:    VecScatterRemap - Remaps the "from" and "to" indices in a 
1727:    vector scatter context. FOR EXPERTS ONLY!

1729:    Collective on VecScatter

1731:    Input Parameters:
1732: +  scat - vector scatter context
1733: .  from - remapping for "from" indices (may be PETSC_NULL)
1734: -  to   - remapping for "to" indices (may be PETSC_NULL)

1736:    Level: developer

1738:    Notes: In the parallel case the todata is actually the indices
1739:           from which the data is TAKEN! The from stuff is where the 
1740:           data is finally put. This is VERY VERY confusing!

1742:           In the sequential case the todata is the indices where the 
1743:           data is put and the fromdata is where it is taken from.
1744:           This is backwards from the paralllel case! CRY! CRY! CRY!

1746: @*/
1747: PetscErrorCode  VecScatterRemap(VecScatter scat,PetscInt *rto,PetscInt *rfrom)
1748: {
1749:   VecScatter_Seq_General *to,*from;
1750:   VecScatter_MPI_General *mto;
1751:   PetscInt               i;


1758:   from = (VecScatter_Seq_General *)scat->fromdata;
1759:   mto  = (VecScatter_MPI_General *)scat->todata;

1761:   if (mto->type == VEC_SCATTER_MPI_TOALL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Not for to all scatter");

1763:   if (rto) {
1764:     if (mto->type == VEC_SCATTER_MPI_GENERAL) {
1765:       /* handle off processor parts */
1766:       for (i=0; i<mto->starts[mto->n]; i++) {
1767:         mto->indices[i] = rto[mto->indices[i]];
1768:       }
1769:       /* handle local part */
1770:       to = &mto->local;
1771:       for (i=0; i<to->n; i++) {
1772:         to->vslots[i] = rto[to->vslots[i]];
1773:       }
1774:     } else if (from->type == VEC_SCATTER_SEQ_GENERAL) {
1775:       for (i=0; i<from->n; i++) {
1776:         from->vslots[i] = rto[from->vslots[i]];
1777:       }
1778:     } else if (from->type == VEC_SCATTER_SEQ_STRIDE) {
1779:       VecScatter_Seq_Stride *sto = (VecScatter_Seq_Stride*)from;
1780: 
1781:       /* if the remapping is the identity and stride is identity then skip remap */
1782:       if (sto->step == 1 && sto->first == 0) {
1783:         for (i=0; i<sto->n; i++) {
1784:           if (rto[i] != i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1785:         }
1786:       } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1787:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Unable to remap such scatters");
1788:   }

1790:   if (rfrom) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unable to remap the FROM in scatters yet");

1792:   /*
1793:      Mark then vector lengths as unknown because we do not know the 
1794:    lengths of the remapped vectors
1795:   */
1796:   scat->from_n = -1;
1797:   scat->to_n   = -1;

1799:   return(0);
1800: }