Actual source code: pbvec.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:    This file contains routines for Parallel vector operations.
  4:  */
  5:  #include <petscsys.h>
  6:  #include <../src/vec/vec/impls/mpi/pvecimpl.h>

  8: PetscErrorCode VecDot_MPI(Vec xin,Vec yin,PetscScalar *z)
  9: {
 10:   PetscScalar    sum,work;

 14:   VecDot_Seq(xin,yin,&work);
 15:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 16:   *z   = sum;
 17:   return(0);
 18: }

 20: PetscErrorCode VecTDot_MPI(Vec xin,Vec yin,PetscScalar *z)
 21: {
 22:   PetscScalar    sum,work;

 26:   VecTDot_Seq(xin,yin,&work);
 27:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 28:   *z   = sum;
 29:   return(0);
 30: }

 32: extern PetscErrorCode VecView_MPI_Draw(Vec,PetscViewer);

 34: static PetscErrorCode VecPlaceArray_MPI(Vec vin,const PetscScalar *a)
 35: {
 37:   Vec_MPI        *v = (Vec_MPI*)vin->data;

 40:   if (v->unplacedarray) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"VecPlaceArray() was already called on this vector, without a call to VecResetArray()");
 41:   v->unplacedarray = v->array;  /* save previous array so reset can bring it back */
 42:   v->array         = (PetscScalar*)a;
 43:   if (v->localrep) {
 44:     VecPlaceArray(v->localrep,a);
 45:   }
 46:   return(0);
 47: }

 49: static PetscErrorCode VecDuplicate_MPI(Vec win,Vec *v)
 50: {
 52:   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
 53:   PetscScalar    *array;

 56:   VecCreate(PetscObjectComm((PetscObject)win),v);
 57:   PetscLayoutReference(win->map,&(*v)->map);

 59:   VecCreate_MPI_Private(*v,PETSC_TRUE,w->nghost,0);
 60:   vw   = (Vec_MPI*)(*v)->data;
 61:   PetscMemcpy((*v)->ops,win->ops,sizeof(struct _VecOps));

 63:   /* save local representation of the parallel vector (and scatter) if it exists */
 64:   if (w->localrep) {
 65:     VecGetArray(*v,&array);
 66:     VecCreateSeqWithArray(PETSC_COMM_SELF,PetscAbs(win->map->bs),win->map->n+w->nghost,array,&vw->localrep);
 67:     PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));
 68:     VecRestoreArray(*v,&array);
 69:     PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);

 71:     vw->localupdate = w->localupdate;
 72:     if (vw->localupdate) {
 73:       PetscObjectReference((PetscObject)vw->localupdate);
 74:     }
 75:   }

 77:   /* New vector should inherit stashing property of parent */
 78:   (*v)->stash.donotstash   = win->stash.donotstash;
 79:   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;

 81:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
 82:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);

 84:   (*v)->map->bs   = PetscAbs(win->map->bs);
 85:   (*v)->bstash.bs = win->bstash.bs;
 86:   return(0);
 87: }


 90: static PetscErrorCode VecSetOption_MPI(Vec V,VecOption op,PetscBool flag)
 91: {
 92:   Vec_MPI        *v = (Vec_MPI*)V->data;

 96:   switch (op) {
 97:   case VEC_IGNORE_OFF_PROC_ENTRIES: V->stash.donotstash = flag;
 98:     break;
 99:   case VEC_IGNORE_NEGATIVE_INDICES: V->stash.ignorenegidx = flag;
100:     break;
101:   case VEC_SUBSET_OFF_PROC_ENTRIES:
102:     v->assembly_subset = flag; /* See the same logic in MatAssembly wrt MAT_SUBSET_OFF_PROC_ENTRIES */
103:     if (!v->assembly_subset) { /* User indicates "do not reuse the communication pattern" */
104:       VecAssemblyReset_MPI(V); /* Reset existing pattern to free memory */
105:       v->first_assembly_done = PETSC_FALSE; /* Mark the first assembly is not done */
106:     }
107:     break;
108:   }

110:   return(0);
111: }


114: static PetscErrorCode VecResetArray_MPI(Vec vin)
115: {
116:   Vec_MPI        *v = (Vec_MPI*)vin->data;

120:   v->array         = v->unplacedarray;
121:   v->unplacedarray = 0;
122:   if (v->localrep) {
123:     VecResetArray(v->localrep);
124:   }
125:   return(0);
126: }

