Actual source code: petscvec.h

  1: /*
  2:     Defines the vector component of PETSc. Vectors generally represent
  3:   degrees of freedom for finite element/finite difference functions
  4:   on a grid. They have more mathematical structure then simple arrays.
  5: */
  6: #pragma once

  8: #include <petscsys.h>
  9: #include <petscsftypes.h>
 10: #include <petscis.h>
 11: #include <petscdevicetypes.h>
 12: #include <petscviewer.h>

 14: /* SUBMANSEC = Vec */

 16: /*S
 17:    Vec - Abstract PETSc vector object. Used for holding solutions and right-hand sides for (non) linear systems and integrators

 19:    Level: beginner

 21:    Note:
 22:    Internally the actual vector representation is generally a simple array but most PETSc code can work on other representations through this abstraction

 24: .seealso: [](doc_vector), [](ch_vectors), `VecCreate()`, `VecType`, `VecSetType()`
 25: S*/
 26: typedef struct _p_Vec *Vec;

 28: /*E
 29:   ScatterMode - Determines the direction of a scatter

 31:   Values:
 32: +  `SCATTER_FORWARD`       - Scatters the values as dictated by the `VecScatterCreate()` call
 33: .  `SCATTER_REVERSE`       - Moves the values in the opposite direction than the directions indicated in the `VecScatterCreate()` call
 34: .  `SCATTER_FORWARD_LOCAL` - Scatters the values as dictated by the `VecScatterCreate()` call except NO MPI communication is done
 35: -  `SCATTER_REVERSE_LOCAL` - Moves the values in the opposite direction than the directions indicated in the `VecScatterCreate()` call
 36:                              except NO MPI communication is done

 38:   Level: beginner

 40: .seealso: [](ch_vectors), `VecScatter`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_FORWARD`, `SCATTER_REVERSE`, `SCATTER_FORWARD_LOCAL`, `SCATTER_REVERSE_LOCAL`
 41: E*/
 42: typedef enum {
 43:   SCATTER_FORWARD       = 0,
 44:   SCATTER_REVERSE       = 1,
 45:   SCATTER_FORWARD_LOCAL = 2,
 46:   SCATTER_REVERSE_LOCAL = 3
 47: } ScatterMode;

 49: /*MC
 50:     SCATTER_FORWARD - Scatters the values as dictated by the `VecScatterCreate()` call

 52:     Level: beginner

 54: .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_REVERSE`, `SCATTER_FORWARD_LOCAL`,
 55:           `SCATTER_REVERSE_LOCAL`
 56: M*/

 58: /*MC
 59:     SCATTER_REVERSE - Moves the values in the opposite direction then the directions indicated in
 60:                       in the `VecScatterCreate()`

 62:     Level: beginner

 64: .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_FORWARD`, `SCATTER_FORWARD_LOCAL`,
 65:           `SCATTER_REVERSE_LOCAL`
 66: M*/

 68: /*MC
 69:     SCATTER_FORWARD_LOCAL - Scatters the values as dictated by the `VecScatterCreate()` call except NO parallel communication
 70:                             is done. Any variables that have be moved between processes are ignored

 72:     Level: developer

 74: .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_REVERSE`, `SCATTER_FORWARD`,
 75:           `SCATTER_REVERSE_LOCAL`
 76: M*/

 78: /*MC
 79:     SCATTER_REVERSE_LOCAL - Moves the values in the opposite direction then the directions indicated in
 80:                            in the `VecScatterCreate()`  except NO parallel communication
 81:                            is done. Any variables that have be moved between processes are ignored

 83:     Level: developer

 85: .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_FORWARD`, `SCATTER_FORWARD_LOCAL`,
 86:           `SCATTER_REVERSE`
 87: M*/

 89: /*J
 90:    VecType - String with the name of a PETSc vector type

 92:    Level: beginner

 94: .seealso: [](doc_vector), [](ch_vectors), `VecSetType()`, `Vec`, `VecCreate()`, `VecDestroy()`
 95: J*/
 96: typedef const char *VecType;
 97: #define VECSEQ         "seq"
 98: #define VECMPI         "mpi"
 99: #define VECSTANDARD    "standard" /* seq on one process and mpi on multiple */
100: #define VECSHARED      "shared"
101: #define VECSEQVIENNACL "seqviennacl"
102: #define VECMPIVIENNACL "mpiviennacl"
103: #define VECVIENNACL    "viennacl" /* seqviennacl on one process and mpiviennacl on multiple */
104: #define VECSEQCUDA     "seqcuda"
105: #define VECMPICUDA     "mpicuda"
106: #define VECCUDA        "cuda" /* seqcuda on one process and mpicuda on multiple */
107: #define VECSEQHIP      "seqhip"
108: #define VECMPIHIP      "mpihip"
109: #define VECHIP         "hip" /* seqhip on one process and mpihip on multiple */
110: #define VECNEST        "nest"
111: #define VECSEQKOKKOS   "seqkokkos"
112: #define VECMPIKOKKOS   "mpikokkos"
113: #define VECKOKKOS      "kokkos" /* seqkokkos on one process and mpikokkos on multiple */

115: /* Dynamic creation and loading functions */
116: PETSC_EXTERN PetscErrorCode VecScatterSetType(VecScatter, VecScatterType);
117: PETSC_EXTERN PetscErrorCode VecScatterGetType(VecScatter, VecScatterType *);
118: PETSC_EXTERN PetscErrorCode VecScatterSetFromOptions(VecScatter);
119: PETSC_EXTERN PetscErrorCode VecScatterRegister(const char[], PetscErrorCode (*)(VecScatter));
120: PETSC_EXTERN PetscErrorCode VecScatterCreate(Vec, IS, Vec, IS, VecScatter *);

122: /* Logging support */
123: #define REAL_FILE_CLASSID 1211213
124: #define VEC_FILE_CLASSID  1211214
125: PETSC_EXTERN PetscClassId VEC_CLASSID;
126: PETSC_EXTERN PetscClassId PETSCSF_CLASSID;

128: PETSC_EXTERN PetscErrorCode VecInitializePackage(void);
129: PETSC_EXTERN PetscErrorCode VecFinalizePackage(void);

131: PETSC_EXTERN PetscErrorCode VecCreate(MPI_Comm, Vec *);
132: PETSC_EXTERN PetscErrorCode VecCreateFromOptions(MPI_Comm, const char *, PetscInt, PetscInt, PetscInt, Vec *);
133: PETSC_EXTERN PetscErrorCode VecCreateSeq(MPI_Comm, PetscInt, Vec *);
134: PETSC_EXTERN PetscErrorCode VecCreateMPI(MPI_Comm, PetscInt, PetscInt, Vec *);
135: PETSC_EXTERN PetscErrorCode VecCreateSeqWithArray(MPI_Comm, PetscInt, PetscInt, const PetscScalar[], Vec *);
136: PETSC_EXTERN PetscErrorCode VecCreateMPIWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar[], Vec *);
137: PETSC_EXTERN PetscErrorCode VecCreateShared(MPI_Comm, PetscInt, PetscInt, Vec *);

139: PETSC_EXTERN PetscErrorCode VecSetFromOptions(Vec);
140: PETSC_EXTERN PetscErrorCode VecViewFromOptions(Vec, PetscObject, const char[]);

