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