128: static PetscErrorCode VecAssemblySend_MPI_Private(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rankid,PetscMPIInt rank,void *sdata,MPI_Request req[],void *ctx)
129: {
130:   Vec X = (Vec)ctx;
131:   Vec_MPI *x = (Vec_MPI*)X->data;
132:   VecAssemblyHeader *hdr = (VecAssemblyHeader*)sdata;
133:   PetscInt bs = X->map->bs;

137:   /* x->first_assembly_done indicates we are reusing a communication network. In that case, some
138:      messages can be empty, but we have to send them this time if we sent them before because the
139:      receiver is expecting them.
140:    */
141:   if (hdr->count || (x->first_assembly_done && x->sendptrs[rankid].ints)) {
142:     MPI_Isend(x->sendptrs[rankid].ints,hdr->count,MPIU_INT,rank,tag[0],comm,&req[0]);
143:     MPI_Isend(x->sendptrs[rankid].scalars,hdr->count,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
144:   }
145:   if (hdr->bcount || (x->first_assembly_done && x->sendptrs[rankid].intb)) {
146:     MPI_Isend(x->sendptrs[rankid].intb,hdr->bcount,MPIU_INT,rank,tag[2],comm,&req[2]);
147:     MPI_Isend(x->sendptrs[rankid].scalarb,hdr->bcount*bs,MPIU_SCALAR,rank,tag[3],comm,&req[3]);
148:   }
149:   return(0);
150: }

152: static PetscErrorCode VecAssemblyRecv_MPI_Private(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *rdata,MPI_Request req[],void *ctx)
153: {
154:   Vec X = (Vec)ctx;
155:   Vec_MPI *x = (Vec_MPI*)X->data;
156:   VecAssemblyHeader *hdr = (VecAssemblyHeader*)rdata;
158:   PetscInt bs = X->map->bs;
159:   VecAssemblyFrame *frame;

162:   PetscSegBufferGet(x->segrecvframe,1,&frame);

164:   if (hdr->count) {
165:     PetscSegBufferGet(x->segrecvint,hdr->count,&frame->ints);
166:     MPI_Irecv(frame->ints,hdr->count,MPIU_INT,rank,tag[0],comm,&req[0]);
167:     PetscSegBufferGet(x->segrecvscalar,hdr->count,&frame->scalars);
168:     MPI_Irecv(frame->scalars,hdr->count,MPIU_SCALAR,rank,tag[1],comm,&req[1]);
169:     frame->pendings = 2;
170:   } else {
171:     frame->ints = NULL;
172:     frame->scalars = NULL;
173:     frame->pendings = 0;
174:   }

176:   if (hdr->bcount) {
177:     PetscSegBufferGet(x->segrecvint,hdr->bcount,&frame->intb);
178:     MPI_Irecv(frame->intb,hdr->bcount,MPIU_INT,rank,tag[2],comm,&req[2]);
179:     PetscSegBufferGet(x->segrecvscalar,hdr->bcount*bs,&frame->scalarb);
180:     MPI_Irecv(frame->scalarb,hdr->bcount*bs,MPIU_SCALAR,rank,tag[3],comm,&req[3]);
181:     frame->pendingb = 2;
182:   } else {
183:     frame->intb = NULL;
184:     frame->scalarb = NULL;
185:     frame->pendingb = 0;
186:   }
187:   return(0);
188: }

190: static PetscErrorCode VecAssemblyBegin_MPI_BTS(Vec X)
191: {
192:   Vec_MPI        *x = (Vec_MPI*)X->data;
194:   MPI_Comm       comm;
195:   PetscInt       i,j,jb,bs;

198:   if (X->stash.donotstash) return(0);

200:   PetscObjectGetComm((PetscObject)X,&comm);
201:   VecGetBlockSize(X,&bs);
202: #if defined(PETSC_USE_DEBUG)
203:   {
204:     InsertMode addv;
205:     MPIU_Allreduce((PetscEnum*)&X->stash.insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);
206:     if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(comm,PETSC_ERR_ARG_NOTSAMETYPE,"Some processors inserted values while others added");
207:   }
208: #endif
209:   X->bstash.insertmode = X->stash.insertmode; /* Block stash implicitly tracks InsertMode of scalar stash */

211:   VecStashSortCompress_Private(&X->stash);
212:   VecStashSortCompress_Private(&X->bstash);

214:   if (!x->sendranks) {
215:     PetscMPIInt nowners,bnowners,*owners,*bowners;
216:     PetscInt ntmp;
217:     VecStashGetOwnerList_Private(&X->stash,X->map,&nowners,&owners);
218:     VecStashGetOwnerList_Private(&X->bstash,X->map,&bnowners,&bowners);
219:     PetscMergeMPIIntArray(nowners,owners,bnowners,bowners,&ntmp,&x->sendranks);
220:     x->nsendranks = ntmp;
221:     PetscFree(owners);
222:     PetscFree(bowners);
223:     PetscMalloc1(x->nsendranks,&x->sendhdr);
224:     PetscCalloc1(x->nsendranks,&x->sendptrs);
225:   }
226:   for (i=0,j=0,jb=0; i<x->nsendranks; i++) {
227:     PetscMPIInt rank = x->sendranks[i];
228:     x->sendhdr[i].insertmode = X->stash.insertmode;
229:     /* Initialize pointers for non-empty stashes the first time around.  Subsequent assemblies with
230:      * VEC_SUBSET_OFF_PROC_ENTRIES will leave the old pointers (dangling because the stash has been collected) when
231:      * there is nothing new to send, so that size-zero messages get sent instead. */
232:     x->sendhdr[i].count = 0;
233:     if (X->stash.n) {
234:       x->sendptrs[i].ints    = &X->stash.idx[j];
235:       x->sendptrs[i].scalars = &X->stash.array[j];
236:       for ( ; j<X->stash.n && X->stash.idx[j] < X->map->range[rank+1]; j++) x->sendhdr[i].count++;
237:     }
238:     x->sendhdr[i].bcount = 0;
239:     if (X->bstash.n) {
240:       x->sendptrs[i].intb    = &X->bstash.idx[jb];
241:       x->sendptrs[i].scalarb = &X->bstash.array[jb*bs];
242:       for ( ; jb<X->bstash.n && X->bstash.idx[jb]*bs < X->map->range[rank+1]; jb++) x->sendhdr[i].bcount++;
243:     }
244:   }

246:   if (!x->segrecvint) {PetscSegBufferCreate(sizeof(PetscInt),1000,&x->segrecvint);}
247:   if (!x->segrecvscalar) {PetscSegBufferCreate(sizeof(PetscScalar),1000,&x->segrecvscalar);}
248:   if (!x->segrecvframe) {PetscSegBufferCreate(sizeof(VecAssemblyFrame),50,&x->segrecvframe);}
249:   if (x->first_assembly_done) { /* this is not the first assembly */
250:     PetscMPIInt tag[4];
251:     for (i=0; i<4; i++) {PetscCommGetNewTag(comm,&tag[i]);}
252:     for (i=0; i<x->nsendranks; i++) {
253:       VecAssemblySend_MPI_Private(comm,tag,i,x->sendranks[i],x->sendhdr+i,x->sendreqs+4*i,X);
254:     }
255:     for (i=0; i<x->nrecvranks; i++) {
256:       VecAssemblyRecv_MPI_Private(comm,tag,x->recvranks[i],x->recvhdr+i,x->recvreqs+4*i,X);
257:     }
258:     x->use_status = PETSC_TRUE;
259:   } else { /* First time assembly */
260:     PetscCommBuildTwoSidedFReq(comm,3,MPIU_INT,x->nsendranks,x->sendranks,(PetscInt*)x->sendhdr,&x->nrecvranks,&x->recvranks,&x->recvhdr,4,&x->sendreqs,&x->recvreqs,VecAssemblySend_MPI_Private,VecAssemblyRecv_MPI_Private,X);
261:     x->use_status = PETSC_FALSE;
262:   }

264:   /* The first_assembly_done flag is only meaningful when x->assembly_subset is set.
265:      This line says when assembly_subset is set, then we mark that the first assembly is done.
266:    */
267:   x->first_assembly_done = x->assembly_subset;

269:   {
270:     PetscInt nstash,reallocs;
271:     VecStashGetInfo_Private(&X->stash,&nstash,&reallocs);
272:     PetscInfo2(X,"Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
273:     VecStashGetInfo_Private(&X->bstash,&nstash,&reallocs);
274:     PetscInfo2(X,"Block-Stash has %D entries, uses %D mallocs.\n",nstash,reallocs);
275:   }
276:   return(0);
277: }

279: static PetscErrorCode VecAssemblyEnd_MPI_BTS(Vec X)
280: {
281:   Vec_MPI *x = (Vec_MPI*)X->data;
282:   PetscInt bs = X->map->bs;
283:   PetscMPIInt npending,*some_indices,r;
284:   MPI_Status  *some_statuses;
285:   PetscScalar *xarray;
287:   VecAssemblyFrame *frame;

290:   if (X->stash.donotstash) {
291:     X->stash.insertmode = NOT_SET_VALUES;
292:     X->bstash.insertmode = NOT_SET_VALUES;
293:     return(0);
294:   }

296:   if (!x->segrecvframe) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing segrecvframe! Probably you forgot to call VecAssemblyBegin first"); 
297:   VecGetArray(X,&xarray);
298:   PetscSegBufferExtractInPlace(x->segrecvframe,&frame);
299:   PetscMalloc2(4*x->nrecvranks,&some_indices,x->use_status?4*x->nrecvranks:0,&some_statuses);
300:   for (r=0,npending=0; r<x->nrecvranks; r++) npending += frame[r].pendings + frame[r].pendingb;
301:   while (npending>0) {
302:     PetscMPIInt ndone=0,ii;
303:     /* Filling MPI_Status fields requires some resources from the MPI library.  We skip it on the first assembly, or
304:      * when VEC_SUBSET_OFF_PROC_ENTRIES has not been set, because we could exchange exact sizes in the initial
305:      * rendezvous.  When the rendezvous is elided, however, we use MPI_Status to get actual message lengths, so that
306:      * subsequent assembly can set a proper subset of the values. */
307:     MPI_Waitsome(4*x->nrecvranks,x->recvreqs,&ndone,some_indices,x->use_status?some_statuses:MPI_STATUSES_IGNORE);
308:     for (ii=0; ii<ndone; ii++) {
309:       PetscInt i = some_indices[ii]/4,j,k;
310:       InsertMode imode = (InsertMode)x->recvhdr[i].insertmode;
311:       PetscInt *recvint;
312:       PetscScalar *recvscalar;
313:       PetscBool intmsg = (PetscBool)(some_indices[ii]%2 == 0);
314:       PetscBool blockmsg = (PetscBool)((some_indices[ii]%4)/2 == 1);
315:       npending--;
316:       if (!blockmsg) { /* Scalar stash */
317:         PetscMPIInt count;
318:         if (--frame[i].pendings > 0) continue;
319:         if (x->use_status) {
320:           MPI_Get_count(&some_statuses[ii],intmsg ? MPIU_INT : MPIU_SCALAR,&count);
321:         } else count = x->recvhdr[i].count;
322:         for (j=0,recvint=frame[i].ints,recvscalar=frame[i].scalars; j<count; j++,recvint++) {
323:           PetscInt loc = *recvint - X->map->rstart;
324:           if (*recvint < X->map->rstart || X->map->rend <= *recvint) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Received vector entry %D out of local range [%D,%D)]",*recvint,X->map->rstart,X->map->rend);
325:           switch (imode) {
326:           case ADD_VALUES:
327:             xarray[loc] += *recvscalar++;
328:             break;
329:           case INSERT_VALUES:
330:             xarray[loc] = *recvscalar++;
331:             break;
332:           default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Insert mode not supported 0x%x",imode);
333:           }
334:         }
335:       } else {                  /* Block stash */
336:         PetscMPIInt count;
337:         if (--frame[i].pendingb > 0) continue;
338:         if (x->use_status) {
339:           MPI_Get_count(&some_statuses[ii],intmsg ? MPIU_INT : MPIU_SCALAR,&count);
340:           if (!intmsg) count /= bs; /* Convert from number of scalars to number of blocks */
341:         } else count = x->recvhdr[i].bcount;
342:         for (j=0,recvint=frame[i].intb,recvscalar=frame[i].scalarb; j<count; j++,recvint++) {
343:           PetscInt loc = (*recvint)*bs - X->map->rstart;
344:           switch (imode) {
345:           case ADD_VALUES:
346:             for (k=loc; k<loc+bs; k++) xarray[k] += *recvscalar++;
347:             break;
348:           case INSERT_VALUES:
349:             for (k=loc; k<loc+bs; k++) xarray[k] = *recvscalar++;
350:             break;
351:           default: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Insert mode not supported 0x%x",imode);
352:           }
353:         }
354:       }
355:     }
356:   }
357:   VecRestoreArray(X,&xarray);
358:   MPI_Waitall(4*x->nsendranks,x->sendreqs,MPI_STATUSES_IGNORE);
359:   PetscFree2(some_indices,some_statuses);
360:   if (x->assembly_subset) {
361:     void *dummy;                /* reset segbuffers */
362:     PetscSegBufferExtractInPlace(x->segrecvint,&dummy);
363:     PetscSegBufferExtractInPlace(x->segrecvscalar,&dummy);
364:   } else {
365:     VecAssemblyReset_MPI(X);
366:   }

368:   X->stash.insertmode = NOT_SET_VALUES;
369:   X->bstash.insertmode = NOT_SET_VALUES;
370:   VecStashScatterEnd_Private(&X->stash);
371:   VecStashScatterEnd_Private(&X->bstash);
372:   return(0);
373: }

375: PetscErrorCode VecAssemblyReset_MPI(Vec X)
376: {
377:   Vec_MPI *x = (Vec_MPI*)X->data;

381:   PetscFree(x->sendreqs);
382:   PetscFree(x->recvreqs);
383:   PetscFree(x->sendranks);
384:   PetscFree(x->recvranks);
385:   PetscFree(x->sendhdr);
386:   PetscFree(x->recvhdr);
387:   PetscFree(x->sendptrs);
388:   PetscSegBufferDestroy(&x->segrecvint);
389:   PetscSegBufferDestroy(&x->segrecvscalar);
390:   PetscSegBufferDestroy(&x->segrecvframe);
391:   return(0);
392: }


395: static PetscErrorCode VecSetFromOptions_MPI(PetscOptionItems *PetscOptionsObject,Vec X)
396: {
397: #if !defined(PETSC_HAVE_MPIUNI)
399:   PetscBool      flg = PETSC_FALSE,set;

402:   PetscOptionsHead(PetscOptionsObject,"VecMPI Options");
403:   PetscOptionsBool("-vec_assembly_legacy","Use MPI 1 version of assembly","",flg,&flg,&set);
404:   if (set) {
405:     X->ops->assemblybegin = flg ? VecAssemblyBegin_MPI : VecAssemblyBegin_MPI_BTS;
406:     X->ops->assemblyend   = flg ? VecAssemblyEnd_MPI   : VecAssemblyEnd_MPI_BTS;
407:   }
408:   PetscOptionsTail();
409: #else
410:   X->ops->assemblybegin =  VecAssemblyBegin_MPI;
411:   X->ops->assemblyend   =  VecAssemblyEnd_MPI;
412: #endif
413:   return(0);
414: }


417: static struct _VecOps DvOps = { VecDuplicate_MPI, /* 1 */
418:                                 VecDuplicateVecs_Default,
419:                                 VecDestroyVecs_Default,
420:                                 VecDot_MPI,
421:                                 VecMDot_MPI,
422:                                 VecNorm_MPI,
423:                                 VecTDot_MPI,
424:                                 VecMTDot_MPI,
425:                                 VecScale_Seq,
426:                                 VecCopy_Seq, /* 10 */
427:                                 VecSet_Seq,
428:                                 VecSwap_Seq,
429:                                 VecAXPY_Seq,
430:                                 VecAXPBY_Seq,
431:                                 VecMAXPY_Seq,
432:                                 VecAYPX_Seq,
433:                                 VecWAXPY_Seq,
434:                                 VecAXPBYPCZ_Seq,
435:                                 VecPointwiseMult_Seq,
436:                                 VecPointwiseDivide_Seq,
437:                                 VecSetValues_MPI, /* 20 */
438:                                 VecAssemblyBegin_MPI_BTS,
439:                                 VecAssemblyEnd_MPI_BTS,
440:                                 0,
441:                                 VecGetSize_MPI,
442:                                 VecGetSize_Seq,
443:                                 0,
444:                                 VecMax_MPI,
445:                                 VecMin_MPI,
446:                                 VecSetRandom_Seq,
447:                                 VecSetOption_MPI,
448:                                 VecSetValuesBlocked_MPI,
449:                                 VecDestroy_MPI,
450:                                 VecView_MPI,
451:                                 VecPlaceArray_MPI,
452:                                 VecReplaceArray_Seq,
453:                                 VecDot_Seq,
454:                                 VecTDot_Seq,
455:                                 VecNorm_Seq,
456:                                 VecMDot_Seq,
457:                                 VecMTDot_Seq,
458:                                 VecLoad_Default,
459:                                 VecReciprocal_Default,
460:                                 VecConjugate_Seq,
461:                                 0,
462:                                 0,
463:                                 VecResetArray_MPI,
464:                                 VecSetFromOptions_MPI,/*set from options */
465:                                 VecMaxPointwiseDivide_Seq,
466:                                 VecPointwiseMax_Seq,
467:                                 VecPointwiseMaxAbs_Seq,
468:                                 VecPointwiseMin_Seq,
469:                                 VecGetValues_MPI,
470:                                 0,
471:                                 0,
472:                                 0,
473:                                 0,
474:                                 0,
475:                                 0,
476:                                 VecStrideGather_Default,
477:                                 VecStrideScatter_Default,
478:                                 0,
479:                                 0,
480:                                 0,
481:                                 0,
482:                                 0,
483:                                 VecStrideSubSetGather_Default,
484:                                 VecStrideSubSetScatter_Default,
485:                                 0,
486:                                 0
487: };

489: /*
490:     VecCreate_MPI_Private - Basic create routine called by VecCreate_MPI() (i.e. VecCreateMPI()),
491:     VecCreateMPIWithArray(), VecCreate_Shared() (i.e. VecCreateShared()), VecCreateGhost(),
492:     VecDuplicate_MPI(), VecCreateGhostWithArray(), VecDuplicate_MPI(), and VecDuplicate_Shared()

494:     If alloc is true and array is NULL then this routine allocates the space, otherwise
495:     no space is allocated.
496: */
497: PetscErrorCode VecCreate_MPI_Private(Vec v,PetscBool alloc,PetscInt nghost,const PetscScalar array[])
498: {
499:   Vec_MPI        *s;

503:   PetscNewLog(v,&s);
504:   v->data        = (void*)s;
505:   PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));
506:   s->nghost      = nghost;
507:   v->petscnative = PETSC_TRUE;