142: PETSC_EXTERN PetscErrorCode VecSetUp(Vec);
143: PETSC_EXTERN PetscErrorCode VecDestroy(Vec *);
144: PETSC_EXTERN PetscErrorCode VecZeroEntries(Vec);
145: PETSC_EXTERN PetscErrorCode VecSetOptionsPrefix(Vec, const char[]);
146: PETSC_EXTERN PetscErrorCode VecAppendOptionsPrefix(Vec, const char[]);
147: PETSC_EXTERN PetscErrorCode VecGetOptionsPrefix(Vec, const char *[]);

149: PETSC_EXTERN PetscErrorCode VecSetSizes(Vec, PetscInt, PetscInt);

151: PETSC_EXTERN PetscErrorCode VecDotNorm2(Vec, Vec, PetscScalar *, PetscReal *);
152: PETSC_EXTERN PetscErrorCode VecDot(Vec, Vec, PetscScalar *);
153: PETSC_EXTERN PetscErrorCode VecDotRealPart(Vec, Vec, PetscReal *);
154: PETSC_EXTERN PetscErrorCode VecTDot(Vec, Vec, PetscScalar *);
155: PETSC_EXTERN PetscErrorCode VecMDot(Vec, PetscInt, const Vec[], PetscScalar[]);
156: PETSC_EXTERN PetscErrorCode VecMTDot(Vec, PetscInt, const Vec[], PetscScalar[]);
157: PETSC_EXTERN PetscErrorCode VecGetSubVector(Vec, IS, Vec *);
158: PETSC_EXTERN PetscErrorCode VecRestoreSubVector(Vec, IS, Vec *);
159: PETSC_EXTERN PetscErrorCode VecConcatenate(PetscInt, const Vec[], Vec *, IS *[]);

161: /*E
162:     NormType - determines what type of norm to compute

164:     Values:
165: +    `NORM_1`         - the one norm, $||v|| = \sum_i | v_i |$. $||A|| = \max_j || v_*j ||$, maximum column sum
166: .    `NORM_2`         - the two norm, $||v|| = sqrt(\sum_i |v_i|^2)$ (vectors only)
167: .    `NORM_FROBENIUS` - $||A|| = sqrt(\sum_ij |A_ij|^2)$, same as `NORM_2` for vectors
168: .    `NORM_INFINITY`  - $||v|| = \max_i |v_i|$. $||A|| = \max_i || v_i* ||$, maximum row sum
169: -    `NORM_1_AND_2`   - computes both the 1 and 2 norm of a vector. The values are stored in two adjacent `PetscReal` memory locations

171:     Level: beginner

173:     Note:
174:     The `v` above represents a `Vec` while the `A` represents a `Mat`

176: .seealso: [](ch_vectors), `Vec`, `Mat`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `MatNorm()`, `NORM_1`,
177:           `NORM_2`, `NORM_FROBENIUS`, `NORM_INFINITY`, `NORM_1_AND_2`
178: E*/
179: typedef enum {
180:   NORM_1         = 0,
181:   NORM_2         = 1,
182:   NORM_FROBENIUS = 2,
183:   NORM_INFINITY  = 3,
184:   NORM_1_AND_2   = 4
185: } NormType;
186: PETSC_EXTERN const char *const NormTypes[];
187: #define NORM_MAX NORM_INFINITY

189: /*MC
190:    NORM_1 - the one norm, $||v|| = \sum_i | v_i |$. $||A|| = \max_j || v_{*,j} ||$, maximum column sum

192:    Level: beginner

194: .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_2`, `NORM_FROBENIUS`,
195:           `NORM_INFINITY`, `NORM_1_AND_2`
196: M*/

198: /*MC
199:    NORM_2 - the two norm, $||v|| = \sqrt{\sum_i |v_i|^2}$ (vectors only)

201:    Level: beginner

203: .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_FROBENIUS`,
204:           `NORM_INFINITY`, `NORM_1_AND_2`
205: M*/

207: /*MC
208:    NORM_FROBENIUS - $||A|| = \sqrt{\sum_{i,j} |A_{i,j}|^2}$, same as `NORM_2` for vectors

210:    Level: beginner

212: .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
213:           `NORM_INFINITY`, `NORM_1_AND_2`
214: M*/

216: /*MC
217:    NORM_INFINITY - $||v|| = \max_i |v_i|$. $||A|| = \max_i || v_{i,*} ||$, maximum row sum

219:    Level: beginner

221: .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
222:           `NORM_FROBENIUS`, `NORM_1_AND_2`
223: M*/

225: /*MC
226:    NORM_1_AND_2 - computes both the 1 and 2 norm of a vector. The values are stored in two adjacent `PetscReal` memory locations

228:    Level: beginner

230: .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
231:           `NORM_FROBENIUS`, `NORM_INFINITY`
232: M*/

234: /*MC
235:    NORM_MAX - see `NORM_INFINITY`

237:    Level: beginner
238: M*/

240: /*E
241:     ReductionType - determines what type of column reduction (one that is not a type of norm defined in `NormType`) to compute

243:     Values:
244: +  `REDUCTION_SUM_REALPART`       - sum of real part of each matrix column
245: .  `REDUCTION_SUM_IMAGINARYPART`  - sum of imaginary part of each matrix column
246: .  `REDUCTION_MEAN_REALPART`      - arithmetic mean of real part of each matrix column
247: -  `REDUCTION_MEAN_IMAGINARYPART` - arithmetic mean of imaginary part of each matrix column

249:     Level: beginner

251:     Developer Note:
252:     The constants defined in `ReductionType` MUST BE DISTINCT from those defined in `NormType`.
253:    This is because `MatGetColumnReductions()` is used to compute both norms and other types of reductions,
254:    and the constants defined in both `NormType` and `ReductionType` are used to designate the desired operation.

256: .seealso: [](ch_vectors), `MatGetColumnReductions()`, `MatGetColumnNorms()`, `NormType`, `REDUCTION_SUM_REALPART`,
257:           `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_NORM_1`, `REDUCTION_NORM_2`, `REDUCTION_NORM_FROBENIUS`, `REDUCTION_NORM_INFINITY`
258: E*/
259: typedef enum {
260:   REDUCTION_SUM_REALPART       = 10,
261:   REDUCTION_MEAN_REALPART      = 11,
262:   REDUCTION_SUM_IMAGINARYPART  = 12,
263:   REDUCTION_MEAN_IMAGINARYPART = 13
264: } ReductionType;

266: /*MC
267:    REDUCTION_SUM_REALPART - sum of real part of matrix column

269:    Level: beginner

271: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_NORM_1`,
272:           `REDUCTION_NORM_2`, `REDUCTION_NORM_FROBENIUS`, `REDUCTION_NORM_INFINITY`
273: M*/

275: /*MC
276:    REDUCTION_SUM_IMAGINARYPART - sum of imaginary part of matrix column

278:    Level: beginner

280: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_SUM_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`, `REDUCTION_NORM_1`,
281:           `REDUCTION_NORM_2`, `REDUCTION_NORM_FROBENIUS`, `REDUCTION_NORM_INFINITY`
282: M*/

284: /*MC
285:    REDUCTION_MEAN_REALPART - arithmetic mean of real part of matrix column

287:    Level: beginner

289: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_IMAGINARYPART`, `REDUCTION_SUM_REALPART`, `REDUCTION_NORM_1`,
290:           `REDUCTION_NORM_2`, `REDUCTION_NORM_FROBENIUS`, `REDUCTION_NORM_INFINITY`
291: M*/

