xref: /petsc/doc/manual/vec.md (revision b59054469e8e75e7f1217f4456448261033ebd2d)
1(ch_vectors)=
2
3# Vectors and Parallel Data
4
5Vectors (denoted by `Vec`) are used to store discrete PDE solutions, right-hand sides for
6linear systems, etc. Users can create and manipulate entries in vectors directly with a basic, low-level interface or
7they can use the PETSc `DM` objects to connect actions on vectors to the type of discretization and grid that they are
8working with. These higher-level interfaces handle much of the details of the interactions with vectors and hence, are preferred
9in most situations. This chapter is organized as follows:
10
11- {any}`sec_veccreate`
12
13  - User managed
14  - {any}`sec_struct`
15  - {any}`sec_stag`
16  - {any}`sec_unstruct`
17  - {any}`sec_network`
18
19- Setting vector values
20
21  - For generic vectors
22  - {any}`sec_struct_set`
23  - {any}`sec_stag_set`
24  - {any}`sec_unstruct_set`
25  - {any}`sec_network_set`
26
27- {any}`sec_vecbasic`
28
29- {any}`sec_localglobal`
30
31- {any}`sec_scatter`
32
33  - {any}`sec_islocaltoglobalmap`
34  - {any}`sec_vecghost`
35
36- {any}`sec_ao`
37
38(sec_veccreate)=
39
40## Creating Vectors
41
42PETSc provides many ways to create vectors. The most basic, where the user is responsible for managing the
43parallel distribution of the vector entries, and a variety of higher-level approaches, based on `DM`, for classes of problems such
44as structured grids, staggered grids, unstructured grids, networks, and particles.
45
46The most basic way to create a vector with a local size of `m` and a global size of `M`, is to
47use
48
49```
50VecCreate(MPI_Comm comm, Vec *v);
51VecSetSizes(Vec v, PetscInt m, PetscInt M);
52VecSetFromOptions(Vec v);
53```
54
55which automatically generates the appropriate vector type (sequential or
56parallel) over all processes in `comm`. The option `-vec_type <type>`
57can be used in conjunction with
58`VecSetFromOptions()` to specify the use of a particular type of vector. For example, for NVIDIA GPU CUDA, use `cuda`.
59The GPU-based vectors allow
60one to set values on either the CPU or GPU but do their computations on the GPU.
61
62We emphasize that all processes in `comm` *must* call the vector
63creation routines since these routines are collective on all
64processes in the communicator. If you are unfamiliar with MPI
65communicators, see the discussion in {any}`sec_writing`. In addition, if a sequence of creation routines is
66used, they must be called in the same order for each process in the
67communicator.
68
69Instead of, or before calling `VecSetFromOptions()`, one can call
70
71```
72VecSetType(Vec v, VecType <VECCUDA, VECHIP, VECKOKKOS, etc.>)
73```
74
75One can create vectors whose entries are stored on GPUs using the convenience routine,
76
77```
78VecCreateMPICUDA(MPI_Comm comm, PetscInt m, PetscInt M, Vec *x);
79```
80
81There are convenience creation routines for almost all vector types; we recommend using the more verbose form because it allows
82selecting CPU or GPU simulations at runtime.
83
84For applications running in parallel that involve multi-dimensional structured grids, unstructured grids, networks, etc, it is cumbersome for users to explicitly manage the needed local and global sizes of the vectors.
85Hence, PETSc provides two powerful abstract objects (lower level) `PetscSection` (see {any}`ch_petscsection`) and (higher level) `DM` (see {any}`ch_dmbase`) to help manage the vectors and matrices needed for such applications. Using `DM`, parallel vectors can be created easily with
86
87```
88DMCreateGlobalVector(DM dm, Vec *v)
89```
90
91The `DM` object, see {any}`sec_struct`, {any}`sec_stag`, and {any}`ch_unstructured` for more details on `DM` for structured grids, staggered
92structured grids, and for unstructured grids,
93manages creating the correctly sized parallel vectors efficiently. One controls the type of vector that `DM` creates by calling
94
95```
96DMSetVecType(DM dm, VecType vt)
97```
98
99or by calling `DMSetFromOptions(DM dm)` and using the option `-dm_vec_type <standard or cuda or kokkos etc>`
100
101(sec_struct)=
102
103### DMDA - Creating vectors for structured grids
104
105Each `DM` type is suitable for a family of problems. The first of these, `DMDA`
106are intended for use with *logically structured rectangular grids*
107when communication of nonlocal data is needed before certain local
108computations can occur. `DMDA` is designed only for
109the case in which data can be thought of as being stored in a standard
110multidimensional array; thus, `DMDA` are *not* intended for
111parallelizing staggered arrays/grids, `DMSTAG` -- {any}`ch_stag`, or unstructured grid problems, `DMPLEX` -- {any}`ch_unstructured`, etc.
112
113For example, a typical situation one encounters in solving PDEs in
114parallel is that, to evaluate a local function, `f(x)`, each process
115requires its local portion of the vector `x` as well as its ghost
116points (the bordering portions of the vector that are owned by
117neighboring processes). Figure {any}`fig_ghosts` illustrates the
118ghost points for the seventh process of a two-dimensional, structured
119parallel grid. Each box represents a process; the ghost points for the
120seventh process’s local part of a parallel array are shown in gray.
121
122:::{figure} /images/manual/ghost.*
123:alt: Ghost Points for Two Stencil Types on the Seventh Process
124:name: fig_ghosts
125
126Ghost Points for Two Stencil Types on the Seventh Process
127:::
128
129The `DMDA` object
130contains parallel data layout information and communication
131information and is used to create vectors and matrices with
132the proper layout.
133
134One creates a `DMDA` two
135dimensions with the convenience routine
136
137```
138DMDACreate2d(MPI_Comm comm, DMBoundaryType xperiod, DMBoundaryType yperiod, DMDAStencilType st, PetscInt M, PetscInt N, PetscInt m, PetscInt n, PetscInt dof, PetscInt s, PetscInt *lx, PetscInt *ly, DM *da);
139```
140
141The arguments `M` and `N` indicate the global numbers of grid points
142in each direction, while `m` and `n` denote the process partition in
143each direction; `m*n` must equal the number of processes in the MPI
144communicator, `comm`. Instead of specifying the process layout, one
145may use `PETSC_DECIDE` for `m` and `n` so that PETSc will
146select the partition. The type of periodicity of the array
147is specified by `xperiod` and `yperiod`, which can be
148`DM_BOUNDARY_NONE` (no periodicity), `DM_BOUNDARY_PERIODIC`
149(periodic in that direction), `DM_BOUNDARY_TWIST` (periodic in that
150direction, but identified in reverse order), `DM_BOUNDARY_GHOSTED` ,
151or `DM_BOUNDARY_MIRROR`. The argument `dof` indicates the number of
152degrees of freedom at each array point, and `s` is the stencil width
153(i.e., the width of the ghost point region). The optional arrays `lx`
154and `ly` may contain the number of nodes along the x and y axis for
155each cell, i.e. the dimension of `lx` is `m` and the dimension of
156`ly` is `n`; alternately, `NULL` may be passed in.
157
158Two types of `DMDA` communication data structures can be
159created, as specified by `st`. Star-type stencils that radiate outward
160only in the coordinate directions are indicated by
161`DMDA_STENCIL_STAR`, while box-type stencils are specified by
162`DMDA_STENCIL_BOX`. For example, for the two-dimensional case,
163`DMDA_STENCIL_STAR` with width 1 corresponds to the standard 5-point
164stencil, while `DMDA_STENCIL_BOX` with width 1 denotes the standard
1659-point stencil. In both instances, the ghost points are identical, the
166only difference being that with star-type stencils, certain ghost points
167are ignored, substantially decreasing the number of messages sent. Note
168that the `DMDA_STENCIL_STAR` stencils can save interprocess
169communication in two and three dimensions.
170
171These `DMDA` stencils have nothing directly to do with a specific finite
172difference stencil one might choose to use for discretization; they
173only ensure that the correct values are in place for the application of a
174user-defined finite difference stencil (or any other discretization
175technique).
176
177The commands for creating `DMDA`
178in one and three dimensions are analogous:
179
180```
181DMDACreate1d(MPI_Comm comm, DMBoundaryType xperiod, PetscInt M, PetscInt w, PetscInt s, PetscInt *lc, DM *inra);
182```
183
184```
185DMDACreate3d(MPI_Comm comm, DMBoundaryType xperiod, DMBoundaryType yperiod, DMBoundaryType zperiod, DMDAStencilType stencil_type, PetscInt M, PetscInt N, PetscInt P, PetscInt m, PetscInt n, PetscInt p, PetscInt w, PetscInt s, PetscInt *lx, PetscInt *ly, PetscInt *lz, DM *inra);
186```
187
188The routines to create a `DM` are collective so that all
189processes in the communicator `comm` must call the same creation routines in the same order.
190
191A `DM` may be created, and its type set with
192
193```
194DMCreate(comm,&dm);
195DMSetType(dm,"Typename");  // for example, "DMDA"
196```
197
198Then `DMType` specific operations can be performed to provide information from which the specifics of the
199`DM` will be provided. For example,
200
201```
202DMSetDimension(dm, 1);
203DMDASetSizes(dm, M, 1, 1));
204DMDASetDof(dm, 1));
205DMSetUp(dm);
206```
207
208We now very briefly introduce a few more `DMType`.
209
210(sec_stag)=
211
212### DMSTAG - Creating vectors for staggered grids
213
214For structured grids with staggered data (living on elements, faces, edges,
215and/or vertices), the `DMSTAG` object is available. It behaves much
216like `DMDA`.
217See {any}`ch_stag` for discussion of creating vectors with `DMSTAG`.
218
219(sec_unstruct)=
220
221### DMPLEX - Creating vectors for unstructured grids
222
223See {any}`ch_unstructured` for a discussion of creating vectors with `DMPLEX`.
224
225(sec_network)=
226
227### DMNETWORK - Creating vectors for networks
228
229See {any}`ch_network` for discussion of creating vectors with `DMNETWORK`.
230
231## Common vector functions and operations
232
233One can examine (print out) a vector with the command
234
235```
236VecView(Vec x, PetscViewer v);
237```
238
239To print the vector to the screen, one can use the viewer
240`PETSC_VIEWER_STDOUT_WORLD`, which ensures that parallel vectors are
241printed correctly to `stdout`. To display the vector in an X-window,
242one can use the default X-windows viewer `PETSC_VIEWER_DRAW_WORLD`, or
243one can create a viewer with the routine `PetscViewerDrawOpen()`. A
244variety of viewers are discussed further in
245{any}`sec_viewers`.
246
247To create a new vector of the same format and parallel layout as an existing vector,
248use
249
250```
251VecDuplicate(Vec old, Vec *new);
252```
253
254To create several new vectors of the same format as an existing vector,
255use
256
257```
258VecDuplicateVecs(Vec old, PetscInt n, Vec **new);
259```
260
261This routine creates an array of pointers to vectors. The two routines
262are useful because they allow one to write library code that does
263not depend on the particular format of the vectors being used. Instead,
264the subroutines can automatically create work vectors based on
265the specified existing vector.
266
267When a vector is no longer needed, it should be destroyed with the
268command
269
270```
271VecDestroy(Vec *x);
272```
273
274To destroy an array of vectors, use the command
275
276```
277VecDestroyVecs(PetscInt n,Vec **vecs);
278```
279
280It is also possible to create vectors that use an array the user provides rather than having PETSc internally allocate the array space. Such
281vectors can be created with the routines such as
282
283```
284VecCreateSeqWithArray(PETSC_COMM_SELF, PetscInt bs, PetscInt n, PetscScalar *array, Vec *V);
285```
286
287```
288VecCreateMPIWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, PetscInt N, PetscScalar *array, Vec *V);
289```
290
291```
292VecCreateMPICUDAWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, PetscInt N, PetscScalar *array, Vec *V);
293```
294
295The `array` pointer should be a GPU memory location for GPU vectors.
296
297Note that here, one must provide the value `n`; it cannot be
298`PETSC_DECIDE` and the user is responsible for providing enough space
299in the array; `n*sizeof(PetscScalar)`.
300
301## Assembling (putting values in) vectors
302
303One can assign a single value to all components of a vector with
304
305```
306VecSet(Vec x, PetscScalar value);
307```
308
309Assigning values to individual vector components is more
310complicated to make it possible to write efficient parallel
311code. Assigning a set of components on a CPU is a two-step process: one first
312calls
313
314```
315VecSetValues(Vec x, PetscInt n, PetscInt *indices, PetscScalar *values, INSERT_VALUES);
316```
317
318any number of times on any or all of the processes. The argument `n`
319gives the number of components being set in this insertion. The integer
320array `indices` contains the *global component indices*, and
321`values` is the array of values to be inserted at those global component index locations. Any process can set
322any vector components; PETSc ensures that they are automatically
323stored in the correct location. Once all of the values have been
324inserted with `VecSetValues()`, one must call
325
326```
327VecAssemblyBegin(Vec x);
328```
329
330followed by
331
332```
333VecAssemblyEnd(Vec x);
334```
335
336to perform any needed message passing of nonlocal components. In order
337to allow the overlap of communication and calculation, the user’s code
338can perform any series of other actions between these two calls while
339the messages are in transition.
340
341Example usage of `VecSetValues()` may be found in
342<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/tutorials/ex2.c.html">src/vec/vec/tutorials/ex2.c</a>
343or <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/vec/vec/tutorials/ex2f.F90.html">src/vec/vec/tutorials/exf.F90</a>.
344
345Rather than inserting elements in a vector, one may wish to add
346values. This process is also done with the command
347
348```
349VecSetValues(Vec x, PetscInt n, PetscInt *indices, PetscScalar *values, ADD_VALUES);
350```
351
352Again, one must call the assembly routines `VecAssemblyBegin()` and
353`VecAssemblyEnd()` after all of the values have been added. Note that
354addition and insertion calls to `VecSetValues()` *cannot* be mixed.
355Instead, one must add and insert vector elements in phases, with
356intervening calls to the assembly routines. This phased assembly
357procedure overcomes the nondeterministic behavior that would occur if
358two different processes generated values for the same location, with one
359process adding while the other is inserting its value. (In this case, the
360addition and insertion actions could be performed in either order, thus
361resulting in different values at the particular location. Since PETSc
362does not allow the simultaneous use of `INSERT_VALUES` and
363`ADD_VALUES` this nondeterministic behavior will not occur in PETSc.)
364
365You can call `VecGetValues()` to pull local values from a vector (but
366not off-process values).
367
368For vectors obtained with `DMCreateGlobalVector()`, one can use `VecSetValuesLocal()` to set values into
369a global vector but using the local (ghosted) vector indexing of the vector entries. See also {any}`sec_islocaltoglobalmap`
370that allows one to provide arbitrary local-to-global mapping when not working with a `DM`.
371
372It is also possible to interact directly with the arrays that the vector values are stored
373in. The routine `VecGetArray()` returns a pointer to the elements local to
374the process:
375
376```
377VecGetArray(Vec v, PetscScalar **array);
378```
379
380When access to the array is no longer needed, the user should call
381
382```
383VecRestoreArray(Vec v, PetscScalar **array);
384```
385
386For vectors that may also have the array data in GPU memory, for example, `VECCUDA`, this call ensures the CPU array has the
387most recent array values by copying the data from the GPU memory if needed.
388
389If the values do not need to be modified, the routines
390
391```
392VecGetArrayRead(Vec v, const PetscScalar **array);
393VecRestoreArrayRead(Vec v, const PetscScalar **array);
394```
395
396should be used instead.
397
398:::{admonition} Listing: <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex1.c.html">SNES Tutorial src/snes/tutorials/ex1.c</a>
399```{literalinclude} /../src/snes/tutorials/ex1.c
400:end-at: PetscFunctionReturn(PETSC_SUCCESS);
401:name: snesex1
402:start-at: PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, PetscCtx ctx)
403```
404:::
405
406Minor differences exist in the Fortran interface for `VecGetArray()`
407and `VecRestoreArray()`, as discussed in
408{any}`sec_fortranarrays`. It is important to note that
409`VecGetArray()` and `VecRestoreArray()` do *not* copy the vector
410elements; they merely give users direct access to the vector elements.
411Thus, these routines require essentially no time to call and can be used
412efficiently.
413
414For GPU vectors, one can access either the values on the CPU as described above or one
415can call, for example,
416
417```
418VecCUDAGetArray(Vec v, PetscScalar **array);
419```
420
421:::{admonition} Listing: <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex47cu.cu.html">SNES Tutorial src/snes/tutorials/ex47cu.cu</a>
422```{literalinclude} /../src/snes/tutorials/ex47cu.cu
423:end-at: '}'
424:name: snesex47
425:start-at: PetscCall(VecCUDAGetArrayRead(xlocal, &xarray));
426```
427:::
428
429or
430
431```
432VecGetArrayAndMemType(Vec v, PetscScalar **array, PetscMemType *mtype);
433```
434
435which, in the first case, returns a GPU memory address and, in the second case, returns either a CPU or GPU memory
436address depending on the type of the vector. One can then launch a GPU kernel function that accesses the
437vector's memory for usage with GPUs. When computing on GPUs, `VecSetValues()` is not used! One always accesses the vector's arrays and passes them
438to the GPU code.
439
440It can also be convenient to treat the vector entries as a Kokkos view. One first creates Kokkos vectors and then calls
441
442```
443VecGetKokkosView(Vec v, Kokkos::View<const PetscScalar*, MemorySpace> *kv)
444```
445
446to set or access the vector entries.
447
448Of course, to provide the correct values to a vector, one must know what parts of the vector are owned by each MPI process.
449For parallel vectors, either CPU or GPU-based, it is possible to determine a process’s local range with the
450routine
451
452```
453VecGetOwnershipRange(Vec vec, PetscInt *start, PetscInt *end);
454```
455
456The argument `start` indicates the first component owned by the local
457process, while `end` specifies *one more than* the last owned by the
458local process. This command is useful, for instance, in assembling
459parallel vectors.
460
461If the `Vec` was obtained from a `DM` with `DMCreateGlobalVector()`, then the range values are determined by the specific `DM`.
462If the `Vec` was created directly, the range values are determined by the local size passed to `VecSetSizes()` or `VecCreateMPI()`.
463If `PETSC_DECIDE` was passed as the local size, then the vector uses default values for the range using `PetscSplitOwnership()`.
464For certain `DM`, such as `DMDA`, it is better to use `DM` specific routines, such as `DMDAGetGhostCorners()`, to determine
465the local values in the vector.
466
467Very occasionally, all MPI processes need to know all the range values, these can be obtained with
468
469```
470VecGetOwnershipRanges(Vec vec, PetscInt range[]);
471```
472
473The number of elements stored locally can be accessed with
474
475```
476VecGetLocalSize(Vec v, PetscInt *size);
477```
478
479The global vector length can be determined by
480
481```
482VecGetSize(Vec v, PetscInt *size);
483```
484
485(sec_struct_set)=
486
487### DMDA - Setting vector values
488
489PETSc provides an easy way to set values into the `DMDA` vectors and
490access them using the natural grid indexing. This is done with the
491routines
492
493```
494DMDAVecGetArray(DM dm, Vec l, void *array);
495... use the array indexing it with 1, 2, or 3 dimensions ...
496... depending on the dimension of the DMDA ...
497DMDAVecRestoreArray(DM dm, Vec l, void *array);
498DMDAVecGetArrayRead(DM dm, Vec l, void *array);
499... use the array indexing it with 1, 2, or 3 dimensions ...
500... depending on the dimension of the DMDA ...
501DMDAVecRestoreArrayRead(DM dm, Vec l, void *array);
502```
503
504where `array` is a multidimensional C array with the same dimension as `da`, and
505
506```
507DMDAVecGetArrayDOF(DM dm, Vec l, void *array);
508... use the array indexing it with 2, 3, or 4 dimensions ...
509... depending on the dimension of the DMDA ...
510DMDAVecRestoreArrayDOF(DM dm, Vec l, void *array);
511DMDAVecGetArrayDOFRead(DM dm, Vec l, void *array);
512... use the array indexing it with 2, 3, or 4 dimensions ...
513... depending on the dimension of the DMDA ...
514DMDAVecRestoreArrayDOFRead(DM dm, Vec l, void *array);
515```
516
517where `array` is a multidimensional C array with one more dimension than
518`da`. The vector `l` can be either a global vector or a local
519vector. The `array` is accessed using the usual *global* indexing on
520the entire grid, but the user may *only* refer to this array's local and ghost
521entries as all other entries are undefined. For example,
522for a scalar problem in two dimensions, one could use
523
524```
525PetscScalar **f, **u;
526...
527DMDAVecGetArrayRead(DM dm, Vec local, &u);
528DMDAVecGetArray(DM dm, Vec global, &f);
529...
530  f[i][j] = u[i][j] - ...
531...
532DMDAVecRestoreArrayRead(DM dm, Vec local, &u);
533DMDAVecRestoreArray(DM dm, Vec global, &f);
534```
535
536:::{admonition} Listing: <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex3.c.html">SNES Tutorial src/snes/tutorials/ex3.c</a>
537```{literalinclude} /../src/snes/tutorials/ex3.c
538:end-at: PetscFunctionReturn(PETSC_SUCCESS);
539:name: snesex3
540:start-at: PetscErrorCode FormFunction(SNES snes, Vec x, Vec f, PetscCtx ctx)
541```
542:::
543
544The recommended approach for multi-component PDEs is to declare a
545`struct` representing the fields defined at each node of the grid,
546e.g.
547
548```
549typedef struct {
550  PetscScalar u, v, omega, temperature;
551} Node;
552```
553
554and write the residual evaluation using
555
556```
557Node **f, **u;
558DMDAVecGetArray(DM dm, Vec local, &u);
559DMDAVecGetArray(DM dm, Vec global, &f);
560 ...
561    f[i][j].omega = ...
562 ...
563DMDAVecRestoreArray(DM dm, Vec local, &u);
564DMDAVecRestoreArray(DM dm, Vec global, &f);
565```
566
567The `DMDAVecGetArray` routines are also provided for GPU access with CUDA, HIP, and Kokkos. For example,
568
569```
570DMDAVecGetKokkosOffsetView(DM dm, Vec vec, Kokkos::View<const PetscScalar*XX*, MemorySpace> *ov)
571```
572
573where `*XX*` can contain any number of `*`. This allows one to write very natural Kokkos multi-dimensional parallel for kernels
574that act on the local portion of `DMDA` vectors.
575
576:::{admonition} Listing: <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex3k.kokkos.cxx.html">SNES Tutorial src/snes/tutorials/ex3k.kokkos.cxx</a>
577:name: snes-ex3-kokkos
578
579```{literalinclude} /../src/snes/tutorials/ex3k.kokkos.cxx
580:end-at: PetscFunctionReturn(PETSC_SUCCESS);
581:start-at: PetscErrorCode KokkosFunction(SNES snes, Vec x, Vec r, PetscCtx ctx)
582```
583:::
584
585The global indices of the lower left corner of the local portion of vectors obtained from `DMDA`
586as well as the local array size can be obtained with the commands
587
588```
589DMDAGetCorners(DM dm, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p);
590DMDAGetGhostCorners(DM dm, PetscInt *x, PetscInt *y, PetscInt *z, PetscInt *m, PetscInt *n, PetscInt *p);
591```
592
593These values can then be used as loop bounds for local function evaluations as demonstrated in the function examples above.
594
595The first version excludes ghost points, while the second
596includes them. The routine `DMDAGetGhostCorners()` deals with the fact
597that subarrays along boundaries of the problem domain have ghost points
598only on their interior edges, but not on their boundary edges.
599
600When either type of stencil is used, `DMDA_STENCIL_STAR` or
601`DMDA_STENCIL_BOX`, the local vectors (with the ghost points)
602represent rectangular arrays, including the extra corner elements in the
603`DMDA_STENCIL_STAR` case. This configuration provides simple access to
604the elements by employing two- (or three--) dimensional indexing. The
605only difference between the two cases is that when `DMDA_STENCIL_STAR`
606is used, the extra corner components are *not* scattered between the
607processes and thus contain undefined values that should *not* be used.
608
609(sec_stag_set)=
610
611### DMSTAG - Setting vector values
612
613For structured grids with staggered data (living on elements, faces, edges,
614and/or vertices), the `DMStag` object is available. It behaves
615like `DMDA`; see the `DMSTAG` manual page for more information.
616
617:::{admonition} Listing: <a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/dm/impls/stag/tutorials/ex6.c.html">SNES Tutorial src/dm/impls/stag/tutorials/ex6.c</a>
618```{literalinclude} /../src/dm/impls/stag/tutorials/ex6.c
619:end-at: /* Update x-velocity */
620:name: stagex6
621:start-at: static PetscErrorCode UpdateVelocity_2d(const Ctx *ctx, Vec velocity, Vec
622:  stress, Vec buoyancy)
623```
624:::
625
626(sec_unstruct_set)=
627
628### DMPLEX - Setting vector values
629
630See {any}`ch_unstructured` for a discussion on setting vector values with `DMPLEX`.
631
632(sec_network_set)=
633
634### DMNETWORK - Setting vector values
635
636See {any}`ch_network` for a discussion on setting vector values with `DMNETWORK`.
637
638(sec_vecbasic)=
639
640## Basic Vector Operations
641
642% Should make the table more attractive by using, for example, cloud_sptheme.ext.table_styling and the lines below
643% :column-alignment: left left
644% :widths: 72 28
645
646:::{container}
647:name: fig_vectorops
648
649```{eval-rst}
650.. table:: PETSc Vector Operations
651
652   +--------------------------------------------------------------+-----------------------------------+
653   | **Function Name**                                            | **Operation**                     |
654   +==============================================================+===================================+
655   | ``VecAXPY(Vec y, PetscScalar a, Vec x);``                    | :math:`y = y + a*x`               |
656   +--------------------------------------------------------------+-----------------------------------+
657   | ``VecAYPX(Vec y, PetscScalar a, Vec x);``                    | :math:`y = x + a*y`               |
658   +--------------------------------------------------------------+-----------------------------------+
659   | ``VecWAXPY(Vec  w, PetscScalar a, Vec x, Vec y);``           | :math:`w = a*x + y`               |
660   +--------------------------------------------------------------+-----------------------------------+
661   | ``VecAXPBY(Vec y, PetscScalar a, PetscScalar b, Vec x);``    | :math:`y = a*x + b*y`             |
662   +--------------------------------------------------------------+-----------------------------------+
663   | ``VecAXPBYPCZ(Vec z, PetscScalar a, PetscScalar b,           | :math:`z = a*x + b*y + c*z`       |
664   | PetscScalar c, Vec x, Vec y);``                              |                                   |
665   +--------------------------------------------------------------+-----------------------------------+
666   | ``VecScale(Vec x, PetscScalar a);``                          | :math:`x = a*x`                   |
667   +--------------------------------------------------------------+-----------------------------------+
668   | ``VecDot(Vec x, Vec y, PetscScalar *r);``                    | :math:`r = \bar{x}^T*y`           |
669   +--------------------------------------------------------------+-----------------------------------+
670   | ``VecTDot(Vec x, Vec y, PetscScalar *r);``                   | :math:`r = x'*y`                  |
671   +--------------------------------------------------------------+-----------------------------------+
672   | ``VecNorm(Vec x, NormType type, PetscReal *r);``             | :math:`r = ||x||_{type}`          |
673   +--------------------------------------------------------------+-----------------------------------+
674   | ``VecSum(Vec x, PetscScalar *r);``                           | :math:`r = \sum x_{i}`            |
675   +--------------------------------------------------------------+-----------------------------------+
676   | ``VecCopy(Vec x, Vec y);``                                   | :math:`y = x`                     |
677   +--------------------------------------------------------------+-----------------------------------+
678   | ``VecSwap(Vec x, Vec y);``                                   | :math:`y = x` while               |
679   |                                                              | :math:`x = y`                     |
680   +--------------------------------------------------------------+-----------------------------------+
681   | ``VecPointwiseMult(Vec w, Vec x, Vec y);``                   | :math:`w_{i} = x_{i}*y_{i}`       |
682   +--------------------------------------------------------------+-----------------------------------+
683   | ``VecPointwiseDivide(Vec w, Vec x, Vec y);``                 | :math:`w_{i} = x_{i}/y_{i}`       |
684   +--------------------------------------------------------------+-----------------------------------+
685   | ``VecMDot(Vec x, PetscInt n, Vec y[], PetscScalar *r);``     | :math:`r[i] = \bar{x}^T*y[i]`     |
686   +--------------------------------------------------------------+-----------------------------------+
687   | ``VecMTDot(Vec x, PetscInt n, Vec y[], PetscScalar *r);``    | :math:`r[i] = x^T*y[i]`           |
688   +--------------------------------------------------------------+-----------------------------------+
689   | ``VecMAXPY(Vec y, PetscInt n, PetscScalar *a, Vec x[]);``    | :math:`y = y + \sum_i a_{i}*x[i]` |
690   +--------------------------------------------------------------+-----------------------------------+
691   | ``VecMax(Vec x, PetscInt *idx, PetscReal *r);``              | :math:`r = \max x_{i}`            |
692   +--------------------------------------------------------------+-----------------------------------+
693   | ``VecMin(Vec x, PetscInt *idx, PetscReal *r);``              | :math:`r = \min x_{i}`            |
694   +--------------------------------------------------------------+-----------------------------------+
695   | ``VecAbs(Vec x);``                                           | :math:`x_i = |x_{i}|`             |
696   +--------------------------------------------------------------+-----------------------------------+
697   | ``VecReciprocal(Vec x);``                                    | :math:`x_i = 1/x_{i}`             |
698   +--------------------------------------------------------------+-----------------------------------+
699   | ``VecShift(Vec x, PetscScalar s);``                          | :math:`x_i = s + x_{i}`           |
700   +--------------------------------------------------------------+-----------------------------------+
701   | ``VecSet(Vec x, PetscScalar alpha);``                        | :math:`x_i = \alpha`              |
702   +--------------------------------------------------------------+-----------------------------------+
703```
704:::
705
706As the table lists, we have chosen certain
707basic vector operations to support within the PETSc vector library.
708These operations were selected because they often arise in application
709codes. The `NormType` argument to `VecNorm()` is one of `NORM_1`,
710`NORM_2`, or `NORM_INFINITY`. The 1-norm is $\sum_i |x_{i}|$,
711the 2-norm is $( \sum_{i} x_{i}^{2})^{1/2}$ and the infinity norm
712is $\max_{i} |x_{i}|$.
713
714In addition to `VecDot()` and `VecMDot()` and `VecNorm()`, PETSc
715provides split phase versions of this functionality that allow several independent
716inner products and/or norms to share the same communication (thus
717improving parallel efficiency). For example, one may have code such as
718
719```
720VecDot(Vec x, Vec y, PetscScalar *dot);
721VecMDot(Vec x, PetscInt nv, Vec y[], PetscScalar *dot);
722VecNorm(Vec x, NormType NORM_2, PetscReal *norm2);
723VecNorm(Vec x, NormType NORM_1, PetscReal *norm1);
724```
725
726This code works fine, but it performs four separate parallel
727communication operations. Instead, one can write
728
729```
730VecDotBegin(Vec x, Vec y, PetscScalar *dot);
731VecMDotBegin(Vec x, PetscInt nv, Vec y[], PetscScalar *dot);
732VecNormBegin(Vec x, NormType NORM_2, PetscReal *norm2);
733VecNormBegin(Vec x, NormType NORM_1, PetscReal *norm1);
734VecDotEnd(Vec x, Vec y, PetscScalar *dot);
735VecMDotEnd(Vec x,  PetscInt nv, Vec y[], PetscScalar *dot);
736VecNormEnd(Vec x, NormType NORM_2, PetscReal *norm2);
737VecNormEnd(Vec x, NormType NORM_1, PetscReal *norm1);
738```
739
740With this code, the communication is delayed until the first call to
741`VecxxxEnd()` at which a single MPI reduction is used to communicate
742all the values. It is required that the calls to the
743`VecxxxEnd()` are performed in the same order as the calls to the
744`VecxxxBegin()`; however, if you mistakenly make the calls in the
745wrong order, PETSc will generate an error informing you of this. There
746are additional routines `VecTDotBegin()` and `VecTDotEnd()`,
747`VecMTDotBegin()`, `VecMTDotEnd()`.
748
749For GPU vectors (like CUDA), the numerical computations will, by default, run on the GPU. Any
750scalar output, like the result of a `VecDot()` are placed in CPU memory.
751
752(sec_localglobal)=
753
754## Local/global vectors and communicating between vectors
755
756Many PDE problems require ghost (or halo) values in each MPI process or even more general parallel communication
757of vector values. These values are needed
758to perform function evaluation on that MPI process. The exact structure of the ghost values needed
759depends on the type of grid being used. `DM` provides a uniform API for communicating the needed
760values. We introduce the concept in detail for `DMDA`.
761
762Each `DM` object defines the layout of two vectors: a distributed
763global vector and a local vector that includes room for the appropriate
764ghost points. The `DM` object provides information about the size
765and layout of these vectors. The user can create
766vector objects that use the `DM` layout information with the
767routines
768
769```
770DMCreateGlobalVector(DM dm, Vec *g);
771DMCreateLocalVector(DM dm, Vec *l);
772```
773
774These vectors will generally serve as the building blocks for local and
775global PDE solutions, etc. If additional vectors with such layout
776information are needed in a code, they can be obtained by duplicating
777`l` or `g` via `VecDuplicate()` or `VecDuplicateVecs()`.
778
779We emphasize that a `DM` provides the information needed to
780communicate the ghost value information between processes. In most
781cases, several different vectors can share the same communication
782information (or, in other words, can share a given `DM`). The design
783of the `DM` object makes this easy, as each `DM` operation may
784operate on vectors of the appropriate size, as obtained via
785`DMCreateLocalVector()` and `DMCreateGlobalVector()` or as produced
786by `VecDuplicate()`.
787
788At certain stages of many applications, there is a need to work on a
789local portion of the vector that includes the ghost points. This may be
790done by scattering a global vector into its local parts by using the
791two-stage commands
792
793```
794DMGlobalToLocalBegin(DM dm, Vec g, InsertMode iora, Vec l);
795DMGlobalToLocalEnd(DM dm, Vec g, InsertMode iora, Vec l);
796```
797
798which allows the overlap of communication and computation. Since the
799global and local vectors, given by `g` and `l`, respectively, must
800be compatible with the `DM`, `da`, they should be
801generated by `DMCreateGlobalVector()` and `DMCreateLocalVector()`
802(or be duplicates of such a vector obtained via `VecDuplicate()`). The
803`InsertMode` can be `ADD_VALUES` or `INSERT_VALUES` among other possible values.
804
805One can scatter the local vectors into the distributed global vector with the
806command
807
808```
809DMLocalToGlobal(DM dm, Vec l, InsertMode mode, Vec g);
810```
811
812or the commands
813
814```
815DMLocalToGlobalBegin(DM dm, Vec l, InsertMode mode, Vec g);
816/* (Computation to overlap with communication) */
817DMLocalToGlobalEnd(DM dm, Vec l, InsertMode mode, Vec g);
818```
819
820In general this is used with an `InsertMode` of `ADD_VALUES`,
821because if one wishes to insert values into the global vector, they
822should access the global vector directly and put in the values.
823
824A third type of `DM` scatter is from a local vector
825(including ghost points that contain irrelevant values) to a local
826vector with correct ghost point values. This scatter may be done with
827the commands
828
829```
830DMLocalToLocalBegin(DM dm, Vec l1, InsertMode iora, Vec l2);
831DMLocalToLocalEnd(DM dm, Vec l1, InsertMode iora, Vec l2);
832```
833
834Since both local vectors, `l1` and `l2`, must be compatible with `da`, they should be generated by
835`DMCreateLocalVector()` (or be duplicates of such vectors obtained via
836`VecDuplicate()`). The `InsertMode` can be either `ADD_VALUES` or
837`INSERT_VALUES`.
838
839In most applications, the local ghosted vectors are only needed temporarily during
840user “function evaluations”. PETSc provides an easy, light-weight
841(requiring essentially no CPU time) way to temporarily obtain these work vectors and
842return them when no longer needed. This is done with the
843routines
844
845```
846DMGetLocalVector(DM dm, Vec *l);
847... use the local vector l ...
848DMRestoreLocalVector(DM dm, Vec *l);
849```
850
851(sec_scatter)=
852
853## Low-level Vector Communication
854
855Most users of PETSc who can utilize a `DM` will not need to utilize the lower-level routines discussed in the rest of this section
856and should skip ahead to {any}`ch_matrices`.
857
858To facilitate creating general vector scatters and gathers used, for example, in
859updating ghost points for problems for which no `DM` currently exists
860PETSc employs the concept of an *index set*, via the `IS` class. An
861index set, a generalization of a set of integer indices, is
862used to define scatters, gathers, and similar operations on vectors and
863matrices. Much of the underlying code that implements `DMGlobalToLocal` communication is built
864on the infrastructure discussed below.
865
866The following command creates an index set based on a list of integers:
867
868```
869ISCreateGeneral(MPI_Comm comm, PetscInt n, PetscInt *indices, PetscCopyMode mode, IS *is);
870```
871
872When `mode` is `PETSC_COPY_VALUES`, this routine copies the `n`
873indices passed to it by the integer array `indices`. Thus, the user
874should be sure to free the integer array `indices` when it is no
875longer needed, perhaps directly after the call to `ISCreateGeneral()`.
876The communicator, `comm`, should include all processes
877using the `IS`.
878
879Another standard index set is defined by a starting point (`first`)
880and a stride (`step`), and can be created with the command
881
882```
883ISCreateStride(MPI_Comm comm, PetscInt n, PetscInt first, PetscInt step, IS *is);
884```
885
886The meaning of `n`, `first`, and `step` correspond to the MATLAB notation
887`first:step:first+n*step`.
888
889Index sets can be destroyed with the command
890
891```
892ISDestroy(IS &is);
893```
894
895On rare occasions, the user may need to access information directly from
896an index set. Several commands assist in this process:
897
898```
899ISGetSize(IS is, PetscInt *size);
900ISStrideGetInfo(IS is, PetscInt *first, PetscInt *stride);
901ISGetIndices(IS is, PetscInt **indices);
902```
903
904The function `ISGetIndices()` returns a pointer to a list of the
905indices in the index set. For certain index sets, this may be a
906temporary array of indices created specifically for the call.
907Thus, once the user finishes using the array of indices, the routine
908
909```
910ISRestoreIndices(IS is, PetscInt **indices);
911```
912
913should be called to ensure that the system can free the space it may
914have used to generate the list of indices.
915
916A blocked version of index sets can be created with the command
917
918```
919ISCreateBlock(MPI_Comm comm, PetscInt bs, PetscInt n, PetscInt *indices, PetscCopyMode mode, IS *is);
920```
921
922This version is used for defining operations in which each element of
923the index set refers to a block of `bs` vector entries. Related
924routines analogous to those described above exist as well, including
925`ISBlockGetIndices()`, `ISBlockGetSize()`,
926`ISBlockGetLocalSize()`, `ISGetBlockSize()`.
927
928Most PETSc applications use a particular `DM` object to manage the communication details needed for their grids.
929In some rare cases, however, codes need to directly set up their required communication patterns. This is done using
930PETSc's `VecScatter` and `PetscSF` (for more general data than vectors). One
931can select any subset of the components of a vector to insert or add to
932any subset of the components of another vector. We refer to these
933operations as *generalized scatters*, though they are a
934combination of scatters and gathers.
935
936To copy selected components from one vector to another, one uses the
937following set of commands:
938
939```
940VecScatterCreate(Vec x, IS ix, Vec y, IS iy, VecScatter *ctx);
941VecScatterBegin(VecScatter ctx, Vec x, Vec y, INSERT_VALUES, SCATTER_FORWARD);
942VecScatterEnd(VecScatter ctx, Vec x, Vec y, INSERT_VALUES, SCATTER_FORWARD);
943VecScatterDestroy(VecScatter *ctx);
944```
945
946Here `ix` denotes the index set of the first vector, while `iy`
947indicates the index set of the destination vector. The vectors can be
948parallel or sequential. The only requirements are that the number of
949entries in the index set of the first vector, `ix`, equals the number
950in the destination index set, `iy`, and that the vectors be long
951enough to contain all the indices referred to in the index sets. If both
952`x` and `y` are parallel, their communicator must have the same set
953of processes, but their process order can differ. The argument
954`INSERT_VALUES` specifies that the vector elements will be inserted
955into the specified locations of the destination vector, overwriting any
956existing values. To add the components, rather than insert them, the
957user should select the option `ADD_VALUES` instead of
958`INSERT_VALUES`. One can also use `MAX_VALUES` or `MIN_VALUES` to
959replace the destination with the maximal or minimal of its current value and
960the scattered values.
961
962To perform a conventional gather operation, the user makes the
963destination index set, `iy`, be a stride index set with a stride of
964one. Similarly, a conventional scatter can be done with an initial
965(sending) index set consisting of a stride. The scatter routines are
966collective operations (i.e. all processes that own a parallel vector
967*must* call the scatter routines). When scattering from a parallel
968vector to sequential vectors, each process has its own sequential vector
969that receives values from locations as indicated in its own index set.
970Similarly, in scattering from sequential vectors to a parallel vector,
971each process has its own sequential vector that contributes to
972the parallel vector.
973
974*Caution*: When `INSERT_VALUES` is used, if two different processes
975contribute different values to the same component in a parallel vector,
976either value may be inserted. When `ADD_VALUES` is used, the
977correct sum is added to the correct location.
978
979In some cases, one may wish to “undo” a scatter, that is, perform the
980scatter backward, switching the roles of the sender and receiver. This
981is done by using
982
983```
984VecScatterBegin(VecScatter ctx, Vec y, Vec x, INSERT_VALUES, SCATTER_REVERSE);
985VecScatterEnd(VecScatter ctx, Vec y, Vec x, INSERT_VALUES, SCATTER_REVERSE);
986```
987
988Note that the roles of the first two arguments to these routines must be
989swapped whenever the `SCATTER_REVERSE` option is used.
990
991Once a `VecScatter` object has been created, it may be used with any
992vectors that have the same parallel data layout. That is, one can
993call `VecScatterBegin()` and `VecScatterEnd()` with different
994vectors than used in the call to `VecScatterCreate()` as long as they
995have the same parallel layout (the number of elements on each process are
996the same). Usually, these “different” vectors would have been obtained
997via calls to `VecDuplicate()` from the original vectors used in the
998call to `VecScatterCreate()`.
999
1000`VecGetValues()` can only access
1001local values from the vector. To get off-process values, the user should
1002create a new vector where the components will be stored and then
1003perform the appropriate vector scatter. For example, if one desires to
1004obtain the values of the 100th and 200th entries of a parallel vector,
1005`p`, one could use a code such as that below. In this example, the
1006values of the 100th and 200th components are placed in the array values.
1007In this example, each process now has the 100th and 200th component, but
1008obviously, each process could gather any elements it needed, or none by
1009creating an index set with no entries.
1010
1011```
1012Vec         p, x;         /* initial vector, destination vector */
1013VecScatter  scatter;      /* scatter context */
1014IS          from, to;     /* index sets that define the scatter */
1015PetscScalar *values;
1016PetscInt    idx_from[] = {100, 200}, idx_to[] = {0, 1};
1017
1018VecCreateSeq(PETSC_COMM_SELF, 2, &x);
1019ISCreateGeneral(PETSC_COMM_SELF, 2, idx_from, PETSC_COPY_VALUES, &from);
1020ISCreateGeneral(PETSC_COMM_SELF, 2, idx_to, PETSC_COPY_VALUES, &to);
1021VecScatterCreate(p, from, x, to, &scatter);
1022VecScatterBegin(scatter, p, x, INSERT_VALUES, SCATTER_FORWARD);
1023VecScatterEnd(scatter, p, x, INSERT_VALUES, SCATTER_FORWARD);
1024VecGetArray(x, &values);
1025ISDestroy(&from);
1026ISDestroy(&to);
1027VecScatterDestroy(&scatter);
1028```
1029
1030The scatter comprises two stages to allow for the overlap of
1031communication and computation. The introduction of the `VecScatter`
1032context allows the communication patterns for the scatter to be computed
1033once and reused repeatedly. Generally, even setting up the
1034communication for a scatter requires communication; hence, it is best to
1035reuse such information when possible.
1036
1037Scatters provide a very general method for managing the
1038communication of required ghost values for unstructured grid
1039computations. One scatters the global vector into a local “ghosted” work
1040vector, performs the computation on the local work vectors, and then
1041scatters back into the global solution vector. In the simplest case, this
1042may be written as
1043
1044```
1045VecScatterBegin(VecScatter scatter, Vec globalin, Vec localin, InsertMode INSERT_VALUES, ScatterMode SCATTER_FORWARD);
1046VecScatterEnd(VecScatter scatter, Vec globalin, Vec localin, InsertMode INSERT_VALUES, ScatterMode SCATTER_FORWARD);
1047/* For example, do local calculations from localin to localout */
1048...
1049VecScatterBegin(VecScatter scatter, Vec localout, Vec globalout, InsertMode ADD_VALUES, ScatterMode SCATTER_REVERSE);
1050VecScatterEnd(VecScatter scatter, Vec localout, Vec globalout, InsertMode ADD_VALUES, ScatterMode SCATTER_REVERSE);
1051```
1052
1053In this case, the scatter is used in a way similar to the usage of `DMGlobalToLocal()` and `DMLocalToGlobal()` discussed above.
1054
1055(sec_islocaltoglobalmap)=
1056
1057### Local to global mappings
1058
1059When working with a global representation of a vector
1060(usually on a vector obtained with `DMCreateGlobalVector()`) and a local
1061representation of the same vector that includes ghost points required
1062for local computation (obtained with `DMCreateLocalVector()`). PETSc provides routines to help map indices from
1063a local numbering scheme to the PETSc global numbering scheme, recall their use above for the routine `VecSetValuesLocal()` introduced above.
1064This is done via the following routines
1065
1066```
1067ISLocalToGlobalMappingCreate(MPI_Comm comm, PetscInt bs, PetscInt N, PetscInt* globalnum, PetscCopyMode mode, ISLocalToGlobalMapping* ctx);
1068ISLocalToGlobalMappingApply(ISLocalToGlobalMapping ctx, PetscInt n, PetscInt *in, PetscInt *out);
1069ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping ctx, IS isin, IS* isout);
1070ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *ctx);
1071```
1072
1073Here `N` denotes the number of local indices, `globalnum` contains
1074the global number of each local number, and `ISLocalToGlobalMapping`
1075is the resulting PETSc object that contains the information needed to
1076apply the mapping with either `ISLocalToGlobalMappingApply()` or
1077`ISLocalToGlobalMappingApplyIS()`.
1078
1079Note that the `ISLocalToGlobalMapping` routines serve a different
1080purpose than the `AO` routines. In the former case, they provide a
1081mapping from a local numbering scheme (including ghost points) to a
1082global numbering scheme, while in the latter, they provide a mapping
1083between two global numbering schemes. Many applications may use
1084both `AO` and `ISLocalToGlobalMapping` routines. The `AO` routines
1085are first used to map from an application global ordering (that has no
1086relationship to parallel processing, etc.) to the PETSc ordering scheme
1087(where each process has a contiguous set of indices in the numbering).
1088Then, to perform function or Jacobian evaluations locally on
1089each process, one works with a local numbering scheme that includes
1090ghost points. The mapping from this local numbering scheme back to the
1091global PETSc numbering can be handled with the
1092`ISLocalToGlobalMapping` routines.
1093
1094If one is given a list of block indices in a global numbering, the
1095routine
1096
1097```
1098ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping ctx, ISGlobalToLocalMappingMode type, PetscInt nin, PetscInt idxin[], PetscInt *nout, PetscInt idxout[]);
1099```
1100
1101will provide a new list of indices in the local numbering. Again,
1102negative values in `idxin` are left unmapped. But in addition, if
1103`type` is set to `IS_GTOLM_MASK` , then `nout` is set to `nin`
1104and all global values in `idxin` that are not represented in the local
1105to global mapping are replaced by -1. When `type` is set to
1106`IS_GTOLM_DROP`, the values in `idxin` that are not represented
1107locally in the mapping are not included in `idxout`, so that
1108potentially `nout` is smaller than `nin`. One must pass in an array
1109long enough to hold all the indices. One can call
1110`ISGlobalToLocalMappingApplyBlock()` with `idxout` equal to `NULL`
1111to determine the required length (returned in `nout`) and then
1112allocate the required space and call
1113`ISGlobalToLocalMappingApplyBlock()` a second time to set the values.
1114
1115Often it is convenient to set elements into a vector using the local
1116node numbering rather than the global node numbering (e.g., each process
1117may maintain its own sublist of vertices and elements and number them
1118locally). To set values into a vector with the local numbering, one must
1119first call
1120
1121```
1122VecSetLocalToGlobalMapping(Vec v, ISLocalToGlobalMapping ctx);
1123```
1124
1125and then call
1126
1127```
1128VecSetValuesLocal(Vec x, PetscInt n, const PetscInt indices[], const PetscScalar values[], INSERT_VALUES);
1129```
1130
1131Now the `indices` use the local numbering rather than the global,
1132meaning the entries lie in $[0,n)$ where $n$ is the local
1133size of the vector. Global vectors obtained from a `DM` already have the global to local mapping
1134provided by the `DM`.
1135
1136One can use global indices
1137with `MatSetValues()` or `MatSetValuesStencil()` to assemble global stiffness matrices. Alternately, the
1138global node number of each local node, including the ghost nodes, can be
1139obtained by calling
1140
1141```
1142DMGetLocalToGlobalMapping(DM dm, ISLocalToGlobalMapping *map);
1143```
1144
1145followed by
1146
1147```
1148VecSetLocalToGlobalMapping(Vec v, ISLocalToGlobalMapping map);
1149MatSetLocalToGlobalMapping(Mat A, ISLocalToGlobalMapping rmapping, ISLocalToGlobalMapping cmapping);
1150```
1151
1152Now, entries may be added to the vector and matrix using the local
1153numbering and `VecSetValuesLocal()` and `MatSetValuesLocal()`.
1154
1155The example
1156<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex5.c.html">SNES Tutorial ex5</a>
1157illustrates the use of a `DMDA` in the solution of a
1158nonlinear problem. The analogous Fortran program is
1159<a href="PETSC_DOC_OUT_ROOT_PLACEHOLDER/src/snes/tutorials/ex5f90.F90.html">SNES Tutorial ex5f90</a>;
1160see {any}`ch_snes` for a discussion of the
1161nonlinear solvers.
1162
1163(sec_vecghost)=
1164
1165### Global Vectors with locations for ghost values
1166
1167There are two minor drawbacks to the basic approach described above for unstructured grids:
1168
1169- the extra memory requirement for the local work vector, `localin`,
1170  which duplicates the local values in the memory in `globalin`, and
1171- the extra time required to copy the local values from `localin` to
1172  `globalin`.
1173
1174An alternative approach is to allocate global vectors with space
1175preallocated for the ghost values.
1176
1177```
1178VecCreateGhost(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, PetscInt *ghosts, Vec *vv)
1179```
1180
1181or
1182
1183```
1184VecCreateGhostWithArray(MPI_Comm comm, PetscInt n, PetscInt N, PetscInt nghost, PetscInt *ghosts, PetscScalar *array, Vec *vv)
1185```
1186
1187Here `n` is the number of local vector entries, `N` is the number of
1188global entries (or `NULL`), and `nghost` is the number of ghost
1189entries. The array `ghosts` is of size `nghost` and contains the
1190global vector location for each local ghost location. Using
1191`VecDuplicate()` or `VecDuplicateVecs()` on a ghosted vector will
1192generate additional ghosted vectors.
1193
1194In many ways, a ghosted vector behaves like any other MPI vector
1195created by `VecCreateMPI()`. The difference is that the ghosted vector
1196has an additional “local” representation that allows one to access the
1197ghost locations. This is done through the call to
1198
1199```
1200VecGhostGetLocalForm(Vec g,Vec *l);
1201```
1202
1203The vector `l` is a sequential representation of the parallel vector
1204`g` that shares the same array space (and hence numerical values); but
1205allows one to access the “ghost” values past “the end of the” array.
1206Note that one accesses the entries in `l` using the local numbering of
1207elements and ghosts, while they are accessed in `g` using the global
1208numbering.
1209
1210A common usage of a ghosted vector is given by
1211
1212```
1213VecGhostUpdateBegin(Vec globalin, InsertMode INSERT_VALUES, ScatterMode SCATTER_FORWARD);
1214VecGhostUpdateEnd(Vec globalin, InsertMode INSERT_VALUES, ScatterMode SCATTER_FORWARD);
1215VecGhostGetLocalForm(Vec globalin, Vec *localin);
1216VecGhostGetLocalForm(Vec globalout, Vec *localout);
1217...  Do local calculations from localin to localout ...
1218VecGhostRestoreLocalForm(Vec globalin, Vec *localin);
1219VecGhostRestoreLocalForm(Vec globalout, Vec *localout);
1220VecGhostUpdateBegin(Vec globalout, InsertMode ADD_VALUES, ScatterMode SCATTER_REVERSE);
1221VecGhostUpdateEnd(Vec globalout, InsertMode ADD_VALUES, ScatterMode SCATTER_REVERSE);
1222```
1223
1224The routines `VecGhostUpdateBegin()` and `VecGhostUpdateEnd()` are
1225equivalent to the routines `VecScatterBegin()` and `VecScatterEnd()`
1226above, except that since they are scattering into the ghost locations,
1227they do not need to copy the local vector values, which are already in
1228place. In addition, the user does not have to allocate the local work
1229vector since the ghosted vector already has allocated slots to contain
1230the ghost values.
1231
1232The input arguments `INSERT_VALUES` and `SCATTER_FORWARD` cause the
1233ghost values to be correctly updated from the appropriate process. The
1234arguments `ADD_VALUES` and `SCATTER_REVERSE` update the “local”
1235portions of the vector from all the other processes’ ghost values. This
1236would be appropriate, for example, when performing a finite element
1237assembly of a load vector. One can also use `MAX_VALUES` or
1238`MIN_VALUES` with `SCATTER_REVERSE`.
1239
1240`DMPLEX` does not yet support ghosted vectors sharing memory with the global representation.
1241This is a work in progress; if you are interested in this feature, please contact the PETSc community members.
1242
1243{any}`sec_partitioning` discusses the important topic of
1244partitioning an unstructured grid.
1245
1246(sec_ao)=
1247
1248## Application Orderings
1249
1250When writing parallel PDE codes, there is extra complexity caused by
1251having multiple ways of indexing (numbering) and ordering objects such
1252as vertices and degrees of freedom. For example, a grid generator or
1253partitioner may renumber the nodes, requiring adjustment of the other
1254data structures that refer to these objects; see Figure
1255{any}`fig_daao`.
1256PETSc provides various tools to help manage the mapping amongst
1257the various numbering systems. The most basic is the `AO`
1258(application ordering), which enables mapping between different global
1259(cross-process) numbering schemes.
1260
1261In many applications, it is desirable to work with one or more
1262“orderings” (or numberings) of degrees of freedom, cells, nodes, etc.
1263Doing so in a parallel environment is complicated by the fact that each
1264process cannot keep complete lists of the mappings between different
1265orderings. In addition, the orderings used in the PETSc linear algebra
1266routines (often contiguous ranges) may not correspond to the “natural”
1267orderings for the application.
1268
1269PETSc provides certain utility routines that allow one to deal cleanly
1270and efficiently with the various orderings. To define a new application
1271ordering (called an `AO` in PETSc), one can call the routine
1272
1273```
1274AOCreateBasic(MPI_Comm comm, PetscInt n, const PetscInt apordering[], const PetscInt petscordering[], AO *ao);
1275```
1276
1277The arrays `apordering` and `petscordering`, respectively, contain a
1278list of integers in the application ordering and their corresponding
1279mapped values in the PETSc ordering. Each process can provide whatever
1280subset of the ordering it chooses, but multiple processes should never
1281contribute duplicate values. The argument `n` indicates the number of
1282local contributed values.
1283
1284For example, consider a vector of length 5, where node 0 in the
1285application ordering corresponds to node 3 in the PETSc ordering. In
1286addition, nodes 1, 2, 3, and 4 of the application ordering correspond,
1287respectively, to nodes 2, 1, 4, and 0 of the PETSc ordering. We can
1288write this correspondence as
1289
1290$$
1291\{ 0, 1, 2, 3, 4 \}  \to  \{ 3, 2, 1, 4, 0 \}.
1292$$
1293
1294The user can create the PETSc `AO` mappings in several ways. For
1295example, if using two processes, one could call
1296
1297```
1298AOCreateBasic(PETSC_COMM_WORLD, 2,{0, 3}, {3, 4}, &ao);
1299```
1300
1301on the first process and
1302
1303```
1304AOCreateBasic(PETSC_COMM_WORLD, 3, {1, 2, 4}, {2, 1, 0}, &ao);
1305```
1306
1307on the other process.
1308
1309Once the application ordering has been created, it can be used with
1310either of the commands
1311
1312```
1313AOPetscToApplication(AO ao, PetscInt n, PetscInt *indices);
1314AOApplicationToPetsc(AO ao, PetscInt n, PetscInt *indices);
1315```
1316
1317Upon input, the `n`-dimensional array `indices` specifies the
1318indices to be mapped, while upon output, `indices` contains the mapped
1319values. Since we, in general, employ a parallel database for the `AO`
1320mappings, it is crucial that all processes that called
1321`AOCreateBasic()` also call these routines; these routines *cannot* be
1322called by just a subset of processes in the MPI communicator that was
1323used in the call to `AOCreateBasic()`.
1324
1325An alternative routine to create the application ordering, `AO`, is
1326
1327```
1328AOCreateBasicIS(IS apordering, IS petscordering, AO *ao);
1329```
1330
1331where index sets are used
1332instead of integer arrays.
1333
1334The mapping routines
1335
1336```
1337AOPetscToApplicationIS(AO ao, IS indices);
1338AOApplicationToPetscIS(AO ao, IS indices);
1339```
1340
1341will map index sets (`IS` objects) between orderings. Both the
1342`AOXxxToYyy()` and `AOXxxToYyyIS()` routines can be used regardless
1343of whether the `AO` was created with a `AOCreateBasic()` or
1344`AOCreateBasicIS()`.
1345
1346The `AO` context should be destroyed with `AODestroy(AO *ao)` and
1347viewed with `AOView(AO ao,PetscViewer viewer)`.
1348
1349Although we refer to the two orderings as “PETSc” and “application”
1350orderings, the user is free to use them both for application orderings
1351and to maintain relationships among a variety of orderings by employing
1352several `AO` contexts.
1353
1354The `AOxxToxx()` routines allow negative entries in the input integer
1355array. These entries are not mapped; they remain unchanged. This
1356functionality enables, for example, mapping neighbor lists that use
1357negative numbers to indicate nonexistent neighbors due to boundary
1358conditions, etc.
1359
1360Since the global ordering that PETSc uses to manage its parallel vectors
1361(and matrices) does not usually correspond to the “natural” ordering of
1362a two- or three-dimensional array, the `DMDA` structure provides an
1363application ordering `AO` (see {any}`sec_ao`) that maps
1364between the natural ordering on a rectangular grid and the ordering
1365PETSc uses to parallelize. This ordering context can be obtained with
1366the command
1367
1368```
1369DMDAGetAO(DM dm, AO *ao);
1370```
1371
1372In Figure {any}`fig_daao`, we indicate the orderings for a
1373two-dimensional `DMDA`, divided among four processes.
1374
1375:::{figure} /images/manual/danumbering.*
1376:alt: Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four Processes)
1377:name: fig_daao
1378
1379Natural Ordering and PETSc Ordering for a 2D Distributed Array (Four
1380Processes)
1381:::
1382
1383
1384(sec_petscsf)=
1385
1386# PetscSF - an alternative to low-level MPI calls for data communication
1387
1388As discussed above, the `VecScatter` object allows one to define parallel communication between vectors by listing, with `IS` objects, which vector entries from one vector
1389are to be communicated to another vector and where in the second vector they are to be inserted. `PetscSF` provides a similar more general functionality for arrays of any MPI datatype.
1390
1391`PetscSF` communicates between `rootdata` and `leafdata` arrays. `rootdata` is distributed across the MPI processes and its entries are indicated by a `PetscSFNode` pair consisting of the MPI `rank` the
1392entry is located on and the `index` in the array on that MPI process.
1393
1394```
1395typedef struct {
1396  PetscInt rank;  /* MPI rank of owner */
1397  PetscInt index; /* Index of node on rank */
1398} PetscSFNode;
1399```
1400
1401Each entry is uniquely owned at that location; in the same way a PETSc global vector has unique MPI process ownership of each entry.
1402
1403`leafdata` is similar to PETSc local vectors; each MPI process's `leafdata` array can contain "ghost values" that match values in other locations of the `leafdata` (on the same or different
1404MPI processes). All these matching ghost values share a common root value in `rootdata`.
1405
1406
1407We begin to explain the use of `PetscSF` with an example. First we construct an array that tells for each leaf entry on that MPI process where its root entry is:
1408
1409```
1410  PetscInt    nroots, nleaves;
1411  PetscSFNode *roots;
1412
1413  if (rank == 0) {
1414    // the number of entries in rootdata on this MPI process
1415    nroots = 2;
1416    // provide the matching location in rootdata for each entry in leafdata
1417    nleaves = 3;
1418    roots[0].rank = 0; roots[0].index = 0;
1419    roots[1].rank = 1; roots[1].index = 0;
1420    roots[2].rank = 0; roots[2].index = 1;
1421  } else {
1422    nroots = 1;
1423    nleaves = 3;
1424    roots[0].rank = 0; roots[0].index = 0;
1425    roots[1].rank = 0; roots[1].index = 1;
1426    roots[2].rank = 1; roots[1].index = 0;
1427  }
1428```
1429
1430Next, we construct the `PetscSF` that encapsulates this information needed for communication:
1431
1432```
1433  PetscSF sf;
1434
1435  PetscSFCreate(PETSC_COMM_WORLD, &sf);
1436  PetscSFSetFromOptions(sf);
1437  PetscSFSetGraph(sf, nroots, nleaves, NULL, PETSC_OWN_POINTER, roots, PETSC_OWN_POINTER);
1438  PetscSFSetUp(sf);
1439```
1440
1441Next we fill `rootdata`:
1442
1443```
1444  PetscInt    *rootdata, *leafdata;
1445
1446  if (rank == 0) {
1447    rootdata[0] = 1;
1448    rootdata[1] = 2;
1449  } else {
1450    rootdata[0] = 3;
1451  }
1452```
1453
1454Finally, we use the `PetscSF` to communicate `rootdata` to `leafdata`:
1455
1456```
1457  PetscSFBcastBegin(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE);
1458  PetscSFBcastEnd(sf, MPIU_INT, rootdata, leafdata, MPI_REPLACE);
1459```
1460
1461Now `leafdata` on MPI rank 0 contains (1, 3, 2) and on MPI rank 1 contains (1, 2, 3).
1462
1463It is also possible to move `leafdata` to `rootdata` using
1464
1465```
1466  PetscSFReduceBegin(sf, MPIU_INT, leafdata, rootdata, MPIU_SUM);
1467  PetscSFReduceEnd(sf, MPIU_INT, leafdata, rootdata, MPIU_SUM);
1468```
1469
1470In this case, since the reduction operation performed (the final argument of `PetscSFReduceBegin()`), is `MPIU_SUM` the final result in each entry of `rootdata` is the sum
1471of the previous value at that location plus all the values it that entries leafs. So `rootdata` on MPI rank 0 contains (3, 6) while on MPI rank
14721 it contains (9).
1473
1474As shown in the example above, `PetscSFBcastBegin()` and `PetscSFBcastEnd()` (as well as other `PetscSF` functions) also take an `MPI_Op` reduction argument, though that is almost always `MPI_REPLACE`.
1475
1476## Non-contiguous storage of leafdata
1477
1478In the example above we treated the `leafdata` as sitting in a contiguous array with entries from 0 to one less than `nleaves`. This is indicated by the
1479`NULL` argument in the call to `PetscSFSetGraph()`. More generally the `leafdata` array can have entries in it that are not accessed by the `PetscSF` operations. For example,
1480
1481```
1482 PetscInt *leaves;
1483
1484 if (rank == 0) {
1485   leaves[0] = 1;
1486   leaves[1] = 2;
1487   leaves[2] = 4;
1488 } else {
1489   leaves[0] = 2;
1490   leaves[1] = 1;
1491   leaves[2] = 0;
1492 }
1493  PetscSFSetGraph(sf, nroots, nleaves, leaves, PETSC_OWN_POINTER, roots, PETSC_OWN_POINTER);
1494```
1495
1496means that the three entries of `leafdata` affected by `PetscSF` communication on MPI rank 0 are the array locations (1, 2, 4); meaning also that `leafdata` must be of length at least 5.
1497On MPI rank 1, the arriving values from the three roots listed in `roots` are placed backwards in `leafdata`. Note that providing the `leaves` permutation array on MPI rank 1 is
1498equivalent to listing the three values in `roots` in the opposite order.
1499
1500If we reran the initial communication with `PetscSFBcastBegin()` and `PetscSFBcastEnd()` using the modified `sf` the resulting values in `leavedata` would be on MPI rank 0 (x, 1, 3, x, 2) and on
1501MPI rank 1 (3, 2, 1) where x indicates the previous value in `leafdata` that was unchanged.
1502
1503
1504
1505## GPU usage
1506
1507`rootdata` and `leafdata` can live either on CPU memory or GPU memory. The `PetscSF` routines automatically detect the memory type. But the time for the calls to the `CUDA` or `HIP`
1508routines for doing this determination (`cudaPointerGetAttributes()` or `hipPointerGetAttributes()`) is not trivial. To avoid the cost of the check,
1509`PetscSF` provides the routines `PetscSFBcastWithMemTypeBegin()` and `PetscSFReduceWithMemTypeBegin()` where the user provides the memory type information.
1510
1511## Gathering leafdata but not reducing it
1512
1513One may wish to gather the entries of the `leafdata` for each root but not reduce them to a single value. This is done with
1514
1515```
1516  PetscSFGatherBegin(sf, MPIU_INT, leafdata, multirootdata);
1517  PetscSFGatherEnd(sf, MPIU_INT, leafdata, multirootdata);
1518```
1519
1520Here `multirootdata` is (generally) an array larger than `rootdata` that has enough locations to store the value of each `leaf` of each local root. The values are stored contiguously
1521for each root; that is `multirootdata` will contain
1522
1523```
1524(first leaf of first root, second leaf of first root, third leaf of first root, ..., last leaf of first root, first leaf of second root, second leaf of second root, ...)
1525```
1526
1527The number of leaves for each local root (sometimes called the degree of the root) can be obtained with calls to `PetscSFComputeDegreeBegin()` and `PetscSFComputeDegreeEnd()`.
1528
1529The data in `multirootdata` can be communicated to `leafdata` using
1530
1531```
1532  PetscSFScatterBegin(sf, MPIU_INT, multirootdata, leafdata);
1533  PetscSFScatterEnd(sf, MPIU_INT, multirootdata, leafdata);
1534```
1535
1536## Optimized communication patterns
1537
1538
1539A performance drawback to using `PetscSFSetGraph()` is that it requires explicitly listing in arrays all the entries of `roots`.
1540`PetscSFSetGraphWithPattern()` provides an alternative way to indicate the communication graph for specific communication patterns.
1541
1542
1543
1544
1545
1546
1547
1548
1549