509:   PetscLayoutSetUp(v->map);

511:   s->array           = (PetscScalar*)array;
512:   s->array_allocated = 0;
513:   if (alloc && !array) {
514:     PetscInt n = v->map->n+nghost;
515:     PetscCalloc1(n,&s->array);
516:     PetscLogObjectMemory((PetscObject)v,n*sizeof(PetscScalar));
517:     s->array_allocated = s->array;
518:   }

520:   /* By default parallel vectors do not have local representation */
521:   s->localrep    = 0;
522:   s->localupdate = 0;

524:   v->stash.insertmode = NOT_SET_VALUES;
525:   v->bstash.insertmode = NOT_SET_VALUES;
526:   /* create the stashes. The block-size for bstash is set later when
527:      VecSetValuesBlocked is called.
528:   */
529:   VecStashCreate_Private(PetscObjectComm((PetscObject)v),1,&v->stash);
530:   VecStashCreate_Private(PetscObjectComm((PetscObject)v),PetscAbs(v->map->bs),&v->bstash);

532: #if defined(PETSC_HAVE_MATLAB_ENGINE)
533:   PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEnginePut_C",VecMatlabEnginePut_Default);
534:   PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEngineGet_C",VecMatlabEngineGet_Default);
535: #endif
536:   PetscObjectChangeTypeName((PetscObject)v,VECMPI);
537:   return(0);
538: }