293: /*MC
294:    REDUCTION_MEAN_IMAGINARYPART - arithmetic mean of imaginary part of matrix column

296:    Level: beginner

298: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_REALPART`, `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_NORM_1`,
299:           `REDUCTION_NORM_2`, `REDUCTION_NORM_FROBENIUS`, `REDUCTION_NORM_INFINITY`
300: M*/

302: PETSC_EXTERN PetscErrorCode VecNorm(Vec, NormType, PetscReal *);
303: PETSC_EXTERN PetscErrorCode VecNormAvailable(Vec, NormType, PetscBool *, PetscReal *);
304: PETSC_EXTERN PetscErrorCode VecNormalize(Vec, PetscReal *);
305: PETSC_EXTERN PetscErrorCode VecSum(Vec, PetscScalar *);
306: PETSC_EXTERN PetscErrorCode VecMean(Vec, PetscScalar *);
307: PETSC_EXTERN PetscErrorCode VecMax(Vec, PetscInt *, PetscReal *);
308: PETSC_EXTERN PetscErrorCode VecMin(Vec, PetscInt *, PetscReal *);
309: PETSC_EXTERN PetscErrorCode VecScale(Vec, PetscScalar);
310: PETSC_EXTERN PetscErrorCode VecCopy(Vec, Vec);
311: PETSC_EXTERN PetscErrorCode VecSetRandom(Vec, PetscRandom);
312: PETSC_EXTERN PetscErrorCode VecSet(Vec, PetscScalar);
313: PETSC_EXTERN PetscErrorCode VecSetInf(Vec);
314: PETSC_EXTERN PetscErrorCode VecSwap(Vec, Vec);
315: PETSC_EXTERN PetscErrorCode VecAXPY(Vec, PetscScalar, Vec);
316: PETSC_EXTERN PetscErrorCode VecAXPBY(Vec, PetscScalar, PetscScalar, Vec);
317: PETSC_EXTERN PetscErrorCode VecMAXPY(Vec, PetscInt, const PetscScalar[], Vec[]);
318: PETSC_EXTERN PetscErrorCode VecMAXPBY(Vec, PetscInt, const PetscScalar[], PetscScalar, Vec[]);
319: PETSC_EXTERN PetscErrorCode VecAYPX(Vec, PetscScalar, Vec);
320: PETSC_EXTERN PetscErrorCode VecWAXPY(Vec, PetscScalar, Vec, Vec);
321: PETSC_EXTERN PetscErrorCode VecAXPBYPCZ(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec);
322: PETSC_EXTERN PetscErrorCode VecPointwiseMax(Vec, Vec, Vec);
323: PETSC_EXTERN PetscErrorCode VecPointwiseMaxAbs(Vec, Vec, Vec);
324: PETSC_EXTERN PetscErrorCode VecPointwiseMin(Vec, Vec, Vec);
325: PETSC_EXTERN PetscErrorCode VecPointwiseMult(Vec, Vec, Vec);
326: PETSC_EXTERN PetscErrorCode VecPointwiseDivide(Vec, Vec, Vec);
327: PETSC_EXTERN PetscErrorCode VecMaxPointwiseDivide(Vec, Vec, PetscReal *);
328: PETSC_EXTERN PetscErrorCode VecShift(Vec, PetscScalar);
329: PETSC_EXTERN PetscErrorCode VecReciprocal(Vec);
330: PETSC_EXTERN PetscErrorCode VecPermute(Vec, IS, PetscBool);
331: PETSC_EXTERN PetscErrorCode VecSqrtAbs(Vec);
332: PETSC_EXTERN PetscErrorCode VecLog(Vec);
333: PETSC_EXTERN PetscErrorCode VecExp(Vec);
334: PETSC_EXTERN PetscErrorCode VecAbs(Vec);
335: PETSC_EXTERN PetscErrorCode VecDuplicate(Vec, Vec *);
336: PETSC_EXTERN PetscErrorCode VecDuplicateVecs(Vec, PetscInt, Vec *[]);
337: PETSC_EXTERN PetscErrorCode VecDestroyVecs(PetscInt, Vec *[]);
338: PETSC_EXTERN PetscErrorCode VecStrideNormAll(Vec, NormType, PetscReal[]);
339: PETSC_EXTERN PetscErrorCode VecStrideMaxAll(Vec, PetscInt[], PetscReal[]);
340: PETSC_EXTERN PetscErrorCode VecStrideMinAll(Vec, PetscInt[], PetscReal[]);
341: PETSC_EXTERN PetscErrorCode VecStrideScaleAll(Vec, const PetscScalar[]);
342: PETSC_EXTERN PetscErrorCode VecStrideSumAll(Vec, PetscScalar *);
343: PETSC_EXTERN PetscErrorCode VecUniqueEntries(Vec, PetscInt *, PetscScalar **);

345: PETSC_EXTERN PetscErrorCode VecStrideNorm(Vec, PetscInt, NormType, PetscReal *);
346: PETSC_EXTERN PetscErrorCode VecStrideMax(Vec, PetscInt, PetscInt *, PetscReal *);
347: PETSC_EXTERN PetscErrorCode VecStrideMin(Vec, PetscInt, PetscInt *, PetscReal *);
348: PETSC_EXTERN PetscErrorCode VecStrideScale(Vec, PetscInt, PetscScalar);
349: PETSC_EXTERN PetscErrorCode VecStrideSum(Vec, PetscInt, PetscScalar *);
350: PETSC_EXTERN PetscErrorCode VecStrideSet(Vec, PetscInt, PetscScalar);

352: PETSC_EXTERN PetscErrorCode VecStrideGather(Vec, PetscInt, Vec, InsertMode);
353: PETSC_EXTERN PetscErrorCode VecStrideScatter(Vec, PetscInt, Vec, InsertMode);
354: PETSC_EXTERN PetscErrorCode VecStrideGatherAll(Vec, Vec[], InsertMode);
355: PETSC_EXTERN PetscErrorCode VecStrideScatterAll(Vec[], Vec, InsertMode);

