xref: /petsc/src/vec/vec/impls/mpi/pbvec.c (revision e607c864f3adee196c0614801d5954446e31f8d0)
1 
2 /*
3    This file contains routines for Parallel vector operations.
4  */
5 #include <petscsys.h>
6 #include <../src/vec/vec/impls/mpi/pvecimpl.h>   /*I  "petscvec.h"   I*/
7 
8 PetscErrorCode VecDot_MPI(Vec xin,Vec yin,PetscScalar *z)
9 {
10   PetscScalar    sum,work;
11   PetscErrorCode ierr;
12 
13   PetscFunctionBegin;
14   ierr = VecDot_Seq(xin,yin,&work);CHKERRQ(ierr);
15   ierr = MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRMPI(ierr);
16   *z   = sum;
17   PetscFunctionReturn(0);
18 }
19 
20 PetscErrorCode VecTDot_MPI(Vec xin,Vec yin,PetscScalar *z)
21 {
22   PetscScalar    sum,work;
23   PetscErrorCode ierr;
24 
25   PetscFunctionBegin;
26   ierr = VecTDot_Seq(xin,yin,&work);CHKERRQ(ierr);
27   ierr = MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRMPI(ierr);
28   *z   = sum;
29   PetscFunctionReturn(0);
30 }
31 
32 extern PetscErrorCode VecView_MPI_Draw(Vec,PetscViewer);
33 
34 static PetscErrorCode VecPlaceArray_MPI(Vec vin,const PetscScalar *a)
35 {
36   PetscErrorCode ierr;
37   Vec_MPI        *v = (Vec_MPI*)vin->data;
38 
39   PetscFunctionBegin;
40   PetscCheckFalse(v->unplacedarray,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     ierr = VecPlaceArray(v->localrep,a);CHKERRQ(ierr);
45   }
46   PetscFunctionReturn(0);
47 }
48 
49 PetscErrorCode VecDuplicate_MPI(Vec win,Vec *v)
50 {
51   PetscErrorCode ierr;
52   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
53   PetscScalar    *array;
54 
55   PetscFunctionBegin;
56   ierr = VecCreate(PetscObjectComm((PetscObject)win),v);CHKERRQ(ierr);
57   ierr = PetscLayoutReference(win->map,&(*v)->map);CHKERRQ(ierr);
58 
59   ierr = VecCreate_MPI_Private(*v,PETSC_TRUE,w->nghost,NULL);CHKERRQ(ierr);
60   vw   = (Vec_MPI*)(*v)->data;
61   ierr = PetscMemcpy((*v)->ops,win->ops,sizeof(struct _VecOps));CHKERRQ(ierr);
62 
63   /* save local representation of the parallel vector (and scatter) if it exists */
64   if (w->localrep) {
65     ierr = VecGetArray(*v,&array);CHKERRQ(ierr);
66     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,PetscAbs(win->map->bs),win->map->n+w->nghost,array,&vw->localrep);CHKERRQ(ierr);
67     ierr = PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));CHKERRQ(ierr);
68     ierr = VecRestoreArray(*v,&array);CHKERRQ(ierr);
69     ierr = PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);CHKERRQ(ierr);
70 
71     vw->localupdate = w->localupdate;
72     if (vw->localupdate) {
73       ierr = PetscObjectReference((PetscObject)vw->localupdate);CHKERRQ(ierr);
74     }
75   }
76 
77   /* New vector should inherit stashing property of parent */
78   (*v)->stash.donotstash   = win->stash.donotstash;
79   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;
80 
81   ierr = PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);CHKERRQ(ierr);
82   ierr = PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);CHKERRQ(ierr);
83 
84   (*v)->map->bs   = PetscAbs(win->map->bs);
85   (*v)->bstash.bs = win->bstash.bs;
86   PetscFunctionReturn(0);
87 }
88 
89 static PetscErrorCode VecSetOption_MPI(Vec V,VecOption op,PetscBool flag)
90 {
91   Vec_MPI        *v = (Vec_MPI*)V->data;
92   PetscErrorCode ierr;
93 
94   PetscFunctionBegin;
95   switch (op) {
96   case VEC_IGNORE_OFF_PROC_ENTRIES: V->stash.donotstash = flag;
97     break;
98   case VEC_IGNORE_NEGATIVE_INDICES: V->stash.ignorenegidx = flag;
99     break;
100   case VEC_SUBSET_OFF_PROC_ENTRIES:
101     v->assembly_subset = flag; /* See the same logic in MatAssembly wrt MAT_SUBSET_OFF_PROC_ENTRIES */
102     if (!v->assembly_subset) { /* User indicates "do not reuse the communication pattern" */
103       ierr = VecAssemblyReset_MPI(V);CHKERRQ(ierr); /* Reset existing pattern to free memory */
104       v->first_assembly_done = PETSC_FALSE; /* Mark the first assembly is not done */
105     }
106     break;
107   }
108 
109   PetscFunctionReturn(0);
110 }
111 
112 static PetscErrorCode VecResetArray_MPI(Vec vin)
113 {
114   Vec_MPI        *v = (Vec_MPI*)vin->data;
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   v->array         = v->unplacedarray;
119   v->unplacedarray = NULL;
120   if (v->localrep) {
121     ierr = VecResetArray(v->localrep);CHKERRQ(ierr);
122   }
123   PetscFunctionReturn(0);
124 }
125 
126 static PetscErrorCode VecAssemblySend_MPI_Private(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rankid,PetscMPIInt rank,void *sdata,MPI_Request req[],void *ctx)
127 {
128   Vec X = (Vec)ctx;
129   Vec_MPI *x = (Vec_MPI*)X->data;
130   VecAssemblyHeader *hdr = (VecAssemblyHeader*)sdata;
131   PetscInt bs = X->map->bs;
132   PetscErrorCode ierr;
133 
134   PetscFunctionBegin;
135   /* x->first_assembly_done indicates we are reusing a communication network. In that case, some
136      messages can be empty, but we have to send them this time if we sent them before because the
137      receiver is expecting them.
138    */
139   if (hdr->count || (x->first_assembly_done && x->sendptrs[rankid].ints)) {
140     ierr = MPI_Isend(x->sendptrs[rankid].ints,hdr->count,MPIU_INT,rank,tag[0],comm,&req[0]);CHKERRMPI(ierr);
141     ierr = MPI_Isend(x->sendptrs[rankid].scalars,hdr->count,MPIU_SCALAR,rank,tag[1],comm,&req[1]);CHKERRMPI(ierr);
142   }
143   if (hdr->bcount || (x->first_assembly_done && x->sendptrs[rankid].intb)) {
144     ierr = MPI_Isend(x->sendptrs[rankid].intb,hdr->bcount,MPIU_INT,rank,tag[2],comm,&req[2]);CHKERRMPI(ierr);
145     ierr = MPI_Isend(x->sendptrs[rankid].scalarb,hdr->bcount*bs,MPIU_SCALAR,rank,tag[3],comm,&req[3]);CHKERRMPI(ierr);
146   }
147   PetscFunctionReturn(0);
148 }
149 
150 static PetscErrorCode VecAssemblyRecv_MPI_Private(MPI_Comm comm,const PetscMPIInt tag[],PetscMPIInt rank,void *rdata,MPI_Request req[],void *ctx)
151 {
152   Vec X = (Vec)ctx;
153   Vec_MPI *x = (Vec_MPI*)X->data;
154   VecAssemblyHeader *hdr = (VecAssemblyHeader*)rdata;
155   PetscErrorCode ierr;
156   PetscInt bs = X->map->bs;
157   VecAssemblyFrame *frame;
158 
159   PetscFunctionBegin;
160   ierr = PetscSegBufferGet(x->segrecvframe,1,&frame);CHKERRQ(ierr);
161 
162   if (hdr->count) {
163     ierr = PetscSegBufferGet(x->segrecvint,hdr->count,&frame->ints);CHKERRQ(ierr);
164     ierr = MPI_Irecv(frame->ints,hdr->count,MPIU_INT,rank,tag[0],comm,&req[0]);CHKERRMPI(ierr);
165     ierr = PetscSegBufferGet(x->segrecvscalar,hdr->count,&frame->scalars);CHKERRQ(ierr);
166     ierr = MPI_Irecv(frame->scalars,hdr->count,MPIU_SCALAR,rank,tag[1],comm,&req[1]);CHKERRMPI(ierr);
167     frame->pendings = 2;
168   } else {
169     frame->ints = NULL;
170     frame->scalars = NULL;
171     frame->pendings = 0;
172   }
173 
174   if (hdr->bcount) {
175     ierr = PetscSegBufferGet(x->segrecvint,hdr->bcount,&frame->intb);CHKERRQ(ierr);
176     ierr = MPI_Irecv(frame->intb,hdr->bcount,MPIU_INT,rank,tag[2],comm,&req[2]);CHKERRMPI(ierr);
177     ierr = PetscSegBufferGet(x->segrecvscalar,hdr->bcount*bs,&frame->scalarb);CHKERRQ(ierr);
178     ierr = MPI_Irecv(frame->scalarb,hdr->bcount*bs,MPIU_SCALAR,rank,tag[3],comm,&req[3]);CHKERRMPI(ierr);
179     frame->pendingb = 2;
180   } else {
181     frame->intb = NULL;
182     frame->scalarb = NULL;
183     frame->pendingb = 0;
184   }
185   PetscFunctionReturn(0);
186 }
187 
188 static PetscErrorCode VecAssemblyBegin_MPI_BTS(Vec X)
189 {
190   Vec_MPI        *x = (Vec_MPI*)X->data;
191   PetscErrorCode ierr;
192   MPI_Comm       comm;
193   PetscInt       i,j,jb,bs;
194 
195   PetscFunctionBegin;
196   if (X->stash.donotstash) PetscFunctionReturn(0);
197 
198   ierr = PetscObjectGetComm((PetscObject)X,&comm);CHKERRQ(ierr);
199   ierr = VecGetBlockSize(X,&bs);CHKERRQ(ierr);
200   if (PetscDefined(USE_DEBUG)) {
201     InsertMode addv;
202     ierr = MPIU_Allreduce((PetscEnum*)&X->stash.insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,comm);CHKERRMPI(ierr);
203     PetscCheckFalse(addv == (ADD_VALUES|INSERT_VALUES),comm,PETSC_ERR_ARG_NOTSAMETYPE,"Some processors inserted values while others added");
204   }
205   X->bstash.insertmode = X->stash.insertmode; /* Block stash implicitly tracks InsertMode of scalar stash */
206 
207   ierr = VecStashSortCompress_Private(&X->stash);CHKERRQ(ierr);
208   ierr = VecStashSortCompress_Private(&X->bstash);CHKERRQ(ierr);
209 
210   if (!x->sendranks) {
211     PetscMPIInt nowners,bnowners,*owners,*bowners;
212     PetscInt ntmp;
213     ierr = VecStashGetOwnerList_Private(&X->stash,X->map,&nowners,&owners);CHKERRQ(ierr);
214     ierr = VecStashGetOwnerList_Private(&X->bstash,X->map,&bnowners,&bowners);CHKERRQ(ierr);
215     ierr = PetscMergeMPIIntArray(nowners,owners,bnowners,bowners,&ntmp,&x->sendranks);CHKERRQ(ierr);
216     x->nsendranks = ntmp;
217     ierr = PetscFree(owners);CHKERRQ(ierr);
218     ierr = PetscFree(bowners);CHKERRQ(ierr);
219     ierr = PetscMalloc1(x->nsendranks,&x->sendhdr);CHKERRQ(ierr);
220     ierr = PetscCalloc1(x->nsendranks,&x->sendptrs);CHKERRQ(ierr);
221   }
222   for (i=0,j=0,jb=0; i<x->nsendranks; i++) {
223     PetscMPIInt rank = x->sendranks[i];
224     x->sendhdr[i].insertmode = X->stash.insertmode;
225     /* Initialize pointers for non-empty stashes the first time around.  Subsequent assemblies with
226      * VEC_SUBSET_OFF_PROC_ENTRIES will leave the old pointers (dangling because the stash has been collected) when
227      * there is nothing new to send, so that size-zero messages get sent instead. */
228     x->sendhdr[i].count = 0;
229     if (X->stash.n) {
230       x->sendptrs[i].ints    = &X->stash.idx[j];
231       x->sendptrs[i].scalars = &X->stash.array[j];
232       for (; j<X->stash.n && X->stash.idx[j] < X->map->range[rank+1]; j++) x->sendhdr[i].count++;
233     }
234     x->sendhdr[i].bcount = 0;
235     if (X->bstash.n) {
236       x->sendptrs[i].intb    = &X->bstash.idx[jb];
237       x->sendptrs[i].scalarb = &X->bstash.array[jb*bs];
238       for (; jb<X->bstash.n && X->bstash.idx[jb]*bs < X->map->range[rank+1]; jb++) x->sendhdr[i].bcount++;
239     }
240   }
241 
242   if (!x->segrecvint) {ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&x->segrecvint);CHKERRQ(ierr);}
243   if (!x->segrecvscalar) {ierr = PetscSegBufferCreate(sizeof(PetscScalar),1000,&x->segrecvscalar);CHKERRQ(ierr);}
244   if (!x->segrecvframe) {ierr = PetscSegBufferCreate(sizeof(VecAssemblyFrame),50,&x->segrecvframe);CHKERRQ(ierr);}
245   if (x->first_assembly_done) { /* this is not the first assembly */
246     PetscMPIInt tag[4];
247     for (i=0; i<4; i++) {ierr = PetscCommGetNewTag(comm,&tag[i]);CHKERRQ(ierr);}
248     for (i=0; i<x->nsendranks; i++) {
249       ierr = VecAssemblySend_MPI_Private(comm,tag,i,x->sendranks[i],x->sendhdr+i,x->sendreqs+4*i,X);CHKERRQ(ierr);
250     }
251     for (i=0; i<x->nrecvranks; i++) {
252       ierr = VecAssemblyRecv_MPI_Private(comm,tag,x->recvranks[i],x->recvhdr+i,x->recvreqs+4*i,X);CHKERRQ(ierr);
253     }
254     x->use_status = PETSC_TRUE;
255   } else { /* First time assembly */
256     ierr = 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);CHKERRQ(ierr);
257     x->use_status = PETSC_FALSE;
258   }
259 
260   /* The first_assembly_done flag is only meaningful when x->assembly_subset is set.
261      This line says when assembly_subset is set, then we mark that the first assembly is done.
262    */
263   x->first_assembly_done = x->assembly_subset;
264 
265   {
266     PetscInt nstash,reallocs;
267     ierr = VecStashGetInfo_Private(&X->stash,&nstash,&reallocs);CHKERRQ(ierr);
268     ierr = PetscInfo(X,"Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
269     ierr = VecStashGetInfo_Private(&X->bstash,&nstash,&reallocs);CHKERRQ(ierr);
270     ierr = PetscInfo(X,"Block-Stash has %" PetscInt_FMT " entries, uses %" PetscInt_FMT " mallocs.\n",nstash,reallocs);CHKERRQ(ierr);
271   }
272   PetscFunctionReturn(0);
273 }
274 
275 static PetscErrorCode VecAssemblyEnd_MPI_BTS(Vec X)
276 {
277   Vec_MPI *x = (Vec_MPI*)X->data;
278   PetscInt bs = X->map->bs;
279   PetscMPIInt npending,*some_indices,r;
280   MPI_Status  *some_statuses;
281   PetscScalar *xarray;
282   PetscErrorCode ierr;
283   VecAssemblyFrame *frame;
284 
285   PetscFunctionBegin;
286   if (X->stash.donotstash) {
287     X->stash.insertmode = NOT_SET_VALUES;
288     X->bstash.insertmode = NOT_SET_VALUES;
289     PetscFunctionReturn(0);
290   }
291 
292   PetscCheckFalse(!x->segrecvframe,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing segrecvframe! Probably you forgot to call VecAssemblyBegin first");
293   ierr = VecGetArray(X,&xarray);CHKERRQ(ierr);
294   ierr = PetscSegBufferExtractInPlace(x->segrecvframe,&frame);CHKERRQ(ierr);
295   ierr = PetscMalloc2(4*x->nrecvranks,&some_indices,x->use_status?4*x->nrecvranks:0,&some_statuses);CHKERRQ(ierr);
296   for (r=0,npending=0; r<x->nrecvranks; r++) npending += frame[r].pendings + frame[r].pendingb;
297   while (npending>0) {
298     PetscMPIInt ndone=0,ii;
299     /* Filling MPI_Status fields requires some resources from the MPI library.  We skip it on the first assembly, or
300      * when VEC_SUBSET_OFF_PROC_ENTRIES has not been set, because we could exchange exact sizes in the initial
301      * rendezvous.  When the rendezvous is elided, however, we use MPI_Status to get actual message lengths, so that
302      * subsequent assembly can set a proper subset of the values. */
303     ierr = MPI_Waitsome(4*x->nrecvranks,x->recvreqs,&ndone,some_indices,x->use_status?some_statuses:MPI_STATUSES_IGNORE);CHKERRMPI(ierr);
304     for (ii=0; ii<ndone; ii++) {
305       PetscInt i = some_indices[ii]/4,j,k;
306       InsertMode imode = (InsertMode)x->recvhdr[i].insertmode;
307       PetscInt *recvint;
308       PetscScalar *recvscalar;
309       PetscBool intmsg = (PetscBool)(some_indices[ii]%2 == 0);
310       PetscBool blockmsg = (PetscBool)((some_indices[ii]%4)/2 == 1);
311       npending--;
312       if (!blockmsg) { /* Scalar stash */
313         PetscMPIInt count;
314         if (--frame[i].pendings > 0) continue;
315         if (x->use_status) {
316           ierr = MPI_Get_count(&some_statuses[ii],intmsg ? MPIU_INT : MPIU_SCALAR,&count);CHKERRMPI(ierr);
317         } else count = x->recvhdr[i].count;
318         for (j=0,recvint=frame[i].ints,recvscalar=frame[i].scalars; j<count; j++,recvint++) {
319           PetscInt loc = *recvint - X->map->rstart;
320           PetscCheckFalse(*recvint < X->map->rstart || X->map->rend <= *recvint,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Received vector entry %" PetscInt_FMT " out of local range [%" PetscInt_FMT ",%" PetscInt_FMT ")]",*recvint,X->map->rstart,X->map->rend);
321           switch (imode) {
322           case ADD_VALUES:
323             xarray[loc] += *recvscalar++;
324             break;
325           case INSERT_VALUES:
326             xarray[loc] = *recvscalar++;
327             break;
328           default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Insert mode not supported 0x%x",imode);
329           }
330         }
331       } else {                  /* Block stash */
332         PetscMPIInt count;
333         if (--frame[i].pendingb > 0) continue;
334         if (x->use_status) {
335           ierr = MPI_Get_count(&some_statuses[ii],intmsg ? MPIU_INT : MPIU_SCALAR,&count);CHKERRMPI(ierr);
336           if (!intmsg) count /= bs; /* Convert from number of scalars to number of blocks */
337         } else count = x->recvhdr[i].bcount;
338         for (j=0,recvint=frame[i].intb,recvscalar=frame[i].scalarb; j<count; j++,recvint++) {
339           PetscInt loc = (*recvint)*bs - X->map->rstart;
340           switch (imode) {
341           case ADD_VALUES:
342             for (k=loc; k<loc+bs; k++) xarray[k] += *recvscalar++;
343             break;
344           case INSERT_VALUES:
345             for (k=loc; k<loc+bs; k++) xarray[k] = *recvscalar++;
346             break;
347           default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Insert mode not supported 0x%x",imode);
348           }
349         }
350       }
351     }
352   }
353   ierr = VecRestoreArray(X,&xarray);CHKERRQ(ierr);
354   ierr = MPI_Waitall(4*x->nsendranks,x->sendreqs,MPI_STATUSES_IGNORE);CHKERRMPI(ierr);
355   ierr = PetscFree2(some_indices,some_statuses);CHKERRQ(ierr);
356   if (x->assembly_subset) {
357     void *dummy;                /* reset segbuffers */
358     ierr = PetscSegBufferExtractInPlace(x->segrecvint,&dummy);CHKERRQ(ierr);
359     ierr = PetscSegBufferExtractInPlace(x->segrecvscalar,&dummy);CHKERRQ(ierr);
360   } else {
361     ierr = VecAssemblyReset_MPI(X);CHKERRQ(ierr);
362   }
363 
364   X->stash.insertmode = NOT_SET_VALUES;
365   X->bstash.insertmode = NOT_SET_VALUES;
366   ierr = VecStashScatterEnd_Private(&X->stash);CHKERRQ(ierr);
367   ierr = VecStashScatterEnd_Private(&X->bstash);CHKERRQ(ierr);
368   PetscFunctionReturn(0);
369 }
370 
371 PetscErrorCode VecAssemblyReset_MPI(Vec X)
372 {
373   Vec_MPI *x = (Vec_MPI*)X->data;
374   PetscErrorCode ierr;
375 
376   PetscFunctionBegin;
377   ierr = PetscFree(x->sendreqs);CHKERRQ(ierr);
378   ierr = PetscFree(x->recvreqs);CHKERRQ(ierr);
379   ierr = PetscFree(x->sendranks);CHKERRQ(ierr);
380   ierr = PetscFree(x->recvranks);CHKERRQ(ierr);
381   ierr = PetscFree(x->sendhdr);CHKERRQ(ierr);
382   ierr = PetscFree(x->recvhdr);CHKERRQ(ierr);
383   ierr = PetscFree(x->sendptrs);CHKERRQ(ierr);
384   ierr = PetscSegBufferDestroy(&x->segrecvint);CHKERRQ(ierr);
385   ierr = PetscSegBufferDestroy(&x->segrecvscalar);CHKERRQ(ierr);
386   ierr = PetscSegBufferDestroy(&x->segrecvframe);CHKERRQ(ierr);
387   PetscFunctionReturn(0);
388 }
389 
390 static PetscErrorCode VecSetFromOptions_MPI(PetscOptionItems *PetscOptionsObject,Vec X)
391 {
392 #if !defined(PETSC_HAVE_MPIUNI)
393   PetscErrorCode ierr;
394   PetscBool      flg = PETSC_FALSE,set;
395 
396   PetscFunctionBegin;
397   ierr = PetscOptionsHead(PetscOptionsObject,"VecMPI Options");CHKERRQ(ierr);
398   ierr = PetscOptionsBool("-vec_assembly_legacy","Use MPI 1 version of assembly","",flg,&flg,&set);CHKERRQ(ierr);
399   if (set) {
400     X->ops->assemblybegin = flg ? VecAssemblyBegin_MPI : VecAssemblyBegin_MPI_BTS;
401     X->ops->assemblyend   = flg ? VecAssemblyEnd_MPI   : VecAssemblyEnd_MPI_BTS;
402   }
403   ierr = PetscOptionsTail();CHKERRQ(ierr);
404 #else
405   PetscFunctionBegin;
406   X->ops->assemblybegin = VecAssemblyBegin_MPI;
407   X->ops->assemblyend   = VecAssemblyEnd_MPI;
408 #endif
409   PetscFunctionReturn(0);
410 }
411 
412 static struct _VecOps DvOps = {
413   PetscDesignatedInitializer(duplicate,VecDuplicate_MPI), /* 1 */
414   PetscDesignatedInitializer(duplicatevecs,VecDuplicateVecs_Default),
415   PetscDesignatedInitializer(destroyvecs,VecDestroyVecs_Default),
416   PetscDesignatedInitializer(dot,VecDot_MPI),
417   PetscDesignatedInitializer(mdot,VecMDot_MPI),
418   PetscDesignatedInitializer(norm,VecNorm_MPI),
419   PetscDesignatedInitializer(tdot,VecTDot_MPI),
420   PetscDesignatedInitializer(mtdot,VecMTDot_MPI),
421   PetscDesignatedInitializer(scale,VecScale_Seq),
422   PetscDesignatedInitializer(copy,VecCopy_Seq), /* 10 */
423   PetscDesignatedInitializer(set,VecSet_Seq),
424   PetscDesignatedInitializer(swap,VecSwap_Seq),
425   PetscDesignatedInitializer(axpy,VecAXPY_Seq),
426   PetscDesignatedInitializer(axpby,VecAXPBY_Seq),
427   PetscDesignatedInitializer(maxpy,VecMAXPY_Seq),
428   PetscDesignatedInitializer(aypx,VecAYPX_Seq),
429   PetscDesignatedInitializer(waxpy,VecWAXPY_Seq),
430   PetscDesignatedInitializer(axpbypcz,VecAXPBYPCZ_Seq),
431   PetscDesignatedInitializer(pointwisemult,VecPointwiseMult_Seq),
432   PetscDesignatedInitializer(pointwisedivide,VecPointwiseDivide_Seq),
433   PetscDesignatedInitializer(setvalues,VecSetValues_MPI), /* 20 */
434   PetscDesignatedInitializer(assemblybegin,VecAssemblyBegin_MPI_BTS),
435   PetscDesignatedInitializer(assemblyend,VecAssemblyEnd_MPI_BTS),
436   PetscDesignatedInitializer(getarray,NULL),
437   PetscDesignatedInitializer(getsize,VecGetSize_MPI),
438   PetscDesignatedInitializer(getlocalsize,VecGetSize_Seq),
439   PetscDesignatedInitializer(restorearray,NULL),
440   PetscDesignatedInitializer(max,VecMax_MPI),
441   PetscDesignatedInitializer(min,VecMin_MPI),
442   PetscDesignatedInitializer(setrandom,VecSetRandom_Seq),
443   PetscDesignatedInitializer(setoption,VecSetOption_MPI),
444   PetscDesignatedInitializer(setvaluesblocked,VecSetValuesBlocked_MPI),
445   PetscDesignatedInitializer(destroy,VecDestroy_MPI),
446   PetscDesignatedInitializer(view,VecView_MPI),
447   PetscDesignatedInitializer(placearray,VecPlaceArray_MPI),
448   PetscDesignatedInitializer(replacearray,VecReplaceArray_Seq),
449   PetscDesignatedInitializer(dot_local,VecDot_Seq),
450   PetscDesignatedInitializer(tdot_local,VecTDot_Seq),
451   PetscDesignatedInitializer(norm_local,VecNorm_Seq),
452   PetscDesignatedInitializer(mdot_local,VecMDot_Seq),
453   PetscDesignatedInitializer(mtdot_local,VecMTDot_Seq),
454   PetscDesignatedInitializer(load,VecLoad_Default),
455   PetscDesignatedInitializer(reciprocal,VecReciprocal_Default),
456   PetscDesignatedInitializer(conjugate,VecConjugate_Seq),
457   PetscDesignatedInitializer(setlocaltoglobalmapping,NULL),
458   PetscDesignatedInitializer(setvalueslocal,NULL),
459   PetscDesignatedInitializer(resetarray,VecResetArray_MPI),
460   PetscDesignatedInitializer(setfromoptions,VecSetFromOptions_MPI),/*set from options */
461   PetscDesignatedInitializer(maxpointwisedivide,VecMaxPointwiseDivide_Seq),
462   PetscDesignatedInitializer(pointwisemax,VecPointwiseMax_Seq),
463   PetscDesignatedInitializer(pointwisemaxabs,VecPointwiseMaxAbs_Seq),
464   PetscDesignatedInitializer(pointwisemin,VecPointwiseMin_Seq),
465   PetscDesignatedInitializer(getvalues,VecGetValues_MPI),
466   PetscDesignatedInitializer(sqrt,NULL),
467   PetscDesignatedInitializer(abs,NULL),
468   PetscDesignatedInitializer(exp,NULL),
469   PetscDesignatedInitializer(log,NULL),
470   PetscDesignatedInitializer(shift,NULL),
471   PetscDesignatedInitializer(create,NULL), /* really? */
472   PetscDesignatedInitializer(stridegather,VecStrideGather_Default),
473   PetscDesignatedInitializer(stridescatter,VecStrideScatter_Default),
474   PetscDesignatedInitializer(dotnorm2,NULL),
475   PetscDesignatedInitializer(getsubvector,NULL),
476   PetscDesignatedInitializer(restoresubvector,NULL),
477   PetscDesignatedInitializer(getarrayread,NULL),
478   PetscDesignatedInitializer(restorearrayread,NULL),
479   PetscDesignatedInitializer(stridesubsetgather,VecStrideSubSetGather_Default),
480   PetscDesignatedInitializer(stridesubsetscatter,VecStrideSubSetScatter_Default),
481   PetscDesignatedInitializer(viewnative,VecView_MPI),
482   PetscDesignatedInitializer(loadnative,NULL),
483   PetscDesignatedInitializer(getlocalvector,NULL),
484 };
485 
486 /*
487     VecCreate_MPI_Private - Basic create routine called by VecCreate_MPI() (i.e. VecCreateMPI()),
488     VecCreateMPIWithArray(), VecCreate_Shared() (i.e. VecCreateShared()), VecCreateGhost(),
489     VecDuplicate_MPI(), VecCreateGhostWithArray(), VecDuplicate_MPI(), and VecDuplicate_Shared()
490 
491     If alloc is true and array is NULL then this routine allocates the space, otherwise
492     no space is allocated.
493 */
494 PetscErrorCode VecCreate_MPI_Private(Vec v,PetscBool alloc,PetscInt nghost,const PetscScalar array[])
495 {
496   Vec_MPI        *s;
497   PetscErrorCode ierr;
498 
499   PetscFunctionBegin;
500   ierr           = PetscNewLog(v,&s);CHKERRQ(ierr);
501   v->data        = (void*)s;
502   ierr           = PetscMemcpy(v->ops,&DvOps,sizeof(DvOps));CHKERRQ(ierr);
503   s->nghost      = nghost;
504   v->petscnative = PETSC_TRUE;
505   if (array) v->offloadmask = PETSC_OFFLOAD_CPU;
506 
507   ierr = PetscLayoutSetUp(v->map);CHKERRQ(ierr);
508 
509   s->array           = (PetscScalar*)array;
510   s->array_allocated = NULL;
511   if (alloc && !array) {
512     PetscInt n = v->map->n+nghost;
513     ierr               = PetscCalloc1(n,&s->array);CHKERRQ(ierr);
514     ierr               = PetscLogObjectMemory((PetscObject)v,n*sizeof(PetscScalar));CHKERRQ(ierr);
515     s->array_allocated = s->array;
516   }
517 
518   /* By default parallel vectors do not have local representation */
519   s->localrep    = NULL;
520   s->localupdate = NULL;
521 
522   v->stash.insertmode = NOT_SET_VALUES;
523   v->bstash.insertmode = NOT_SET_VALUES;
524   /* create the stashes. The block-size for bstash is set later when
525      VecSetValuesBlocked is called.
526   */
527   ierr = VecStashCreate_Private(PetscObjectComm((PetscObject)v),1,&v->stash);CHKERRQ(ierr);
528   ierr = VecStashCreate_Private(PetscObjectComm((PetscObject)v),PetscAbs(v->map->bs),&v->bstash);CHKERRQ(ierr);
529 
530 #if defined(PETSC_HAVE_MATLAB_ENGINE)
531   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEnginePut_C",VecMatlabEnginePut_Default);CHKERRQ(ierr);
532   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscMatlabEngineGet_C",VecMatlabEngineGet_Default);CHKERRQ(ierr);
533 #endif
534   ierr = PetscObjectChangeTypeName((PetscObject)v,VECMPI);CHKERRQ(ierr);
535   PetscFunctionReturn(0);
536 }
537 
538 /*MC
539    VECMPI - VECMPI = "mpi" - The basic parallel vector
540 
541    Options Database Keys:
542 . -vec_type mpi - sets the vector type to VECMPI during a call to VecSetFromOptions()
543 
544   Level: beginner
545 
546 .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECMPI, VecType, VecCreateMPI(), VecCreateMPI()
547 M*/
548 
549 PetscErrorCode VecCreate_MPI(Vec vv)
550 {
551   PetscErrorCode ierr;
552 
553   PetscFunctionBegin;
554   ierr = VecCreate_MPI_Private(vv,PETSC_TRUE,0,NULL);CHKERRQ(ierr);
555   PetscFunctionReturn(0);
556 }
557 
558 /*MC
559    VECSTANDARD = "standard" - A VECSEQ on one process and VECMPI on more than one process
560 
561    Options Database Keys:
562 . -vec_type standard - sets a vector type to standard on calls to VecSetFromOptions()
563 
564   Level: beginner
565 
566 .seealso: VecCreateSeq(), VecCreateMPI()
567 M*/
568 
569 PETSC_EXTERN PetscErrorCode VecCreate_Standard(Vec v)
570 {
571   PetscErrorCode ierr;
572   PetscMPIInt    size;
573 
574   PetscFunctionBegin;
575   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);CHKERRMPI(ierr);
576   if (size == 1) {
577     ierr = VecSetType(v,VECSEQ);CHKERRQ(ierr);
578   } else {
579     ierr = VecSetType(v,VECMPI);CHKERRQ(ierr);
580   }
581   PetscFunctionReturn(0);
582 }
583 
584 /*@C
585    VecCreateMPIWithArray - Creates a parallel, array-style vector,
586    where the user provides the array space to store the vector values.
587 
588    Collective
589 
590    Input Parameters:
591 +  comm  - the MPI communicator to use
592 .  bs    - block size, same meaning as VecSetBlockSize()
593 .  n     - local vector length, cannot be PETSC_DECIDE
594 .  N     - global vector length (or PETSC_DECIDE to have calculated)
595 -  array - the user provided array to store the vector values
596 
597    Output Parameter:
598 .  vv - the vector
599 
600    Notes:
601    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
602    same type as an existing vector.
603 
604    If the user-provided array is NULL, then VecPlaceArray() can be used
605    at a later stage to SET the array for storing the vector values.
606 
607    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
608    The user should not free the array until the vector is destroyed.
609 
610    Level: intermediate
611 
612 .seealso: VecCreateSeqWithArray(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
613           VecCreateMPI(), VecCreateGhostWithArray(), VecPlaceArray()
614 
615 @*/
616 PetscErrorCode  VecCreateMPIWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar array[],Vec *vv)
617 {
618   PetscErrorCode ierr;
619 
620   PetscFunctionBegin;
621   PetscCheckFalse(n == PETSC_DECIDE,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size of vector");
622   ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
623   ierr = VecCreate(comm,vv);CHKERRQ(ierr);
624   ierr = VecSetSizes(*vv,n,N);CHKERRQ(ierr);
625   ierr = VecSetBlockSize(*vv,bs);CHKERRQ(ierr);
626   ierr = VecCreate_MPI_Private(*vv,PETSC_FALSE,0,array);CHKERRQ(ierr);
627   PetscFunctionReturn(0);
628 }
629 
630 /*@C
631    VecCreateGhostWithArray - Creates a parallel vector with ghost padding on each processor;
632    the caller allocates the array space.
633 
634    Collective
635 
636    Input Parameters:
637 +  comm - the MPI communicator to use
638 .  n - local vector length
639 .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
640 .  nghost - number of local ghost points
641 .  ghosts - global indices of ghost points (or NULL if not needed), these do not need to be in increasing order (sorted)
642 -  array - the space to store the vector values (as long as n + nghost)
643 
644    Output Parameter:
645 .  vv - the global vector representation (without ghost points as part of vector)
646 
647    Notes:
648    Use VecGhostGetLocalForm() to access the local, ghosted representation
649    of the vector.
650 
651    This also automatically sets the ISLocalToGlobalMapping() for this vector.
652 
653    Level: advanced
654 
655 .seealso: VecCreate(), VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
656           VecCreateGhost(), VecCreateSeqWithArray(), VecCreateMPIWithArray(),
657           VecCreateGhostBlock(), VecCreateGhostBlockWithArray(), VecMPISetGhost()
658 
659 @*/
660 PetscErrorCode  VecCreateGhostWithArray(MPI_Comm comm,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],const PetscScalar array[],Vec *vv)
661 {
662   PetscErrorCode         ierr;
663   Vec_MPI                *w;
664   PetscScalar            *larray;
665   IS                     from,to;
666   ISLocalToGlobalMapping ltog;
667   PetscInt               rstart,i,*indices;
668 
669   PetscFunctionBegin;
670   *vv = NULL;
671 
672   PetscCheckFalse(n == PETSC_DECIDE,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size");
673   PetscCheckFalse(nghost == PETSC_DECIDE,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local ghost size");
674   PetscCheckFalse(nghost < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Ghost length must be >= 0");
675   ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
676   /* Create global representation */
677   ierr = VecCreate(comm,vv);CHKERRQ(ierr);
678   ierr = VecSetSizes(*vv,n,N);CHKERRQ(ierr);
679   ierr = VecCreate_MPI_Private(*vv,PETSC_TRUE,nghost,array);CHKERRQ(ierr);
680   w    = (Vec_MPI*)(*vv)->data;
681   /* Create local representation */
682   ierr = VecGetArray(*vv,&larray);CHKERRQ(ierr);
683   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,n+nghost,larray,&w->localrep);CHKERRQ(ierr);
684   ierr = PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localrep);CHKERRQ(ierr);
685   ierr = VecRestoreArray(*vv,&larray);CHKERRQ(ierr);
686 
687   /*
688        Create scatter context for scattering (updating) ghost values
689   */
690   ierr = ISCreateGeneral(comm,nghost,ghosts,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
691   ierr = ISCreateStride(PETSC_COMM_SELF,nghost,n,1,&to);CHKERRQ(ierr);
692   ierr = VecScatterCreate(*vv,from,w->localrep,to,&w->localupdate);CHKERRQ(ierr);
693   ierr = PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localupdate);CHKERRQ(ierr);
694   ierr = ISDestroy(&to);CHKERRQ(ierr);
695   ierr = ISDestroy(&from);CHKERRQ(ierr);
696 
697   /* set local to global mapping for ghosted vector */
698   ierr = PetscMalloc1(n+nghost,&indices);CHKERRQ(ierr);
699   ierr = VecGetOwnershipRange(*vv,&rstart,NULL);CHKERRQ(ierr);
700   for (i=0; i<n; i++) {
701     indices[i] = rstart + i;
702   }
703   for (i=0; i<nghost; i++) {
704     indices[n+i] = ghosts[i];
705   }
706   ierr = ISLocalToGlobalMappingCreate(comm,1,n+nghost,indices,PETSC_OWN_POINTER,&ltog);CHKERRQ(ierr);
707   ierr = VecSetLocalToGlobalMapping(*vv,ltog);CHKERRQ(ierr);
708   ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
709   PetscFunctionReturn(0);
710 }
711 
712 /*@
713    VecCreateGhost - Creates a parallel vector with ghost padding on each processor.
714 
715    Collective
716 
717    Input Parameters:
718 +  comm - the MPI communicator to use
719 .  n - local vector length
720 .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
721 .  nghost - number of local ghost points
722 -  ghosts - global indices of ghost points, these do not need to be in increasing order (sorted)
723 
724    Output Parameter:
725 .  vv - the global vector representation (without ghost points as part of vector)
726 
727    Notes:
728    Use VecGhostGetLocalForm() to access the local, ghosted representation
729    of the vector.
730 
731    This also automatically sets the ISLocalToGlobalMapping() for this vector.
732 
733    Level: advanced
734 
735 .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
736           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(), VecGhostUpdateBegin(),
737           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecGhostUpdateEnd(),
738           VecCreateGhostBlock(), VecCreateGhostBlockWithArray(), VecMPISetGhost()
739 
740 @*/
741 PetscErrorCode  VecCreateGhost(MPI_Comm comm,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],Vec *vv)
742 {
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   ierr = VecCreateGhostWithArray(comm,n,N,nghost,ghosts,NULL,vv);CHKERRQ(ierr);
747   PetscFunctionReturn(0);
748 }
749 
750 /*@
751    VecMPISetGhost - Sets the ghost points for an MPI ghost vector
752 
753    Collective on Vec
754 
755    Input Parameters:
756 +  vv - the MPI vector
757 .  nghost - number of local ghost points
758 -  ghosts - global indices of ghost points, these do not need to be in increasing order (sorted)
759 
760    Notes:
761    Use VecGhostGetLocalForm() to access the local, ghosted representation
762    of the vector.
763 
764    This also automatically sets the ISLocalToGlobalMapping() for this vector.
765 
766    You must call this AFTER you have set the type of the vector (with VecSetType()) and the size (with VecSetSizes()).
767 
768    Level: advanced
769 
770 .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
771           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(), VecGhostUpdateBegin(),
772           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecGhostUpdateEnd(),
773           VecCreateGhostBlock(), VecCreateGhostBlockWithArray()
774 
775 @*/
776 PetscErrorCode  VecMPISetGhost(Vec vv,PetscInt nghost,const PetscInt ghosts[])
777 {
778   PetscErrorCode ierr;
779   PetscBool      flg;
780 
781   PetscFunctionBegin;
782   ierr = PetscObjectTypeCompare((PetscObject)vv,VECMPI,&flg);CHKERRQ(ierr);
783   /* if already fully existant VECMPI then basically destroy it and rebuild with ghosting */
784   if (flg) {
785     PetscInt               n,N;
786     Vec_MPI                *w;
787     PetscScalar            *larray;
788     IS                     from,to;
789     ISLocalToGlobalMapping ltog;
790     PetscInt               rstart,i,*indices;
791     MPI_Comm               comm;
792 
793     ierr = PetscObjectGetComm((PetscObject)vv,&comm);CHKERRQ(ierr);
794     n    = vv->map->n;
795     N    = vv->map->N;
796     ierr = (*vv->ops->destroy)(vv);CHKERRQ(ierr);
797     ierr = VecSetSizes(vv,n,N);CHKERRQ(ierr);
798     ierr = VecCreate_MPI_Private(vv,PETSC_TRUE,nghost,NULL);CHKERRQ(ierr);
799     w    = (Vec_MPI*)(vv)->data;
800     /* Create local representation */
801     ierr = VecGetArray(vv,&larray);CHKERRQ(ierr);
802     ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,n+nghost,larray,&w->localrep);CHKERRQ(ierr);
803     ierr = PetscLogObjectParent((PetscObject)vv,(PetscObject)w->localrep);CHKERRQ(ierr);
804     ierr = VecRestoreArray(vv,&larray);CHKERRQ(ierr);
805 
806     /*
807      Create scatter context for scattering (updating) ghost values
808      */
809     ierr = ISCreateGeneral(comm,nghost,ghosts,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
810     ierr = ISCreateStride(PETSC_COMM_SELF,nghost,n,1,&to);CHKERRQ(ierr);
811     ierr = VecScatterCreate(vv,from,w->localrep,to,&w->localupdate);CHKERRQ(ierr);
812     ierr = PetscLogObjectParent((PetscObject)vv,(PetscObject)w->localupdate);CHKERRQ(ierr);
813     ierr = ISDestroy(&to);CHKERRQ(ierr);
814     ierr = ISDestroy(&from);CHKERRQ(ierr);
815 
816     /* set local to global mapping for ghosted vector */
817     ierr = PetscMalloc1(n+nghost,&indices);CHKERRQ(ierr);
818     ierr = VecGetOwnershipRange(vv,&rstart,NULL);CHKERRQ(ierr);
819 
820     for (i=0; i<n; i++)      indices[i]   = rstart + i;
821     for (i=0; i<nghost; i++) indices[n+i] = ghosts[i];
822 
823     ierr = ISLocalToGlobalMappingCreate(comm,1,n+nghost,indices,PETSC_OWN_POINTER,&ltog);CHKERRQ(ierr);
824     ierr = VecSetLocalToGlobalMapping(vv,ltog);CHKERRQ(ierr);
825     ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
826   } else PetscCheckFalse(vv->ops->create == VecCreate_MPI,PetscObjectComm((PetscObject)vv),PETSC_ERR_ARG_WRONGSTATE,"Must set local or global size before setting ghosting");
827   else PetscCheckFalse(!((PetscObject)vv)->type_name,PetscObjectComm((PetscObject)vv),PETSC_ERR_ARG_WRONGSTATE,"Must set type to VECMPI before ghosting");
828   PetscFunctionReturn(0);
829 }
830 
831 /* ------------------------------------------------------------------------------------------*/
832 /*@C
833    VecCreateGhostBlockWithArray - Creates a parallel vector with ghost padding on each processor;
834    the caller allocates the array space. Indices in the ghost region are based on blocks.
835 
836    Collective
837 
838    Input Parameters:
839 +  comm - the MPI communicator to use
840 .  bs - block size
841 .  n - local vector length
842 .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
843 .  nghost - number of local ghost blocks
844 .  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)
845 -  array - the space to store the vector values (as long as n + nghost*bs)
846 
847    Output Parameter:
848 .  vv - the global vector representation (without ghost points as part of vector)
849 
850    Notes:
851    Use VecGhostGetLocalForm() to access the local, ghosted representation
852    of the vector.
853 
854    n is the local vector size (total local size not the number of blocks) while nghost
855    is the number of blocks in the ghost portion, i.e. the number of elements in the ghost
856    portion is bs*nghost
857 
858    Level: advanced
859 
860 .seealso: VecCreate(), VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
861           VecCreateGhost(), VecCreateSeqWithArray(), VecCreateMPIWithArray(),
862           VecCreateGhostWithArray(), VecCreateGhostBlock()
863 
864 @*/
865 PetscErrorCode  VecCreateGhostBlockWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],const PetscScalar array[],Vec *vv)
866 {
867   PetscErrorCode         ierr;
868   Vec_MPI                *w;
869   PetscScalar            *larray;
870   IS                     from,to;
871   ISLocalToGlobalMapping ltog;
872   PetscInt               rstart,i,nb,*indices;
873 
874   PetscFunctionBegin;
875   *vv = NULL;
876 
877   PetscCheckFalse(n == PETSC_DECIDE,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size");
878   PetscCheckFalse(nghost == PETSC_DECIDE,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local ghost size");
879   PetscCheckFalse(nghost < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Ghost length must be >= 0");
880   PetscCheckFalse(n % bs,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local size must be a multiple of block size");
881   ierr = PetscSplitOwnership(comm,&n,&N);CHKERRQ(ierr);
882   /* Create global representation */
883   ierr = VecCreate(comm,vv);CHKERRQ(ierr);
884   ierr = VecSetSizes(*vv,n,N);CHKERRQ(ierr);
885   ierr = VecSetBlockSize(*vv,bs);CHKERRQ(ierr);
886   ierr = VecCreate_MPI_Private(*vv,PETSC_TRUE,nghost*bs,array);CHKERRQ(ierr);
887   w    = (Vec_MPI*)(*vv)->data;
888   /* Create local representation */
889   ierr = VecGetArray(*vv,&larray);CHKERRQ(ierr);
890   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs,n+bs*nghost,larray,&w->localrep);CHKERRQ(ierr);
891   ierr = PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localrep);CHKERRQ(ierr);
892   ierr = VecRestoreArray(*vv,&larray);CHKERRQ(ierr);
893 
894   /*
895        Create scatter context for scattering (updating) ghost values
896   */
897   ierr = ISCreateBlock(comm,bs,nghost,ghosts,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
898   ierr = ISCreateStride(PETSC_COMM_SELF,bs*nghost,n,1,&to);CHKERRQ(ierr);
899   ierr = VecScatterCreate(*vv,from,w->localrep,to,&w->localupdate);CHKERRQ(ierr);
900   ierr = PetscLogObjectParent((PetscObject)*vv,(PetscObject)w->localupdate);CHKERRQ(ierr);
901   ierr = ISDestroy(&to);CHKERRQ(ierr);
902   ierr = ISDestroy(&from);CHKERRQ(ierr);
903 
904   /* set local to global mapping for ghosted vector */
905   nb     = n/bs;
906   ierr   = PetscMalloc1(nb+nghost,&indices);CHKERRQ(ierr);
907   ierr   = VecGetOwnershipRange(*vv,&rstart,NULL);CHKERRQ(ierr);
908   rstart = rstart/bs;
909 
910   for (i=0; i<nb; i++)      indices[i]    = rstart + i;
911   for (i=0; i<nghost; i++)  indices[nb+i] = ghosts[i];
912 
913   ierr = ISLocalToGlobalMappingCreate(comm,bs,nb+nghost,indices,PETSC_OWN_POINTER,&ltog);CHKERRQ(ierr);
914   ierr = VecSetLocalToGlobalMapping(*vv,ltog);CHKERRQ(ierr);
915   ierr = ISLocalToGlobalMappingDestroy(&ltog);CHKERRQ(ierr);
916   PetscFunctionReturn(0);
917 }
918 
919 /*@
920    VecCreateGhostBlock - Creates a parallel vector with ghost padding on each processor.
921         The indicing of the ghost points is done with blocks.
922 
923    Collective
924 
925    Input Parameters:
926 +  comm - the MPI communicator to use
927 .  bs - the block size
928 .  n - local vector length
929 .  N - global vector length (or PETSC_DECIDE to have calculated if n is given)
930 .  nghost - number of local ghost blocks
931 -  ghosts - global indices of ghost blocks, counts are by block, not by individual index, these do not need to be in increasing order (sorted)
932 
933    Output Parameter:
934 .  vv - the global vector representation (without ghost points as part of vector)
935 
936    Notes:
937    Use VecGhostGetLocalForm() to access the local, ghosted representation
938    of the vector.
939 
940    n is the local vector size (total local size not the number of blocks) while nghost
941    is the number of blocks in the ghost portion, i.e. the number of elements in the ghost
942    portion is bs*nghost
943 
944    Level: advanced
945 
946 .seealso: VecCreateSeq(), VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateMPI(),
947           VecGhostGetLocalForm(), VecGhostRestoreLocalForm(),
948           VecCreateGhostWithArray(), VecCreateMPIWithArray(), VecCreateGhostBlockWithArray()
949 
950 @*/
951 PetscErrorCode  VecCreateGhostBlock(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,PetscInt nghost,const PetscInt ghosts[],Vec *vv)
952 {
953   PetscErrorCode ierr;
954 
955   PetscFunctionBegin;
956   ierr = VecCreateGhostBlockWithArray(comm,bs,n,N,nghost,ghosts,NULL,vv);CHKERRQ(ierr);
957   PetscFunctionReturn(0);
958 }
959