540: /*MC
541:    VECMPI - VECMPI = "mpi" - The basic parallel vector

543:    Options Database Keys:
544: . -vec_type mpi - sets the vector type to VECMPI during a call to VecSetFromOptions()

546:   Level: beginner

548: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateMPI()
549: M*/

551: PETSC_EXTERN PetscErrorCode VecCreate_MPI(Vec vv)
552: {

556:   VecCreate_MPI_Private(vv,PETSC_TRUE,0,0);
557:   return(0);
558: }

560: /*MC
561:    VECSTANDARD = "standard" - A VECSEQ on one process and VECMPI on more than one process

563:    Options Database Keys:
564: . -vec_type standard - sets a vector type to standard on calls to VecSetFromOptions()

566:   Level: beginner

568: .seealso: VecCreateSeq(), VecCreateMPI()
569: M*/

571: PETSC_EXTERN PetscErrorCode VecCreate_Standard(Vec v)
572: {
574:   PetscMPIInt    size;

577:   MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
578:   if (size == 1) {
579:     VecSetType(v,VECSEQ);
580:   } else {
581:     VecSetType(v,VECMPI);
582:   }
583:   return(0);
584: }

586: /*@C
587:    VecCreateMPIWithArray - Creates a parallel, array-style vector,
588:    where the user provides the array space to store the vector values.

590:    Collective

592:    Input Parameters:
593: +  comm  - the MPI communicator to use
594: .  bs    - block size, same meaning as VecSetBlockSize()
595: .  n     - local vector length, cannot be PETSC_DECIDE
596: .  N     - global vector length (or PETSC_DECIDE to have calculated)
597: -  array - the user provided array to store the vector values

599:    Output Parameter:
600: .  vv - the vector

602:    Notes:
603:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
604:    same type as an existing vector.

606:    If the user-provided array is NULL, then VecPlaceArray() can be used
607:    at a later stage to SET the array for storing the vector values.

609:    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
610:    The user should not free the array until the vector is destroyed.

612:    Level: intermediate

614: .seealso: VecCreateSeqWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
615:           VecCreateMPI(), VecCreateGhostWithArray(), VecPlaceArray()

617: @*/
618: PetscErrorCode  VecCreateMPIWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar array[],Vec *vv)
619: {

623:   if (n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size of vector");
624:   PetscSplitOwnership(comm,&n,&N);
625:   VecCreate(comm,vv);
626:   VecSetSizes(*vv,n,N);
627:   VecSetBlockSize(*vv,bs);
628:   VecCreate_MPI_Private(*vv,PETSC_FALSE,0,array);
629:   return(0);
630: }

632: /*@C
633:    VecCreateGhostWithArray - Creates a parallel vector with ghost padding on each processor;
634:    the caller allocates the array space.

636:    Collective

638:    Input Parameters:
639: +  comm - the MPI communicator to use
640: .  n - local vector length
641: .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
642: .  nghost - number of local ghost points
643: .  ghosts - global indices of ghost points (or NULL if not needed), these do not need to be in increasing order (sorted)
644: -  array - the space to store the vector values (as long as n + nghost)

646:    Output Parameter:
647: .  vv - the global vector representation (without ghost points as part of vector)

649:    Notes:
650:    Use VecGhostGetLocalForm() to access the local, ghosted representation
651:    of the vector.

653:    This also automatically sets the ISLocalToGlobalMapping() for this vector.

655:    Level: advanced


658: .seealso: VecCreate(), VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
659:           VecCreateGhost(), VecCreateSeqWithArray(), VecCreateMPIWithArray(),
660:           VecCreateGhostBlock(), VecCreateGhostBlockWithArray(), VecMPISetGhost()

662: @*/
663: PetscErrorCode  VecCreateGhostWithArray(MPI_Comm comm,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],const PetscScalar array[],Vec *vv)
664: {
665:   PetscErrorCode         ierr;
666:   Vec_MPI                *w;
667:   PetscScalar            *larray;
668:   IS                     from,to;
669:   ISLocalToGlobalMapping ltog;
670:   PetscInt               rstart,i,*indices;

673:   *vv = 0;

675:   if (n == PETSC_DECIDE)      SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size");
676:   if (nghost == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local ghost size");
677:   if (nghost < 0)             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Ghost length must be >= 0");
678:   PetscSplitOwnership(comm,&n,&N);
679:   /* Create global representation */
680:   VecCreate(comm,vv);
681:   VecSetSizes(*vv,n,N);
682:   VecCreate_MPI_Private(*vv,PETSC_TRUE,nghost,array);
683:   w    = (Vec_MPI*)(*vv)->data;
684:   /* Create local representation */
685:   VecGetArray(*vv,&larray);
686:   VecCreateSeqWithArray(PETSC_COMM_SELF,1,n+nghost,larray,&w->localrep);
687:   PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localrep);
688:   VecRestoreArray(*vv,&larray);

690:   /*
691:        Create scatter context for scattering (updating) ghost values
692:   */
693:   ISCreateGeneral(comm,nghost,ghosts,PETSC_COPY_VALUES,&from);
694:   ISCreateStride(PETSC_COMM_SELF,nghost,n,1,&to);
695:   VecScatterCreate(*vv,from,w->localrep,to,&w->localupdate);
696:   PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localupdate);
697:   ISDestroy(&to);
698:   ISDestroy(&from);