357: PETSC_EXTERN PetscErrorCode VecStrideSubSetScatter(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
358: PETSC_EXTERN PetscErrorCode VecStrideSubSetGather(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);

360: PETSC_EXTERN PetscErrorCode VecSetValues(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
361: PETSC_EXTERN PetscErrorCode VecGetValues(Vec, PetscInt, const PetscInt[], PetscScalar[]);
362: PETSC_EXTERN PetscErrorCode VecAssemblyBegin(Vec);
363: PETSC_EXTERN PetscErrorCode VecAssemblyEnd(Vec);
364: PETSC_EXTERN PetscErrorCode VecStashSetInitialSize(Vec, PetscInt, PetscInt);
365: PETSC_EXTERN PetscErrorCode VecStashView(Vec, PetscViewer);
366: PETSC_EXTERN PetscErrorCode VecStashViewFromOptions(Vec, PetscObject, const char[]);
367: PETSC_EXTERN PetscErrorCode VecStashGetInfo(Vec, PetscInt *, PetscInt *, PetscInt *, PetscInt *);

369: PETSC_EXTERN PetscErrorCode VecSetPreallocationCOO(Vec, PetscCount, const PetscInt[]);
370: PETSC_EXTERN PetscErrorCode VecSetPreallocationCOOLocal(Vec, PetscCount, PetscInt[]);
371: PETSC_EXTERN PetscErrorCode VecSetValuesCOO(Vec, const PetscScalar[], InsertMode);

373: /*MC
374:    VecSetValue - Set a single entry into a vector.

376:    Synopsis:
377: #include <petscvec.h>
378:    PetscErrorCode VecSetValue(Vec v,PetscInt row,PetscScalar value, InsertMode mode);

380:    Not Collective

382:    Input Parameters:
383: +  v     - the vector
384: .  row   - the row location of the entry
385: .  value - the value to insert
386: -  mode  - either `INSERT_VALUES` or `ADD_VALUES`

388:    Level: beginner

390:    Notes:
391:    For efficiency one should use `VecSetValues()` and set several or
392:    many values simultaneously if possible.

394:    These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
395:    MUST be called after all calls to `VecSetValue()` have been completed.

397:    `VecSetValue()` uses 0-based indices in Fortran as well as in C.

399: .seealso: [](ch_vectors), `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValueLocal()`
400: M*/
401: static inline PetscErrorCode VecSetValue(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
402: {
403:   return VecSetValues(v, 1, &i, &va, mode);
404: }

406: PETSC_EXTERN PetscErrorCode VecSetBlockSize(Vec, PetscInt);
407: PETSC_EXTERN PetscErrorCode VecGetBlockSize(Vec, PetscInt *);
408: PETSC_EXTERN PetscErrorCode VecSetValuesBlocked(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);

410: /* Dynamic creation and loading functions */
411: PETSC_EXTERN PetscFunctionList VecList;
412: PETSC_EXTERN PetscErrorCode    VecSetType(Vec, VecType);
413: PETSC_EXTERN PetscErrorCode    VecGetType(Vec, VecType *);
414: PETSC_EXTERN PetscErrorCode    VecRegister(const char[], PetscErrorCode (*)(Vec));
415: PETSC_EXTERN PetscErrorCode    VecRegisterAll(void);

417: PETSC_EXTERN PetscErrorCode VecScatterBegin(VecScatter, Vec, Vec, InsertMode, ScatterMode);
418: PETSC_EXTERN PetscErrorCode VecScatterEnd(VecScatter, Vec, Vec, InsertMode, ScatterMode);
419: PETSC_EXTERN PetscErrorCode VecScatterDestroy(VecScatter *);
420: PETSC_EXTERN PetscErrorCode VecScatterSetUp(VecScatter);
421: PETSC_EXTERN PetscErrorCode VecScatterCopy(VecScatter, VecScatter *);
422: PETSC_EXTERN PetscErrorCode VecScatterView(VecScatter, PetscViewer);
423: PETSC_EXTERN PetscErrorCode VecScatterViewFromOptions(VecScatter, PetscObject, const char[]);
424: PETSC_EXTERN PetscErrorCode VecScatterRemap(VecScatter, PetscInt[], PetscInt[]);
425: PETSC_EXTERN PetscErrorCode VecScatterGetMerged(VecScatter, PetscBool *);

427: PETSC_EXTERN PetscErrorCode VecGetArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
428: PETSC_EXTERN PetscErrorCode VecRestoreArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
429: PETSC_EXTERN PetscErrorCode VecGetArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
430: PETSC_EXTERN PetscErrorCode VecRestoreArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
431: PETSC_EXTERN PetscErrorCode VecGetArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
432: PETSC_EXTERN PetscErrorCode VecRestoreArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
433: PETSC_EXTERN PetscErrorCode VecGetArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);
434: PETSC_EXTERN PetscErrorCode VecRestoreArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);

436: PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
437: PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
438: PETSC_EXTERN PetscErrorCode VecRestoreArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
439: PETSC_EXTERN PetscErrorCode VecGetArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
440: PETSC_EXTERN PetscErrorCode VecRestoreArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
441: PETSC_EXTERN PetscErrorCode VecGetArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
442: PETSC_EXTERN PetscErrorCode VecRestoreArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
443: PETSC_EXTERN PetscErrorCode VecGetArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);
444: PETSC_EXTERN PetscErrorCode VecRestoreArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);

446: PETSC_EXTERN PetscErrorCode VecGetArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
447: PETSC_EXTERN PetscErrorCode VecRestoreArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
448: PETSC_EXTERN PetscErrorCode VecGetArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
449: PETSC_EXTERN PetscErrorCode VecRestoreArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
450: PETSC_EXTERN PetscErrorCode VecGetArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
451: PETSC_EXTERN PetscErrorCode VecRestoreArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
452: PETSC_EXTERN PetscErrorCode VecGetArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);
453: PETSC_EXTERN PetscErrorCode VecRestoreArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);

455: PETSC_EXTERN PetscErrorCode VecPlaceArray(Vec, const PetscScalar[]);
456: PETSC_EXTERN PetscErrorCode VecResetArray(Vec);
457: PETSC_EXTERN PetscErrorCode VecReplaceArray(Vec, const PetscScalar[]);

459: PETSC_EXTERN PetscErrorCode VecGetArrays(const Vec[], PetscInt, PetscScalar **[]);
460: PETSC_EXTERN PetscErrorCode VecRestoreArrays(const Vec[], PetscInt, PetscScalar **[]);

462: PETSC_EXTERN PetscErrorCode VecView(Vec, PetscViewer);
463: PETSC_EXTERN PetscErrorCode VecViewNative(Vec, PetscViewer);
464: PETSC_EXTERN PetscErrorCode VecEqual(Vec, Vec, PetscBool *);
465: PETSC_EXTERN PetscErrorCode VecLoad(Vec, PetscViewer);

467: PETSC_EXTERN PetscErrorCode VecGetSize(Vec, PetscInt *);
468: PETSC_EXTERN PetscErrorCode VecGetLocalSize(Vec, PetscInt *);
469: PETSC_EXTERN PetscErrorCode VecGetOwnershipRange(Vec, PetscInt *, PetscInt *);
470: PETSC_EXTERN PetscErrorCode VecGetOwnershipRanges(Vec, const PetscInt *[]);

472: PETSC_EXTERN PetscErrorCode VecSetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping);
473: PETSC_EXTERN PetscErrorCode VecSetValuesLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);

475: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLContext(Vec, PETSC_UINTPTR_T *);
476: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLQueue(Vec, PETSC_UINTPTR_T *);
477: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemRead(Vec, PETSC_UINTPTR_T *);
478: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemWrite(Vec, PETSC_UINTPTR_T *);
479: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMemWrite(Vec);
480: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMem(Vec, PETSC_UINTPTR_T *);
481: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMem(Vec);