700:   /* set local to global mapping for ghosted vector */
701:   PetscMalloc1(n+nghost,&indices);
702:   VecGetOwnershipRange(*vv,&rstart,NULL);
703:   for (i=0; i<n; i++) {
704:     indices[i] = rstart + i;
705:   }
706:   for (i=0; i<nghost; i++) {
707:     indices[n+i] = ghosts[i];
708:   }
709:   ISLocalToGlobalMappingCreate(comm,1,n+nghost,indices,PETSC_OWN_POINTER,&ltog);
710:   VecSetLocalToGlobalMapping(*vv,ltog);
711:   ISLocalToGlobalMappingDestroy(&ltog);
712:   return(0);
713: }

715: /*@
716:    VecCreateGhost - Creates a parallel vector with ghost padding on each processor.

718:    Collective

720:    Input Parameters:
721: +  comm - the MPI communicator to use
722: .  n - local vector length
723: .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
724: .  nghost - number of local ghost points
725: -  ghosts - global indices of ghost points, these do not need to be in increasing order (sorted)

727:    Output Parameter:
728: .  vv - the global vector representation (without ghost points as part of vector)

730:    Notes:
731:    Use VecGhostGetLocalForm() to access the local, ghosted representation
732:    of the vector.

734:    This also automatically sets the ISLocalToGlobalMapping() for this vector.

736:    Level: advanced

738: .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
739:           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(), VecGhostUpdateBegin(),
740:           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecGhostUpdateEnd(),
741:           VecCreateGhostBlock(), VecCreateGhostBlockWithArray(), VecMPISetGhost()

743: @*/
744: PetscErrorCode  VecCreateGhost(MPI_Comm comm,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],Vec *vv)
745: {

749:   VecCreateGhostWithArray(comm,n,N,nghost,ghosts,0,vv);
750:   return(0);
751: }

753: /*@
754:    VecMPISetGhost - Sets the ghost points for an MPI ghost vector

756:    Collective on Vec

758:    Input Parameters:
759: +  vv - the MPI vector
760: .  nghost - number of local ghost points
761: -  ghosts - global indices of ghost points, these do not need to be in increasing order (sorted)


764:    Notes:
765:    Use VecGhostGetLocalForm() to access the local, ghosted representation
766:    of the vector.

768:    This also automatically sets the ISLocalToGlobalMapping() for this vector.

770:    You must call this AFTER you have set the type of the vector (with VecSetType()) and the size (with VecSetSizes()).

772:    Level: advanced

774: .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
775:           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(), VecGhostUpdateBegin(),
776:           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecGhostUpdateEnd(),
777:           VecCreateGhostBlock(), VecCreateGhostBlockWithArray()

779: @*/
780: PetscErrorCode  VecMPISetGhost(Vec vv,PetscInt nghost,const PetscInt ghosts[])
781: {
783:   PetscBool      flg;

786:   PetscObjectTypeCompare((PetscObject)vv,VECMPI,&flg);
787:   /* if already fully existant VECMPI then basically destroy it and rebuild with ghosting */
788:   if (flg) {
789:     PetscInt               n,N;
790:     Vec_MPI                *w;
791:     PetscScalar            *larray;
792:     IS                     from,to;
793:     ISLocalToGlobalMapping ltog;
794:     PetscInt               rstart,i,*indices;
795:     MPI_Comm               comm;

797:     PetscObjectGetComm((PetscObject)vv,&comm);
798:     n    = vv->map->n;
799:     N    = vv->map->N;
800:     (*vv->ops->destroy)(vv);
801:     VecSetSizes(vv,n,N);
802:     VecCreate_MPI_Private(vv,PETSC_TRUE,nghost,NULL);
803:     w    = (Vec_MPI*)(vv)->data;
804:     /* Create local representation */
805:     VecGetArray(vv,&larray);
806:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,n+nghost,larray,&w->localrep);
807:     PetscLogObjectParent((PetscObject)vv,(PetscObject)w->localrep);
808:     VecRestoreArray(vv,&larray);

810:     /*
811:      Create scatter context for scattering (updating) ghost values
812:      */
813:     ISCreateGeneral(comm,nghost,ghosts,PETSC_COPY_VALUES,&from);
814:     ISCreateStride(PETSC_COMM_SELF,nghost,n,1,&to);
815:     VecScatterCreate(vv,from,w->localrep,to,&w->localupdate);
816:     PetscLogObjectParent((PetscObject)vv,(PetscObject)w->localupdate);
817:     ISDestroy(&to);
818:     ISDestroy(&from);

820:     /* set local to global mapping for ghosted vector */
821:     PetscMalloc1(n+nghost,&indices);
822:     VecGetOwnershipRange(vv,&rstart,NULL);

824:     for (i=0; i<n; i++)      indices[i]   = rstart + i;
825:     for (i=0; i<nghost; i++) indices[n+i] = ghosts[i];

827:     ISLocalToGlobalMappingCreate(comm,1,n+nghost,indices,PETSC_OWN_POINTER,&ltog);
828:     VecSetLocalToGlobalMapping(vv,ltog);
829:     ISLocalToGlobalMappingDestroy(&ltog);
830:   } else if (vv->ops->create == VecCreate_MPI) SETERRQ(PetscObjectComm((PetscObject)vv),PETSC_ERR_ARG_WRONGSTATE,"Must set local or global size before setting ghosting");
831:   else if (!((PetscObject)vv)->type_name) SETERRQ(PetscObjectComm((PetscObject)vv),PETSC_ERR_ARG_WRONGSTATE,"Must set type to VECMPI before ghosting");
832:   return(0);
833: }