483: /*MC
484:    VecSetValueLocal - Set a single entry into a vector using the local numbering, see `VecSetValuesLocal()`

486:    Synopsis:
487: #include <petscvec.h>
488:    PetscErrorCode VecSetValueLocal(Vec v,PetscInt row,PetscScalar value, InsertMode mode);

490:    Not Collective

492:    Input Parameters:
493: +  v     - the vector
494: .  row   - the row location of the entry
495: .  value - the value to insert
496: -  mode  - either `INSERT_VALUES` or `ADD_VALUES`

498:    Level: beginner

500:    Notes:
501:    For efficiency one should use `VecSetValuesLocal()` and set several or
502:    many values simultaneously if possible.

504:    These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
505:    MUST be called after all calls to `VecSetValueLocal()` have been completed.

507:    `VecSetValues()` uses 0-based indices in Fortran as well as in C.

509: .seealso: [](ch_vectors), `VecSetValuesLocal()`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValue()`
510: M*/
511: static inline PetscErrorCode VecSetValueLocal(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
512: {
513:   return VecSetValuesLocal(v, 1, &i, &va, mode);
514: }

516: /*MC
517:    VecCheckAssembled - checks if values have been changed in the vector, by `VecSetValues()` or related routines,  but it has not been assembled

519:    Synopsis:
520: #include <petscvec.h>
521:    VecCheckAssembled(Vec v);

523:    Not Collective

525:    Input Parameter:
526: .  v - the vector to check

528:    Level: developer

530:    Note:
531:    After calls to `VecSetValues()` and related routines on must call ``VecAssemblyBegin()` and `VecAssemblyEnd()` before using the vector

533: .seealso: [](ch_vectors), `Vec`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`
534: M*/
535: #define VecCheckAssembled(v) PetscCheck(v->stash.insertmode == NOT_SET_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled vector, did you call VecAssemblyBegin()/VecAssemblyEnd()?")

537: PETSC_EXTERN PetscErrorCode VecSetValuesBlockedLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
538: PETSC_EXTERN PetscErrorCode VecGetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping *);

540: PETSC_EXTERN PetscErrorCode VecDotBegin(Vec, Vec, PetscScalar *);
541: PETSC_EXTERN PetscErrorCode VecDotEnd(Vec, Vec, PetscScalar *);
542: PETSC_EXTERN PetscErrorCode VecTDotBegin(Vec, Vec, PetscScalar *);
543: PETSC_EXTERN PetscErrorCode VecTDotEnd(Vec, Vec, PetscScalar *);
544: PETSC_EXTERN PetscErrorCode VecNormBegin(Vec, NormType, PetscReal *);
545: PETSC_EXTERN PetscErrorCode VecNormEnd(Vec, NormType, PetscReal *);
546: PETSC_EXTERN PetscErrorCode VecErrorWeightedNorms(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *);

548: PETSC_EXTERN PetscErrorCode VecMDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
549: PETSC_EXTERN PetscErrorCode VecMDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
550: PETSC_EXTERN PetscErrorCode VecMTDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
551: PETSC_EXTERN PetscErrorCode VecMTDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
552: PETSC_EXTERN PetscErrorCode PetscCommSplitReductionBegin(MPI_Comm);

554: PETSC_EXTERN PetscErrorCode VecBindToCPU(Vec, PetscBool);
555: PETSC_DEPRECATED_FUNCTION(3, 13, 0, "VecBindToCPU()", ) static inline PetscErrorCode VecPinToCPU(Vec v, PetscBool flg)
556: {
557:   return VecBindToCPU(v, flg);
558: }
559: PETSC_EXTERN PetscErrorCode VecBoundToCPU(Vec, PetscBool *);
560: PETSC_EXTERN PetscErrorCode VecSetBindingPropagates(Vec, PetscBool);
561: PETSC_EXTERN PetscErrorCode VecGetBindingPropagates(Vec, PetscBool *);
562: PETSC_EXTERN PetscErrorCode VecSetPinnedMemoryMin(Vec, size_t);
563: PETSC_EXTERN PetscErrorCode VecGetPinnedMemoryMin(Vec, size_t *);

565: PETSC_EXTERN PetscErrorCode VecGetOffloadMask(Vec, PetscOffloadMask *);

567: typedef enum {
568:   VEC_IGNORE_OFF_PROC_ENTRIES,
569:   VEC_IGNORE_NEGATIVE_INDICES,
570:   VEC_SUBSET_OFF_PROC_ENTRIES
571: } VecOption;
572: PETSC_EXTERN PetscErrorCode VecSetOption(Vec, VecOption, PetscBool);

574: PETSC_EXTERN PetscErrorCode VecGetArray(Vec, PetscScalar **);
575: PETSC_EXTERN PetscErrorCode VecGetArrayWrite(Vec, PetscScalar **);
576: PETSC_EXTERN PetscErrorCode VecGetArrayRead(Vec, const PetscScalar **);
577: PETSC_EXTERN PetscErrorCode VecRestoreArray(Vec, PetscScalar **);
578: PETSC_EXTERN PetscErrorCode VecRestoreArrayWrite(Vec, PetscScalar **);
579: PETSC_EXTERN PetscErrorCode VecRestoreArrayRead(Vec, const PetscScalar **);
580: PETSC_EXTERN PetscErrorCode VecCreateLocalVector(Vec, Vec *);
581: PETSC_EXTERN PetscErrorCode VecGetLocalVector(Vec, Vec);
582: PETSC_EXTERN PetscErrorCode VecRestoreLocalVector(Vec, Vec);
583: PETSC_EXTERN PetscErrorCode VecGetLocalVectorRead(Vec, Vec);
584: PETSC_EXTERN PetscErrorCode VecRestoreLocalVectorRead(Vec, Vec);
585: PETSC_EXTERN PetscErrorCode VecGetArrayAndMemType(Vec, PetscScalar **, PetscMemType *);
586: PETSC_EXTERN PetscErrorCode VecRestoreArrayAndMemType(Vec, PetscScalar **);
587: PETSC_EXTERN PetscErrorCode VecGetArrayReadAndMemType(Vec, const PetscScalar **, PetscMemType *);
588: PETSC_EXTERN PetscErrorCode VecRestoreArrayReadAndMemType(Vec, const PetscScalar **);
589: PETSC_EXTERN PetscErrorCode VecGetArrayWriteAndMemType(Vec, PetscScalar **, PetscMemType *);
590: PETSC_EXTERN PetscErrorCode VecRestoreArrayWriteAndMemType(Vec, PetscScalar **);

592: /*@C
593:    VecGetArrayPair - Accesses a pair of pointers for two vectors that may be common. When the vectors are not the same the first pointer is read only

595:    Logically Collective; No Fortran Support

597:    Input Parameters:
598: +  x - the vector
599: -  y - the second vector

601:    Output Parameters:
602: +  xv - location to put pointer to the first array
603: -  yv - location to put pointer to the second array

605:    Level: developer

607: .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecRestoreArrayPair()`
608: @*/
609: static inline PetscErrorCode VecGetArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
610: {
611:   PetscFunctionBegin;
612:   PetscCall(VecGetArray(y, yv));
613:   if (x == y) *xv = *yv;
614:   else PetscCall(VecGetArrayRead(x, (const PetscScalar **)xv));
615:   PetscFunctionReturn(PETSC_SUCCESS);
616: }

618: /*@C
619:    VecRestoreArrayPair - Returns a pair of pointers for two vectors that may be common obtained with `VecGetArrayPair()`

621:    Logically Collective; No Fortran Support

623:    Input Parameters:
624: +  x  - the vector
625: .  y  - the second vector
626: .  xv - location to put pointer to the first array
627: -  yv - location to put pointer to the second array

629:    Level: developer

631: .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecGetArrayPair()`
632: @*/
633: static inline PetscErrorCode VecRestoreArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
634: {
635:   PetscFunctionBegin;
636:   PetscCall(VecRestoreArray(y, yv));
637:   if (x != y) PetscCall(VecRestoreArrayRead(x, (const PetscScalar **)xv));
638:   PetscFunctionReturn(PETSC_SUCCESS);
639: }