836: /* ------------------------------------------------------------------------------------------*/
837: /*@C
838:    VecCreateGhostBlockWithArray - Creates a parallel vector with ghost padding on each processor;
839:    the caller allocates the array space. Indices in the ghost region are based on blocks.

841:    Collective

843:    Input Parameters:
844: +  comm - the MPI communicator to use
845: .  bs - block size
846: .  n - local vector length
847: .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
848: .  nghost - number of local ghost blocks
849: .  ghosts - global indices of ghost blocks (or NULL if not needed), counts are by block not by index, these do not need to be in increasing order (sorted)
850: -  array - the space to store the vector values (as long as n + nghost*bs)

852:    Output Parameter:
853: .  vv - the global vector representation (without ghost points as part of vector)

855:    Notes:
856:    Use VecGhostGetLocalForm() to access the local, ghosted representation
857:    of the vector.

859:    n is the local vector size (total local size not the number of blocks) while nghost
860:    is the number of blocks in the ghost portion, i.e. the number of elements in the ghost
861:    portion is bs*nghost

863:    Level: advanced


866: .seealso: VecCreate(), VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
867:           VecCreateGhost(), VecCreateSeqWithArray(), VecCreateMPIWithArray(),
868:           VecCreateGhostWithArray(), VecCreateGhostBlock()

870: @*/
871: PetscErrorCode  VecCreateGhostBlockWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],const PetscScalar array[],Vec *vv)
872: {
873:   PetscErrorCode         ierr;
874:   Vec_MPI                *w;
875:   PetscScalar            *larray;
876:   IS                     from,to;
877:   ISLocalToGlobalMapping ltog;
878:   PetscInt               rstart,i,nb,*indices;

881:   *vv = 0;

883:   if (n == PETSC_DECIDE)      SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size");
884:   if (nghost == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local ghost size");
885:   if (nghost < 0)             SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Ghost length must be >= 0");
886:   if (n % bs)                 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local size must be a multiple of block size");
887:   PetscSplitOwnership(comm,&n,&N);
888:   /* Create global representation */
889:   VecCreate(comm,vv);
890:   VecSetSizes(*vv,n,N);
891:   VecSetBlockSize(*vv,bs);
892:   VecCreate_MPI_Private(*vv,PETSC_TRUE,nghost*bs,array);
893:   w    = (Vec_MPI*)(*vv)->data;
894:   /* Create local representation */
895:   VecGetArray(*vv,&larray);
896:   VecCreateSeqWithArray(PETSC_COMM_SELF,bs,n+bs*nghost,larray,&w->localrep);
897:   PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localrep);
898:   VecRestoreArray(*vv,&larray);

900:   /*
901:        Create scatter context for scattering (updating) ghost values
902:   */
903:   ISCreateBlock(comm,bs,nghost,ghosts,PETSC_COPY_VALUES,&from);
904:   ISCreateStride(PETSC_COMM_SELF,bs*nghost,n,1,&to);
905:   VecScatterCreate(*vv,from,w->localrep,to,&w->localupdate);
906:   PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localupdate);
907:   ISDestroy(&to);
908:   ISDestroy(&from);