641: #if PetscDefined(USE_DEBUG)
642: PETSC_EXTERN PetscErrorCode  VecLockReadPush(Vec);
643: PETSC_EXTERN PetscErrorCode  VecLockReadPop(Vec);
644: PETSC_EXTERN PetscErrorCode  VecLockWriteSet(Vec, PetscBool);
645: PETSC_EXTERN PetscErrorCode  VecLockGet(Vec, PetscInt *);
646: PETSC_EXTERN PetscErrorCode  VecLockGetLocation(Vec, const char *[], const char *[], int *);
647: static inline PetscErrorCode VecSetErrorIfLocked(Vec x, PetscInt arg)
648: {
649:   PetscInt state;

651:   PetscFunctionBegin;
652:   PetscCall(VecLockGet(x, &state));
653:   if (PetscUnlikely(state != 0)) {
654:     const char *file, *func, *name;
655:     int         line;

657:     PetscCall(VecLockGetLocation(x, &file, &func, &line));
658:     PetscCall(PetscObjectGetName((PetscObject)x, &name));
659:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector '%s' (argument #%" PetscInt_FMT ") was locked for %s access in %s() at %s:%d (line numbers only accurate to function begin)", name, arg, state > 0 ? "read-only" : "write-only", func ? func : "unknown_function", file ? file : "unknown file", line);
660:   }
661:   PetscFunctionReturn(PETSC_SUCCESS);
662: }
663: /* The three are deprecated */
664: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPush()", ) PetscErrorCode VecLockPush(Vec);
665: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPop()", ) PetscErrorCode VecLockPop(Vec);
666:   #define VecLocked(x, arg) VecSetErrorIfLocked(x, arg) PETSC_DEPRECATED_MACRO(3, 11, 0, "VecSetErrorIfLocked()", )
667: #else
668:   #define VecLockReadPush(x)          PETSC_SUCCESS
669:   #define VecLockReadPop(x)           PETSC_SUCCESS
670:   #define VecLockGet(x, s)            (*(s) = 0, PETSC_SUCCESS)
671:   #define VecSetErrorIfLocked(x, arg) PETSC_SUCCESS
672:   #define VecLockWriteSet(x, flg)     PETSC_SUCCESS
673:   /* The three are deprecated */
674:   #define VecLockPush(x)              PETSC_SUCCESS
675:   #define VecLockPop(x)               PETSC_SUCCESS
676:   #define VecLocked(x, arg)           PETSC_SUCCESS
677: #endif

679: /*E
680:   VecOperation - Enumeration of overide-able methods in the `Vec` implementation function-table.

682:   Values:
683: + `VECOP_DUPLICATE`  - `VecDuplicate()`
684: . `VECOP_SET`        - `VecSet()`
685: . `VECOP_VIEW`       - `VecView()`
686: . `VECOP_LOAD`       - `VecLoad()`
687: . `VECOP_VIEWNATIVE` - `VecViewNative()`
688: - `VECOP_LOADNATIVE` - `VecLoadNative()`

690:   Level: advanced

692:   Notes:
693:   Some operations may serve as the implementation for other routines not listed above. For
694:   example `VECOP_SET` can be used to simultaneously overriding the implementation used in
695:   `VecSet()`, `VecSetInf()`, and `VecZeroEntries()`.

697:   Entries to `VecOperation` are added as needed so if you do not see the operation listed which
698:   you'd like to replace, please send mail to `petsc-maint@mcs.anl.gov`!

700: .seealso: [](ch_vectors), `Vec`, `VecSetOperation()`
701: E*/
702: typedef enum {
703:   VECOP_DUPLICATE  = 0,
704:   VECOP_SET        = 10,
705:   VECOP_VIEW       = 33,
706:   VECOP_LOAD       = 41,
707:   VECOP_VIEWNATIVE = 69,
708:   VECOP_LOADNATIVE = 70
709: } VecOperation;
710: PETSC_EXTERN PetscErrorCode VecSetOperation(Vec, VecOperation, void (*)(void));

712: /*
713:      Routines for dealing with ghosted vectors:
714:   vectors with ghost elements at the end of the array.
715: */
716: PETSC_EXTERN PetscErrorCode VecMPISetGhost(Vec, PetscInt, const PetscInt[]);
717: PETSC_EXTERN PetscErrorCode VecCreateGhost(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
718: PETSC_EXTERN PetscErrorCode VecCreateGhostWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
719: PETSC_EXTERN PetscErrorCode VecCreateGhostBlock(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
720: PETSC_EXTERN PetscErrorCode VecCreateGhostBlockWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
721: PETSC_EXTERN PetscErrorCode VecGhostGetGhostIS(Vec, IS *);
722: PETSC_EXTERN PetscErrorCode VecGhostGetLocalForm(Vec, Vec *);
723: PETSC_EXTERN PetscErrorCode VecGhostRestoreLocalForm(Vec, Vec *);
724: PETSC_EXTERN PetscErrorCode VecGhostIsLocalForm(Vec, Vec, PetscBool *);
725: PETSC_EXTERN PetscErrorCode VecGhostUpdateBegin(Vec, InsertMode, ScatterMode);
726: PETSC_EXTERN PetscErrorCode VecGhostUpdateEnd(Vec, InsertMode, ScatterMode);

728: PETSC_EXTERN PetscErrorCode VecConjugate(Vec);
729: PETSC_EXTERN PetscErrorCode VecImaginaryPart(Vec);
730: PETSC_EXTERN PetscErrorCode VecRealPart(Vec);

732: PETSC_EXTERN PetscErrorCode VecScatterCreateToAll(Vec, VecScatter *, Vec *);
733: PETSC_EXTERN PetscErrorCode VecScatterCreateToZero(Vec, VecScatter *, Vec *);

735: /* These vector calls were included from TAO. They miss vtable entries and GPU implementation */
736: PETSC_EXTERN PetscErrorCode ISComplementVec(IS, Vec, IS *);
737: PETSC_EXTERN PetscErrorCode VecPow(Vec, PetscScalar);
738: PETSC_EXTERN PetscErrorCode VecMedian(Vec, Vec, Vec, Vec);
739: PETSC_EXTERN PetscErrorCode VecWhichInactive(Vec, Vec, Vec, Vec, PetscBool, IS *);
740: PETSC_EXTERN PetscErrorCode VecWhichBetween(Vec, Vec, Vec, IS *);
741: PETSC_EXTERN PetscErrorCode VecWhichBetweenOrEqual(Vec, Vec, Vec, IS *);
742: PETSC_EXTERN PetscErrorCode VecWhichGreaterThan(Vec, Vec, IS *);
743: PETSC_EXTERN PetscErrorCode VecWhichLessThan(Vec, Vec, IS *);
744: PETSC_EXTERN PetscErrorCode VecWhichEqual(Vec, Vec, IS *);
745: PETSC_EXTERN PetscErrorCode VecISAXPY(Vec, IS, PetscScalar, Vec);
746: PETSC_EXTERN PetscErrorCode VecISCopy(Vec, IS, ScatterMode, Vec);
747: PETSC_EXTERN PetscErrorCode VecISSet(Vec, IS, PetscScalar);
748: PETSC_EXTERN PetscErrorCode VecISShift(Vec, IS, PetscScalar);
749: PETSC_EXTERN PetscErrorCode VecBoundGradientProjection(Vec, Vec, Vec, Vec, Vec);
750: PETSC_EXTERN PetscErrorCode VecStepBoundInfo(Vec, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
751: PETSC_EXTERN PetscErrorCode VecStepMax(Vec, Vec, PetscReal *);
752: PETSC_EXTERN PetscErrorCode VecStepMaxBounded(Vec, Vec, Vec, Vec, PetscReal *);

754: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaGetVector(PetscViewer, Vec);
755: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutVector(PetscViewer, Vec);

757: /*S
758:    Vecs - Collection of vectors where the data for the vectors is stored in
759:           one contiguous memory

761:    Level: advanced

763:    Notes:
764:    Temporary construct for handling multiply right-hand side solves

766:    This is faked by storing a single vector that has enough array space for
767:    n vectors

769: S*/
770: struct _n_Vecs {
771:   PetscInt n;
772:   Vec      v;
773: };
774: typedef struct _n_Vecs     *Vecs;
775: PETSC_EXTERN PetscErrorCode VecsDestroy(Vecs);
776: PETSC_EXTERN PetscErrorCode VecsCreateSeq(MPI_Comm, PetscInt, PetscInt, Vecs *);
777: PETSC_EXTERN PetscErrorCode VecsCreateSeqWithArray(MPI_Comm, PetscInt, PetscInt, PetscScalar *, Vecs *);
778: PETSC_EXTERN PetscErrorCode VecsDuplicate(Vecs, Vecs *);

780: #if PetscDefined(HAVE_VIENNACL)
781: typedef struct _p_PetscViennaCLIndices *PetscViennaCLIndices;
782: PETSC_EXTERN PetscErrorCode             VecViennaCLCopyToGPUSome_Public(Vec, PetscViennaCLIndices);
783: PETSC_EXTERN PetscErrorCode             VecViennaCLCopyFromGPUSome_Public(Vec, PetscViennaCLIndices);
784: PETSC_EXTERN PetscErrorCode             VecCreateSeqViennaCL(MPI_Comm, PetscInt, Vec *);
785: PETSC_EXTERN PetscErrorCode             VecCreateMPIViennaCL(MPI_Comm, PetscInt, PetscInt, Vec *);
786: #endif
787: #if PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP)
788: PETSC_EXTERN PetscErrorCode VecScatterInitializeForGPU(VecScatter, Vec);
789: PETSC_EXTERN PetscErrorCode VecScatterFinalizeForGPU(VecScatter);
790: #endif
791: #if PetscDefined(HAVE_KOKKOS_KERNELS)
792: PETSC_EXTERN PetscErrorCode VecCreateSeqKokkos(MPI_Comm, PetscInt, Vec *);
793: PETSC_EXTERN PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm, PetscInt, PetscInt, const PetscScalar *, Vec *);
794: PETSC_EXTERN PetscErrorCode VecCreateMPIKokkos(MPI_Comm, PetscInt, PetscInt, Vec *);
795: PETSC_EXTERN PetscErrorCode VecCreateMPIKokkosWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar *, Vec *);
796: #endif

798: PETSC_EXTERN PetscErrorCode VecNestGetSubVecs(Vec, PetscInt *, Vec **);
799: PETSC_EXTERN PetscErrorCode VecNestGetSubVec(Vec, PetscInt, Vec *);
800: PETSC_EXTERN PetscErrorCode VecNestSetSubVecs(Vec, PetscInt, PetscInt *, Vec *);
801: PETSC_EXTERN PetscErrorCode VecNestSetSubVec(Vec, PetscInt, Vec);
802: PETSC_EXTERN PetscErrorCode VecCreateNest(MPI_Comm, PetscInt, IS *, Vec *, Vec *);
803: PETSC_EXTERN PetscErrorCode VecNestGetSize(Vec, PetscInt *);

805: PETSC_EXTERN PetscErrorCode PetscOptionsGetVec(PetscOptions, const char[], const char[], Vec, PetscBool *);
806: PETSC_EXTERN PetscErrorCode VecFilter(Vec, PetscReal);
807: PETSC_DEPRECATED_FUNCTION(3, 20, 0, "VecFilter()", ) static inline PetscErrorCode VecChop(Vec v, PetscReal tol)
808: {
809:   return VecFilter(v, tol);
810: }

812: PETSC_EXTERN PetscErrorCode VecGetLayout(Vec, PetscLayout *);
813: PETSC_EXTERN PetscErrorCode VecSetLayout(Vec, PetscLayout);

815: PETSC_EXTERN PetscErrorCode PetscSectionVecView(PetscSection, Vec, PetscViewer);
816: PETSC_EXTERN PetscErrorCode VecGetValuesSection(Vec, PetscSection, PetscInt, PetscScalar **);
817: PETSC_EXTERN PetscErrorCode VecSetValuesSection(Vec, PetscSection, PetscInt, const PetscScalar[], InsertMode);
818: PETSC_EXTERN PetscErrorCode PetscSectionVecNorm(PetscSection, PetscSection, Vec, NormType, PetscReal[]);

820: /*S
821:   VecTagger - Object used to manage the tagging of a subset of indices based on the values of a vector.  The
822:               motivating application is the selection of cells for refinement or coarsening based on vector containing
823:               the values in an error indicator metric.

825:   Values:
826: +  `VECTAGGERABSOLUTE` - "absolute" values are in a interval (box for complex values) of explicitly defined values
827: .  `VECTAGGERRELATIVE` - "relative" values are in a interval (box for complex values) of values relative to the set of all values in the vector
828: .  `VECTAGGERCDF`      - "cdf" values are in a relative range of the *cumulative distribution* of values in the vector
829: .  `VECTAGGEROR`       - "or" values are in the union of other tags
830: -  `VECTAGGERAND`      - "and" values are in the intersection of other tags

832:   Level: advanced

834:   Developer Note:
835:   Why not use a `DMLabel` or similar object

837: .seealso: [](ch_vectors), `Vec`, `VecTaggerType`, `VecTaggerCreate()`
838: S*/
839: typedef struct _p_VecTagger *VecTagger;

841: /*J
842:   VecTaggerType - String with the name of a `VecTagger` type

844:   Level: advanced

846: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerCreate()`
847: J*/
848: typedef const char *VecTaggerType;
849: #define VECTAGGERABSOLUTE "absolute"
850: #define VECTAGGERRELATIVE "relative"
851: #define VECTAGGERCDF      "cdf"
852: #define VECTAGGEROR       "or"
853: #define VECTAGGERAND      "and"

855: PETSC_EXTERN PetscClassId      VEC_TAGGER_CLASSID;
856: PETSC_EXTERN PetscFunctionList VecTaggerList;
857: PETSC_EXTERN PetscErrorCode    VecTaggerRegister(const char[], PetscErrorCode (*)(VecTagger));
858: PETSC_EXTERN PetscErrorCode    VecTaggerRegisterAll(void);

860: PETSC_EXTERN PetscErrorCode VecTaggerCreate(MPI_Comm, VecTagger *);
861: PETSC_EXTERN PetscErrorCode VecTaggerSetBlockSize(VecTagger, PetscInt);
862: PETSC_EXTERN PetscErrorCode VecTaggerGetBlockSize(VecTagger, PetscInt *);
863: PETSC_EXTERN PetscErrorCode VecTaggerSetType(VecTagger, VecTaggerType);
864: PETSC_EXTERN PetscErrorCode VecTaggerGetType(VecTagger, VecTaggerType *);
865: PETSC_EXTERN PetscErrorCode VecTaggerSetInvert(VecTagger, PetscBool);
866: PETSC_EXTERN PetscErrorCode VecTaggerGetInvert(VecTagger, PetscBool *);
867: PETSC_EXTERN PetscErrorCode VecTaggerSetFromOptions(VecTagger);
868: PETSC_EXTERN PetscErrorCode VecTaggerSetUp(VecTagger);
869: PETSC_EXTERN PetscErrorCode VecTaggerView(VecTagger, PetscViewer);
870: PETSC_EXTERN PetscErrorCode VecTaggerComputeIS(VecTagger, Vec, IS *, PetscBool *);
871: PETSC_EXTERN PetscErrorCode VecTaggerDestroy(VecTagger *);

873: /*S
874:    VecTaggerBox - A interval (box for complex numbers) range used to tag values.  For real scalars, this is just a closed interval; for complex scalars,
875:    the box is the closed region in the complex plane such that real(min) <= real(z) <= real(max) and imag(min) <= imag(z) <= imag(max).  `INF` is an acceptable endpoint.

877:    Level: beginner

879: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerComputeIntervals()`
880: S*/
881: typedef struct {
882:   PetscScalar min;
883:   PetscScalar max;
884: } VecTaggerBox;
885: PETSC_EXTERN PetscErrorCode VecTaggerComputeBoxes(VecTagger, Vec, PetscInt *, VecTaggerBox **, PetscBool *);

887: PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteSetBox(VecTagger, VecTaggerBox *);
888: PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteGetBox(VecTagger, const VecTaggerBox **);

890: PETSC_EXTERN PetscErrorCode VecTaggerRelativeSetBox(VecTagger, VecTaggerBox *);
891: PETSC_EXTERN PetscErrorCode VecTaggerRelativeGetBox(VecTagger, const VecTaggerBox **);

893: PETSC_EXTERN PetscErrorCode VecTaggerCDFSetBox(VecTagger, VecTaggerBox *);
894: PETSC_EXTERN PetscErrorCode VecTaggerCDFGetBox(VecTagger, const VecTaggerBox **);

896: /*E
897:   VecTaggerCDFMethod - Determines what method is used to compute absolute values from cumulative distribution values (e.g., what value is the preimage of .95 in the cdf).

899:    Values:
900: +  `VECTAGGER_CDF_GATHER`    - gather results to MPI rank 0, perform the computation and broadcast the result
901: -  `VECTAGGER_CDF_ITERATIVE` - compute the results on all ranks iteratively using `MPI_Allreduce()`

903:   Level: advanced

905:   Note:
906:   Relevant only in parallel: in serial it is directly computed.

908: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerCDFSetMethod()`, `VecTaggerCDFMethods`
909: E*/
910: typedef enum {
911:   VECTAGGER_CDF_GATHER,
912:   VECTAGGER_CDF_ITERATIVE,
913:   VECTAGGER_CDF_NUM_METHODS
914: } VecTaggerCDFMethod;
915: PETSC_EXTERN const char *const VecTaggerCDFMethods[];

917: PETSC_EXTERN PetscErrorCode VecTaggerCDFSetMethod(VecTagger, VecTaggerCDFMethod);
918: PETSC_EXTERN PetscErrorCode VecTaggerCDFGetMethod(VecTagger, VecTaggerCDFMethod *);
919: PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeSetTolerances(VecTagger, PetscInt, PetscReal, PetscReal);
920: PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeGetTolerances(VecTagger, PetscInt *, PetscReal *, PetscReal *);

922: PETSC_EXTERN PetscErrorCode VecTaggerOrSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
923: PETSC_EXTERN PetscErrorCode VecTaggerOrGetSubs(VecTagger, PetscInt *, VecTagger **);

925: PETSC_EXTERN PetscErrorCode VecTaggerAndSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
926: PETSC_EXTERN PetscErrorCode VecTaggerAndGetSubs(VecTagger, PetscInt *, VecTagger **);

928: PETSC_EXTERN PetscErrorCode VecTaggerInitializePackage(void);
929: PETSC_EXTERN PetscErrorCode VecTaggerFinalizePackage(void);

931: #if PetscDefined(USE_DEBUG)
932: /* This is an internal debug-only routine that should not be used by users */
933: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecValidValues_Internal(Vec, PetscInt, PetscBool);
934: #else
935:   #define VecValidValues_Internal(...) PETSC_SUCCESS
936: #endif /* PETSC_USE_DEBUG */

938: #define VEC_CUPM_NOT_CONFIGURED(impl) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP_SYS, "Must configure PETSc with " PetscStringize(impl) " support to use %s", PETSC_FUNCTION_NAME)

940: #if PetscDefined(HAVE_CUDA)
941:   #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
942: #else
943:   #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) \
944:     static inline __decl__ \
945:     { \
946:       __VA_ARGS__; \
947:       VEC_CUPM_NOT_CONFIGURED(cuda); \
948:     }
949: #endif /* PETSC_HAVE_CUDA */

951: /* extra underscore here to make it line up with the cuda versions */
952: #if PetscDefined(HAVE_HIP)
953:   #define VEC_HIP__DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
954: #else
955:   #define VEC_HIP__DECL_OR_STUB(__decl__, ...) \
956:     static inline __decl__ \
957:     { \
958:       __VA_ARGS__; \
959:       VEC_CUPM_NOT_CONFIGURED(hip); \
960:     }
961: #endif /* PETSC_HAVE_HIP */

963: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateSeqCUDA(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)
964: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateSeqHIP(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)

966: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateSeqCUDAWithArray(MPI_Comm a, PetscInt b, PetscInt c, const PetscScalar *d, Vec *e), (void)a, (void)b, (void)c, (void)d, (void)e)
967: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateSeqHIPWithArray(MPI_Comm a, PetscInt b, PetscInt c, const PetscScalar *d, Vec *e), (void)a, (void)b, (void)c, (void)d, (void)e)