910:   /* set local to global mapping for ghosted vector */
911:   nb     = n/bs;
912:   PetscMalloc1(nb+nghost,&indices);
913:   VecGetOwnershipRange(*vv,&rstart,NULL);
914:   rstart = rstart/bs;

916:   for (i=0; i<nb; i++)      indices[i]    = rstart + i;
917:   for (i=0; i<nghost; i++)  indices[nb+i] = ghosts[i];

919:   ISLocalToGlobalMappingCreate(comm,bs,nb+nghost,indices,PETSC_OWN_POINTER,&ltog);
920:   VecSetLocalToGlobalMapping(*vv,ltog);
921:   ISLocalToGlobalMappingDestroy(&ltog);
922:   return(0);
923: }

925: /*@
926:    VecCreateGhostBlock - Creates a parallel vector with ghost padding on each processor.
927:         The indicing of the ghost points is done with blocks.

929:    Collective

931:    Input Parameters:
932: +  comm - the MPI communicator to use
933: .  bs - the block size
934: .  n - local vector length
935: .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
936: .  nghost - number of local ghost blocks
937: -  ghosts - global indices of ghost blocks, counts are by block, not by individual index, these do not need to be in increasing order (sorted)

939:    Output Parameter:
940: .  vv - the global vector representation (without ghost points as part of vector)

942:    Notes:
943:    Use VecGhostGetLocalForm() to access the local, ghosted representation
944:    of the vector.

946:    n is the local vector size (total local size not the number of blocks) while nghost
947:    is the number of blocks in the ghost portion, i.e. the number of elements in the ghost
948:    portion is bs*nghost

950:    Level: advanced

952: .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
953:           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
954:           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecCreateGhostBlockWithArray()

956: @*/
957: PetscErrorCode  VecCreateGhostBlock(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],Vec *vv)
958: {

962:   VecCreateGhostBlockWithArray(comm,bs,n,N,nghost,ghosts,0,vv);
963:   return(0);
964: }