969: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateSeqCUDAWithArrays(MPI_Comm a, PetscInt b, PetscInt c, const PetscScalar *d, const PetscScalar *e, Vec *f), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f)
970: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateSeqHIPWithArrays(MPI_Comm a, PetscInt b, PetscInt c, const PetscScalar *d, const PetscScalar *e, Vec *f), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f)

972: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateMPICUDA(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)
973: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateMPIHIP(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)

975: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateMPICUDAWithArray(MPI_Comm a, PetscInt b, PetscInt c, PetscInt d, const PetscScalar *e, Vec *f), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f)
976: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateMPIHIPWithArray(MPI_Comm a, PetscInt b, PetscInt c, PetscInt d, const PetscScalar *e, Vec *f), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f)

978: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateMPICUDAWithArrays(MPI_Comm a, PetscInt b, PetscInt c, PetscInt d, const PetscScalar *e, const PetscScalar *f, Vec *g), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f, (void)g)
979: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateMPIHIPWithArrays(MPI_Comm a, PetscInt b, PetscInt c, PetscInt d, const PetscScalar *e, const PetscScalar *f, Vec *g), (void)a, (void)b, (void)c, (void)d, (void)e, (void)f, (void)g)

981: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArray(Vec a, PetscScalar **b), (void)a, (void)b)
982: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArray(Vec a, PetscScalar **b), (void)a, (void)b)

984: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)
985: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)

987: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
988: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)

990: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
991: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)

993: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
994: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)

996: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
997: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)

999: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
1000: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)

1002: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
1003: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)

1005: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAResetArray(Vec a), (void)a)
1006: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPResetArray(Vec a), (void)a)

1008: #undef VEC_CUPM_NOT_CONFIGURED
1009: #undef VEC_CUDA_DECL_OR_STUB
1010: #undef VEC_HIP__DECL_OR_STUB