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`
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`
272: M*/
274: /*MC
275: REDUCTION_SUM_IMAGINARYPART - sum of imaginary part of matrix column
277: Level: beginner
279: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_SUM_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
280: M*/
282: /*MC
283: REDUCTION_MEAN_REALPART - arithmetic mean of real part of matrix column
285: Level: beginner
287: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_IMAGINARYPART`, `REDUCTION_SUM_REALPART`
288: M*/
290: /*MC
291: REDUCTION_MEAN_IMAGINARYPART - arithmetic mean of imaginary part of matrix column
293: Level: beginner
295: .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_REALPART`, `REDUCTION_SUM_IMAGINARYPART`
296: M*/
298: PETSC_EXTERN PetscErrorCode VecNorm(Vec, NormType, PetscReal *);
299: PETSC_EXTERN PetscErrorCode VecNormAvailable(Vec, NormType, PetscBool *, PetscReal *);
300: PETSC_EXTERN PetscErrorCode VecNormalize(Vec, PetscReal *);
301: PETSC_EXTERN PetscErrorCode VecSum(Vec, PetscScalar *);
302: PETSC_EXTERN PetscErrorCode VecMean(Vec, PetscScalar *);
303: PETSC_EXTERN PetscErrorCode VecMax(Vec, PetscInt *, PetscReal *);
304: PETSC_EXTERN PetscErrorCode VecMin(Vec, PetscInt *, PetscReal *);
305: PETSC_EXTERN PetscErrorCode VecScale(Vec, PetscScalar);
306: PETSC_EXTERN PetscErrorCode VecCopy(Vec, Vec);
307: PETSC_EXTERN PetscErrorCode VecSetRandom(Vec, PetscRandom);
308: PETSC_EXTERN PetscErrorCode VecSet(Vec, PetscScalar);
309: PETSC_EXTERN PetscErrorCode VecSetInf(Vec);
310: PETSC_EXTERN PetscErrorCode VecSwap(Vec, Vec);
311: PETSC_EXTERN PetscErrorCode VecAXPY(Vec, PetscScalar, Vec);
312: PETSC_EXTERN PetscErrorCode VecAXPBY(Vec, PetscScalar, PetscScalar, Vec);
313: PETSC_EXTERN PetscErrorCode VecMAXPY(Vec, PetscInt, const PetscScalar[], Vec[]);
314: PETSC_EXTERN PetscErrorCode VecMAXPBY(Vec, PetscInt, const PetscScalar[], PetscScalar, Vec[]);
315: PETSC_EXTERN PetscErrorCode VecAYPX(Vec, PetscScalar, Vec);
316: PETSC_EXTERN PetscErrorCode VecWAXPY(Vec, PetscScalar, Vec, Vec);
317: PETSC_EXTERN PetscErrorCode VecAXPBYPCZ(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec);
318: PETSC_EXTERN PetscErrorCode VecPointwiseMax(Vec, Vec, Vec);
319: PETSC_EXTERN PetscErrorCode VecPointwiseMaxAbs(Vec, Vec, Vec);
320: PETSC_EXTERN PetscErrorCode VecPointwiseMin(Vec, Vec, Vec);
321: PETSC_EXTERN PetscErrorCode VecPointwiseMult(Vec, Vec, Vec);
322: PETSC_EXTERN PetscErrorCode VecPointwiseDivide(Vec, Vec, Vec);
323: PETSC_EXTERN PetscErrorCode VecMaxPointwiseDivide(Vec, Vec, PetscReal *);
324: PETSC_EXTERN PetscErrorCode VecShift(Vec, PetscScalar);
325: PETSC_EXTERN PetscErrorCode VecReciprocal(Vec);
326: PETSC_EXTERN PetscErrorCode VecPermute(Vec, IS, PetscBool);
327: PETSC_EXTERN PetscErrorCode VecSqrtAbs(Vec);
328: PETSC_EXTERN PetscErrorCode VecLog(Vec);
329: PETSC_EXTERN PetscErrorCode VecExp(Vec);
330: PETSC_EXTERN PetscErrorCode VecAbs(Vec);
331: PETSC_EXTERN PetscErrorCode VecDuplicate(Vec, Vec *);
332: PETSC_EXTERN PetscErrorCode VecDuplicateVecs(Vec, PetscInt, Vec *[]);
333: PETSC_EXTERN PetscErrorCode VecDestroyVecs(PetscInt, Vec *[]);
334: PETSC_EXTERN PetscErrorCode VecStrideNormAll(Vec, NormType, PetscReal[]);
335: PETSC_EXTERN PetscErrorCode VecStrideMaxAll(Vec, PetscInt[], PetscReal[]);
336: PETSC_EXTERN PetscErrorCode VecStrideMinAll(Vec, PetscInt[], PetscReal[]);
337: PETSC_EXTERN PetscErrorCode VecStrideScaleAll(Vec, const PetscScalar[]);
338: PETSC_EXTERN PetscErrorCode VecStrideSumAll(Vec, PetscScalar *);
339: PETSC_EXTERN PetscErrorCode VecUniqueEntries(Vec, PetscInt *, PetscScalar **);
341: PETSC_EXTERN PetscErrorCode VecStrideNorm(Vec, PetscInt, NormType, PetscReal *);
342: PETSC_EXTERN PetscErrorCode VecStrideMax(Vec, PetscInt, PetscInt *, PetscReal *);
343: PETSC_EXTERN PetscErrorCode VecStrideMin(Vec, PetscInt, PetscInt *, PetscReal *);
344: PETSC_EXTERN PetscErrorCode VecStrideScale(Vec, PetscInt, PetscScalar);
345: PETSC_EXTERN PetscErrorCode VecStrideSum(Vec, PetscInt, PetscScalar *);
346: PETSC_EXTERN PetscErrorCode VecStrideSet(Vec, PetscInt, PetscScalar);
348: PETSC_EXTERN PetscErrorCode VecStrideGather(Vec, PetscInt, Vec, InsertMode);
349: PETSC_EXTERN PetscErrorCode VecStrideScatter(Vec, PetscInt, Vec, InsertMode);
350: PETSC_EXTERN PetscErrorCode VecStrideGatherAll(Vec, Vec[], InsertMode);
351: PETSC_EXTERN PetscErrorCode VecStrideScatterAll(Vec[], Vec, InsertMode);
353: PETSC_EXTERN PetscErrorCode VecStrideSubSetScatter(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
354: PETSC_EXTERN PetscErrorCode VecStrideSubSetGather(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
356: PETSC_EXTERN PetscErrorCode VecSetValues(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
357: PETSC_EXTERN PetscErrorCode VecGetValues(Vec, PetscInt, const PetscInt[], PetscScalar[]);
358: PETSC_EXTERN PetscErrorCode VecAssemblyBegin(Vec);
359: PETSC_EXTERN PetscErrorCode VecAssemblyEnd(Vec);
360: PETSC_EXTERN PetscErrorCode VecStashSetInitialSize(Vec, PetscInt, PetscInt);
361: PETSC_EXTERN PetscErrorCode VecStashView(Vec, PetscViewer);
362: PETSC_EXTERN PetscErrorCode VecStashViewFromOptions(Vec, PetscObject, const char[]);
363: PETSC_EXTERN PetscErrorCode VecStashGetInfo(Vec, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
365: PETSC_EXTERN PetscErrorCode VecSetPreallocationCOO(Vec, PetscCount, const PetscInt[]);
366: PETSC_EXTERN PetscErrorCode VecSetPreallocationCOOLocal(Vec, PetscCount, PetscInt[]);
367: PETSC_EXTERN PetscErrorCode VecSetValuesCOO(Vec, const PetscScalar[], InsertMode);
369: /*MC
370: VecSetValue - Set a single entry into a vector.
372: Synopsis:
373: #include <petscvec.h>
374: PetscErrorCode VecSetValue(Vec v,PetscInt row,PetscScalar value, InsertMode mode);
376: Not Collective
378: Input Parameters:
379: + v - the vector
380: . row - the row location of the entry
381: . value - the value to insert
382: - mode - either `INSERT_VALUES` or `ADD_VALUES`
384: Level: beginner
386: Notes:
387: For efficiency one should use `VecSetValues()` and set several or
388: many values simultaneously if possible.
390: These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
391: MUST be called after all calls to `VecSetValue()` have been completed.
393: `VecSetValue()` uses 0-based indices in Fortran as well as in C.
395: .seealso: [](ch_vectors), `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValueLocal()`
396: M*/
397: static inline PetscErrorCode VecSetValue(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
398: {
399: return VecSetValues(v, 1, &i, &va, mode);
400: }
402: PETSC_EXTERN PetscErrorCode VecSetBlockSize(Vec, PetscInt);
403: PETSC_EXTERN PetscErrorCode VecGetBlockSize(Vec, PetscInt *);
404: PETSC_EXTERN PetscErrorCode VecSetValuesBlocked(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
406: /* Dynamic creation and loading functions */
407: PETSC_EXTERN PetscFunctionList VecList;
408: PETSC_EXTERN PetscErrorCode VecSetType(Vec, VecType);
409: PETSC_EXTERN PetscErrorCode VecGetType(Vec, VecType *);
410: PETSC_EXTERN PetscErrorCode VecRegister(const char[], PetscErrorCode (*)(Vec));
411: PETSC_EXTERN PetscErrorCode VecRegisterAll(void);
413: PETSC_EXTERN PetscErrorCode VecScatterBegin(VecScatter, Vec, Vec, InsertMode, ScatterMode);
414: PETSC_EXTERN PetscErrorCode VecScatterEnd(VecScatter, Vec, Vec, InsertMode, ScatterMode);
415: PETSC_EXTERN PetscErrorCode VecScatterDestroy(VecScatter *);
416: PETSC_EXTERN PetscErrorCode VecScatterSetUp(VecScatter);
417: PETSC_EXTERN PetscErrorCode VecScatterCopy(VecScatter, VecScatter *);
418: PETSC_EXTERN PetscErrorCode VecScatterView(VecScatter, PetscViewer);
419: PETSC_EXTERN PetscErrorCode VecScatterViewFromOptions(VecScatter, PetscObject, const char[]);
420: PETSC_EXTERN PetscErrorCode VecScatterRemap(VecScatter, PetscInt[], PetscInt[]);
421: PETSC_EXTERN PetscErrorCode VecScatterGetMerged(VecScatter, PetscBool *);
423: PETSC_EXTERN PetscErrorCode VecGetArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
424: PETSC_EXTERN PetscErrorCode VecRestoreArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
425: PETSC_EXTERN PetscErrorCode VecGetArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
426: PETSC_EXTERN PetscErrorCode VecRestoreArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
427: PETSC_EXTERN PetscErrorCode VecGetArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
428: PETSC_EXTERN PetscErrorCode VecRestoreArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
429: PETSC_EXTERN PetscErrorCode VecGetArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);
430: PETSC_EXTERN PetscErrorCode VecRestoreArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);
432: PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
433: PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
434: PETSC_EXTERN PetscErrorCode VecRestoreArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
435: PETSC_EXTERN PetscErrorCode VecGetArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
436: PETSC_EXTERN PetscErrorCode VecRestoreArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
437: PETSC_EXTERN PetscErrorCode VecGetArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
438: PETSC_EXTERN PetscErrorCode VecRestoreArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
439: PETSC_EXTERN PetscErrorCode VecGetArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);
440: PETSC_EXTERN PetscErrorCode VecRestoreArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);
442: PETSC_EXTERN PetscErrorCode VecGetArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
443: PETSC_EXTERN PetscErrorCode VecRestoreArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
444: PETSC_EXTERN PetscErrorCode VecGetArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
445: PETSC_EXTERN PetscErrorCode VecRestoreArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
446: PETSC_EXTERN PetscErrorCode VecGetArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
447: PETSC_EXTERN PetscErrorCode VecRestoreArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
448: PETSC_EXTERN PetscErrorCode VecGetArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);
449: PETSC_EXTERN PetscErrorCode VecRestoreArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);
451: PETSC_EXTERN PetscErrorCode VecPlaceArray(Vec, const PetscScalar[]);
452: PETSC_EXTERN PetscErrorCode VecResetArray(Vec);
453: PETSC_EXTERN PetscErrorCode VecReplaceArray(Vec, const PetscScalar[]);
455: PETSC_EXTERN PetscErrorCode VecGetArrays(const Vec[], PetscInt, PetscScalar **[]);
456: PETSC_EXTERN PetscErrorCode VecRestoreArrays(const Vec[], PetscInt, PetscScalar **[]);
458: PETSC_EXTERN PetscErrorCode VecView(Vec, PetscViewer);
459: PETSC_EXTERN PetscErrorCode VecViewNative(Vec, PetscViewer);
460: PETSC_EXTERN PetscErrorCode VecEqual(Vec, Vec, PetscBool *);
461: PETSC_EXTERN PetscErrorCode VecLoad(Vec, PetscViewer);
463: PETSC_EXTERN PetscErrorCode VecGetSize(Vec, PetscInt *);
464: PETSC_EXTERN PetscErrorCode VecGetLocalSize(Vec, PetscInt *);
465: PETSC_EXTERN PetscErrorCode VecGetOwnershipRange(Vec, PetscInt *, PetscInt *);
466: PETSC_EXTERN PetscErrorCode VecGetOwnershipRanges(Vec, const PetscInt *[]);
468: PETSC_EXTERN PetscErrorCode VecSetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping);
469: PETSC_EXTERN PetscErrorCode VecSetValuesLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
471: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLContext(Vec, PETSC_UINTPTR_T *);
472: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLQueue(Vec, PETSC_UINTPTR_T *);
473: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemRead(Vec, PETSC_UINTPTR_T *);
474: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemWrite(Vec, PETSC_UINTPTR_T *);
475: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMemWrite(Vec);
476: PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMem(Vec, PETSC_UINTPTR_T *);
477: PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMem(Vec);
479: /*MC
480: VecSetValueLocal - Set a single entry into a vector using the local numbering, see `VecSetValuesLocal()`
482: Synopsis:
483: #include <petscvec.h>
484: PetscErrorCode VecSetValueLocal(Vec v,PetscInt row,PetscScalar value, InsertMode mode);
486: Not Collective
488: Input Parameters:
489: + v - the vector
490: . row - the row location of the entry
491: . value - the value to insert
492: - mode - either `INSERT_VALUES` or `ADD_VALUES`
494: Level: beginner
496: Notes:
497: For efficiency one should use `VecSetValuesLocal()` and set several or
498: many values simultaneously if possible.
500: These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
501: MUST be called after all calls to `VecSetValueLocal()` have been completed.
503: `VecSetValueLocal()` uses 0-based indices in Fortran as well as in C.
505: .seealso: [](ch_vectors), `VecSetValuesLocal()`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValue()`
506: M*/
507: static inline PetscErrorCode VecSetValueLocal(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
508: {
509: return VecSetValuesLocal(v, 1, &i, &va, mode);
510: }
512: /*MC
513: VecCheckAssembled - checks if values have been changed in the vector, by `VecSetValues()` or related routines, but it has not been assembled
515: Synopsis:
516: #include <petscvec.h>
517: VecCheckAssembled(Vec v);
519: Not Collective
521: Input Parameter:
522: . v - the vector to check
524: Level: developer
526: Note:
527: After calls to `VecSetValues()` and related routines on must call ``VecAssemblyBegin()` and `VecAssemblyEnd()` before using the vector
529: .seealso: [](ch_vectors), `Vec`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`
530: M*/
531: #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()?")
533: PETSC_EXTERN PetscErrorCode VecSetValuesBlockedLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
534: PETSC_EXTERN PetscErrorCode VecGetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping *);
536: PETSC_EXTERN PetscErrorCode VecDotBegin(Vec, Vec, PetscScalar *);
537: PETSC_EXTERN PetscErrorCode VecDotEnd(Vec, Vec, PetscScalar *);
538: PETSC_EXTERN PetscErrorCode VecTDotBegin(Vec, Vec, PetscScalar *);
539: PETSC_EXTERN PetscErrorCode VecTDotEnd(Vec, Vec, PetscScalar *);
540: PETSC_EXTERN PetscErrorCode VecNormBegin(Vec, NormType, PetscReal *);
541: PETSC_EXTERN PetscErrorCode VecNormEnd(Vec, NormType, PetscReal *);
542: PETSC_EXTERN PetscErrorCode VecErrorWeightedNorms(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *);
544: PETSC_EXTERN PetscErrorCode VecMDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
545: PETSC_EXTERN PetscErrorCode VecMDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
546: PETSC_EXTERN PetscErrorCode VecMTDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
547: PETSC_EXTERN PetscErrorCode VecMTDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
548: PETSC_EXTERN PetscErrorCode PetscCommSplitReductionBegin(MPI_Comm);
550: PETSC_EXTERN PetscErrorCode VecBindToCPU(Vec, PetscBool);
551: PETSC_DEPRECATED_FUNCTION(3, 13, 0, "VecBindToCPU()", ) static inline PetscErrorCode VecPinToCPU(Vec v, PetscBool flg)
552: {
553: return VecBindToCPU(v, flg);
554: }
555: PETSC_EXTERN PetscErrorCode VecBoundToCPU(Vec, PetscBool *);
556: PETSC_EXTERN PetscErrorCode VecSetBindingPropagates(Vec, PetscBool);
557: PETSC_EXTERN PetscErrorCode VecGetBindingPropagates(Vec, PetscBool *);
558: PETSC_EXTERN PetscErrorCode VecSetPinnedMemoryMin(Vec, size_t);
559: PETSC_EXTERN PetscErrorCode VecGetPinnedMemoryMin(Vec, size_t *);
561: PETSC_EXTERN PetscErrorCode VecGetOffloadMask(Vec, PetscOffloadMask *);
563: typedef enum {
564: VEC_IGNORE_OFF_PROC_ENTRIES,
565: VEC_IGNORE_NEGATIVE_INDICES,
566: VEC_SUBSET_OFF_PROC_ENTRIES
567: } VecOption;
568: PETSC_EXTERN PetscErrorCode VecSetOption(Vec, VecOption, PetscBool);
570: PETSC_EXTERN PetscErrorCode VecGetArray(Vec, PetscScalar **);
571: PETSC_EXTERN PetscErrorCode VecGetArrayWrite(Vec, PetscScalar **);
572: PETSC_EXTERN PetscErrorCode VecGetArrayRead(Vec, const PetscScalar **);
573: PETSC_EXTERN PetscErrorCode VecRestoreArray(Vec, PetscScalar **);
574: PETSC_EXTERN PetscErrorCode VecRestoreArrayWrite(Vec, PetscScalar **);
575: PETSC_EXTERN PetscErrorCode VecRestoreArrayRead(Vec, const PetscScalar **);
576: PETSC_EXTERN PetscErrorCode VecCreateLocalVector(Vec, Vec *);
577: PETSC_EXTERN PetscErrorCode VecGetLocalVector(Vec, Vec);
578: PETSC_EXTERN PetscErrorCode VecRestoreLocalVector(Vec, Vec);
579: PETSC_EXTERN PetscErrorCode VecGetLocalVectorRead(Vec, Vec);
580: PETSC_EXTERN PetscErrorCode VecRestoreLocalVectorRead(Vec, Vec);
581: PETSC_EXTERN PetscErrorCode VecGetArrayAndMemType(Vec, PetscScalar **, PetscMemType *);
582: PETSC_EXTERN PetscErrorCode VecRestoreArrayAndMemType(Vec, PetscScalar **);
583: PETSC_EXTERN PetscErrorCode VecGetArrayReadAndMemType(Vec, const PetscScalar **, PetscMemType *);
584: PETSC_EXTERN PetscErrorCode VecRestoreArrayReadAndMemType(Vec, const PetscScalar **);
585: PETSC_EXTERN PetscErrorCode VecGetArrayWriteAndMemType(Vec, PetscScalar **, PetscMemType *);
586: PETSC_EXTERN PetscErrorCode VecRestoreArrayWriteAndMemType(Vec, PetscScalar **);
588: /*@C
589: 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
591: Logically Collective; No Fortran Support
593: Input Parameters:
594: + x - the vector
595: - y - the second vector
597: Output Parameters:
598: + xv - location to put pointer to the first array
599: - yv - location to put pointer to the second array
601: Level: developer
603: .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecRestoreArrayPair()`
604: @*/
605: static inline PetscErrorCode VecGetArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
606: {
607: PetscFunctionBegin;
608: PetscCall(VecGetArray(y, yv));
609: if (x == y) *xv = *yv;
610: else PetscCall(VecGetArrayRead(x, (const PetscScalar **)xv));
611: PetscFunctionReturn(PETSC_SUCCESS);
612: }
614: /*@C
615: VecRestoreArrayPair - Returns a pair of pointers for two vectors that may be common obtained with `VecGetArrayPair()`
617: Logically Collective; No Fortran Support
619: Input Parameters:
620: + x - the vector
621: . y - the second vector
622: . xv - location to put pointer to the first array
623: - yv - location to put pointer to the second array
625: Level: developer
627: .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecGetArrayPair()`
628: @*/
629: static inline PetscErrorCode VecRestoreArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
630: {
631: PetscFunctionBegin;
632: PetscCall(VecRestoreArray(y, yv));
633: if (x != y) PetscCall(VecRestoreArrayRead(x, (const PetscScalar **)xv));
634: PetscFunctionReturn(PETSC_SUCCESS);
635: }
637: #if PetscDefined(USE_DEBUG)
638: PETSC_EXTERN PetscErrorCode VecLockReadPush(Vec);
639: PETSC_EXTERN PetscErrorCode VecLockReadPop(Vec);
640: PETSC_EXTERN PetscErrorCode VecLockWriteSet(Vec, PetscBool);
641: PETSC_EXTERN PetscErrorCode VecLockGet(Vec, PetscInt *);
642: PETSC_EXTERN PetscErrorCode VecLockGetLocation(Vec, const char *[], const char *[], int *);
643: static inline PetscErrorCode VecSetErrorIfLocked(Vec x, PetscInt arg)
644: {
645: PetscInt state;
647: PetscFunctionBegin;
648: PetscCall(VecLockGet(x, &state));
649: if (PetscUnlikely(state != 0)) {
650: const char *file, *func, *name;
651: int line;
653: PetscCall(VecLockGetLocation(x, &file, &func, &line));
654: PetscCall(PetscObjectGetName((PetscObject)x, &name));
655: 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);
656: }
657: PetscFunctionReturn(PETSC_SUCCESS);
658: }
659: /* The three are deprecated */
660: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPush()", ) PetscErrorCode VecLockPush(Vec);
661: PETSC_EXTERN PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPop()", ) PetscErrorCode VecLockPop(Vec);
662: #define VecLocked(x, arg) VecSetErrorIfLocked(x, arg) PETSC_DEPRECATED_MACRO(3, 11, 0, "VecSetErrorIfLocked()", )
663: #else
664: #define VecLockReadPush(x) PETSC_SUCCESS
665: #define VecLockReadPop(x) PETSC_SUCCESS
666: #define VecLockGet(x, s) (*(s) = 0, PETSC_SUCCESS)
667: #define VecSetErrorIfLocked(x, arg) PETSC_SUCCESS
668: #define VecLockWriteSet(x, flg) PETSC_SUCCESS
669: /* The three are deprecated */
670: #define VecLockPush(x) PETSC_SUCCESS
671: #define VecLockPop(x) PETSC_SUCCESS
672: #define VecLocked(x, arg) PETSC_SUCCESS
673: #endif
675: /*E
676: VecOperation - Enumeration of overide-able methods in the `Vec` implementation function-table.
678: Values:
679: + `VECOP_DUPLICATE` - `VecDuplicate()`
680: . `VECOP_SET` - `VecSet()`
681: . `VECOP_VIEW` - `VecView()`
682: . `VECOP_LOAD` - `VecLoad()`
683: . `VECOP_VIEWNATIVE` - `VecViewNative()`
684: - `VECOP_LOADNATIVE` - `VecLoadNative()`
686: Level: advanced
688: Notes:
689: Some operations may serve as the implementation for other routines not listed above. For
690: example `VECOP_SET` can be used to simultaneously overriding the implementation used in
691: `VecSet()`, `VecSetInf()`, and `VecZeroEntries()`.
693: Entries to `VecOperation` are added as needed so if you do not see the operation listed which
694: you'd like to replace, please send mail to `petsc-maint@mcs.anl.gov`!
696: .seealso: [](ch_vectors), `Vec`, `VecSetOperation()`
697: E*/
698: typedef enum {
699: VECOP_DUPLICATE = 0,
700: VECOP_SET = 10,
701: VECOP_VIEW = 33,
702: VECOP_LOAD = 41,
703: VECOP_VIEWNATIVE = 69,
704: VECOP_LOADNATIVE = 70
705: } VecOperation;
706: PETSC_EXTERN PetscErrorCode VecSetOperation(Vec, VecOperation, void (*)(void));
708: /*
709: Routines for dealing with ghosted vectors:
710: vectors with ghost elements at the end of the array.
711: */
712: PETSC_EXTERN PetscErrorCode VecMPISetGhost(Vec, PetscInt, const PetscInt[]);
713: PETSC_EXTERN PetscErrorCode VecCreateGhost(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
714: PETSC_EXTERN PetscErrorCode VecCreateGhostWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
715: PETSC_EXTERN PetscErrorCode VecCreateGhostBlock(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
716: PETSC_EXTERN PetscErrorCode VecCreateGhostBlockWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
717: PETSC_EXTERN PetscErrorCode VecGhostGetGhostIS(Vec, IS *);
718: PETSC_EXTERN PetscErrorCode VecGhostGetLocalForm(Vec, Vec *);
719: PETSC_EXTERN PetscErrorCode VecGhostRestoreLocalForm(Vec, Vec *);
720: PETSC_EXTERN PetscErrorCode VecGhostIsLocalForm(Vec, Vec, PetscBool *);
721: PETSC_EXTERN PetscErrorCode VecGhostUpdateBegin(Vec, InsertMode, ScatterMode);
722: PETSC_EXTERN PetscErrorCode VecGhostUpdateEnd(Vec, InsertMode, ScatterMode);
724: PETSC_EXTERN PetscErrorCode VecConjugate(Vec);
725: PETSC_EXTERN PetscErrorCode VecImaginaryPart(Vec);
726: PETSC_EXTERN PetscErrorCode VecRealPart(Vec);
728: PETSC_EXTERN PetscErrorCode VecScatterCreateToAll(Vec, VecScatter *, Vec *);
729: PETSC_EXTERN PetscErrorCode VecScatterCreateToZero(Vec, VecScatter *, Vec *);
731: /* These vector calls were included from TAO. They miss vtable entries and GPU implementation */
732: PETSC_EXTERN PetscErrorCode ISComplementVec(IS, Vec, IS *);
733: PETSC_EXTERN PetscErrorCode VecPow(Vec, PetscScalar);
734: PETSC_EXTERN PetscErrorCode VecMedian(Vec, Vec, Vec, Vec);
735: PETSC_EXTERN PetscErrorCode VecWhichInactive(Vec, Vec, Vec, Vec, PetscBool, IS *);
736: PETSC_EXTERN PetscErrorCode VecWhichBetween(Vec, Vec, Vec, IS *);
737: PETSC_EXTERN PetscErrorCode VecWhichBetweenOrEqual(Vec, Vec, Vec, IS *);
738: PETSC_EXTERN PetscErrorCode VecWhichGreaterThan(Vec, Vec, IS *);
739: PETSC_EXTERN PetscErrorCode VecWhichLessThan(Vec, Vec, IS *);
740: PETSC_EXTERN PetscErrorCode VecWhichEqual(Vec, Vec, IS *);
741: PETSC_EXTERN PetscErrorCode VecISAXPY(Vec, IS, PetscScalar, Vec);
742: PETSC_EXTERN PetscErrorCode VecISCopy(Vec, IS, ScatterMode, Vec);
743: PETSC_EXTERN PetscErrorCode VecISSet(Vec, IS, PetscScalar);
744: PETSC_EXTERN PetscErrorCode VecISShift(Vec, IS, PetscScalar);
745: PETSC_EXTERN PetscErrorCode VecBoundGradientProjection(Vec, Vec, Vec, Vec, Vec);
746: PETSC_EXTERN PetscErrorCode VecStepBoundInfo(Vec, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
747: PETSC_EXTERN PetscErrorCode VecStepMax(Vec, Vec, PetscReal *);
748: PETSC_EXTERN PetscErrorCode VecStepMaxBounded(Vec, Vec, Vec, Vec, PetscReal *);
750: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaGetVector(PetscViewer, Vec);
751: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutVector(PetscViewer, Vec);
753: /*S
754: Vecs - Collection of vectors where the data for the vectors is stored in
755: one contiguous memory
757: Level: advanced
759: Notes:
760: Temporary construct for handling multiply right-hand side solves
762: This is faked by storing a single vector that has enough array space for
763: n vectors
765: S*/
766: struct _n_Vecs {
767: PetscInt n;
768: Vec v;
769: };
770: typedef struct _n_Vecs *Vecs;
771: PETSC_EXTERN PetscErrorCode VecsDestroy(Vecs);
772: PETSC_EXTERN PetscErrorCode VecsCreateSeq(MPI_Comm, PetscInt, PetscInt, Vecs *);
773: PETSC_EXTERN PetscErrorCode VecsCreateSeqWithArray(MPI_Comm, PetscInt, PetscInt, PetscScalar *, Vecs *);
774: PETSC_EXTERN PetscErrorCode VecsDuplicate(Vecs, Vecs *);
776: #if PetscDefined(HAVE_VIENNACL)
777: typedef struct _p_PetscViennaCLIndices *PetscViennaCLIndices;
778: PETSC_EXTERN PetscErrorCode VecViennaCLCopyToGPUSome_Public(Vec, PetscViennaCLIndices);
779: PETSC_EXTERN PetscErrorCode VecViennaCLCopyFromGPUSome_Public(Vec, PetscViennaCLIndices);
780: PETSC_EXTERN PetscErrorCode VecCreateSeqViennaCL(MPI_Comm, PetscInt, Vec *);
781: PETSC_EXTERN PetscErrorCode VecCreateMPIViennaCL(MPI_Comm, PetscInt, PetscInt, Vec *);
782: #endif
783: #if PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP)
784: PETSC_EXTERN PetscErrorCode VecScatterInitializeForGPU(VecScatter, Vec);
785: PETSC_EXTERN PetscErrorCode VecScatterFinalizeForGPU(VecScatter);
786: #endif
787: #if PetscDefined(HAVE_KOKKOS_KERNELS)
788: PETSC_EXTERN PetscErrorCode VecCreateSeqKokkos(MPI_Comm, PetscInt, Vec *);
789: PETSC_EXTERN PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm, PetscInt, PetscInt, const PetscScalar *, Vec *);
790: PETSC_EXTERN PetscErrorCode VecCreateMPIKokkos(MPI_Comm, PetscInt, PetscInt, Vec *);
791: PETSC_EXTERN PetscErrorCode VecCreateMPIKokkosWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar *, Vec *);
792: #endif
794: PETSC_EXTERN PetscErrorCode VecNestGetSubVecs(Vec, PetscInt *, Vec **);
795: PETSC_EXTERN PetscErrorCode VecNestGetSubVec(Vec, PetscInt, Vec *);
796: PETSC_EXTERN PetscErrorCode VecNestSetSubVecs(Vec, PetscInt, PetscInt *, Vec *);
797: PETSC_EXTERN PetscErrorCode VecNestSetSubVec(Vec, PetscInt, Vec);
798: PETSC_EXTERN PetscErrorCode VecCreateNest(MPI_Comm, PetscInt, IS *, Vec *, Vec *);
799: PETSC_EXTERN PetscErrorCode VecNestGetSize(Vec, PetscInt *);
801: PETSC_EXTERN PetscErrorCode PetscOptionsGetVec(PetscOptions, const char[], const char[], Vec, PetscBool *);
802: PETSC_EXTERN PetscErrorCode VecFilter(Vec, PetscReal);
803: PETSC_DEPRECATED_FUNCTION(3, 20, 0, "VecFilter()", ) static inline PetscErrorCode VecChop(Vec v, PetscReal tol)
804: {
805: return VecFilter(v, tol);
806: }
808: PETSC_EXTERN PetscErrorCode VecGetLayout(Vec, PetscLayout *);
809: PETSC_EXTERN PetscErrorCode VecSetLayout(Vec, PetscLayout);
811: PETSC_EXTERN PetscErrorCode PetscSectionVecView(PetscSection, Vec, PetscViewer);
812: PETSC_EXTERN PetscErrorCode VecGetValuesSection(Vec, PetscSection, PetscInt, PetscScalar **);
813: PETSC_EXTERN PetscErrorCode VecSetValuesSection(Vec, PetscSection, PetscInt, const PetscScalar[], InsertMode);
814: PETSC_EXTERN PetscErrorCode PetscSectionVecNorm(PetscSection, PetscSection, Vec, NormType, PetscReal[]);
816: /*S
817: VecTagger - Object used to manage the tagging of a subset of indices based on the values of a vector. The
818: motivating application is the selection of cells for refinement or coarsening based on vector containing
819: the values in an error indicator metric.
821: Values:
822: + `VECTAGGERABSOLUTE` - "absolute" values are in a interval (box for complex values) of explicitly defined values
823: . `VECTAGGERRELATIVE` - "relative" values are in a interval (box for complex values) of values relative to the set of all values in the vector
824: . `VECTAGGERCDF` - "cdf" values are in a relative range of the *cumulative distribution* of values in the vector
825: . `VECTAGGEROR` - "or" values are in the union of other tags
826: - `VECTAGGERAND` - "and" values are in the intersection of other tags
828: Level: advanced
830: Developer Note:
831: Why not use a `DMLabel` or similar object
833: .seealso: [](ch_vectors), `Vec`, `VecTaggerType`, `VecTaggerCreate()`
834: S*/
835: typedef struct _p_VecTagger *VecTagger;
837: /*J
838: VecTaggerType - String with the name of a `VecTagger` type
840: Level: advanced
842: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerCreate()`
843: J*/
844: typedef const char *VecTaggerType;
845: #define VECTAGGERABSOLUTE "absolute"
846: #define VECTAGGERRELATIVE "relative"
847: #define VECTAGGERCDF "cdf"
848: #define VECTAGGEROR "or"
849: #define VECTAGGERAND "and"
851: PETSC_EXTERN PetscClassId VEC_TAGGER_CLASSID;
852: PETSC_EXTERN PetscFunctionList VecTaggerList;
853: PETSC_EXTERN PetscErrorCode VecTaggerRegister(const char[], PetscErrorCode (*)(VecTagger));
854: PETSC_EXTERN PetscErrorCode VecTaggerRegisterAll(void);
856: PETSC_EXTERN PetscErrorCode VecTaggerCreate(MPI_Comm, VecTagger *);
857: PETSC_EXTERN PetscErrorCode VecTaggerSetBlockSize(VecTagger, PetscInt);
858: PETSC_EXTERN PetscErrorCode VecTaggerGetBlockSize(VecTagger, PetscInt *);
859: PETSC_EXTERN PetscErrorCode VecTaggerSetType(VecTagger, VecTaggerType);
860: PETSC_EXTERN PetscErrorCode VecTaggerGetType(VecTagger, VecTaggerType *);
861: PETSC_EXTERN PetscErrorCode VecTaggerSetInvert(VecTagger, PetscBool);
862: PETSC_EXTERN PetscErrorCode VecTaggerGetInvert(VecTagger, PetscBool *);
863: PETSC_EXTERN PetscErrorCode VecTaggerSetFromOptions(VecTagger);
864: PETSC_EXTERN PetscErrorCode VecTaggerSetUp(VecTagger);
865: PETSC_EXTERN PetscErrorCode VecTaggerView(VecTagger, PetscViewer);
866: PETSC_EXTERN PetscErrorCode VecTaggerComputeIS(VecTagger, Vec, IS *, PetscBool *);
867: PETSC_EXTERN PetscErrorCode VecTaggerDestroy(VecTagger *);
869: /*S
870: VecTaggerBox - A interval (box for complex numbers) range used to tag values. For real scalars, this is just a closed interval; for complex scalars,
871: 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.
873: Level: beginner
875: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerComputeIntervals()`
876: S*/
877: typedef struct {
878: PetscScalar min;
879: PetscScalar max;
880: } VecTaggerBox;
881: PETSC_EXTERN PetscErrorCode VecTaggerComputeBoxes(VecTagger, Vec, PetscInt *, VecTaggerBox **, PetscBool *);
883: PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteSetBox(VecTagger, VecTaggerBox *);
884: PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteGetBox(VecTagger, const VecTaggerBox **);
886: PETSC_EXTERN PetscErrorCode VecTaggerRelativeSetBox(VecTagger, VecTaggerBox *);
887: PETSC_EXTERN PetscErrorCode VecTaggerRelativeGetBox(VecTagger, const VecTaggerBox **);
889: PETSC_EXTERN PetscErrorCode VecTaggerCDFSetBox(VecTagger, VecTaggerBox *);
890: PETSC_EXTERN PetscErrorCode VecTaggerCDFGetBox(VecTagger, const VecTaggerBox **);
892: /*E
893: 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).
895: Values:
896: + `VECTAGGER_CDF_GATHER` - gather results to MPI rank 0, perform the computation and broadcast the result
897: - `VECTAGGER_CDF_ITERATIVE` - compute the results on all ranks iteratively using `MPI_Allreduce()`
899: Level: advanced
901: Note:
902: Relevant only in parallel: in serial it is directly computed.
904: .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerCDFSetMethod()`, `VecTaggerCDFMethods`
905: E*/
906: typedef enum {
907: VECTAGGER_CDF_GATHER,
908: VECTAGGER_CDF_ITERATIVE,
909: VECTAGGER_CDF_NUM_METHODS
910: } VecTaggerCDFMethod;
911: PETSC_EXTERN const char *const VecTaggerCDFMethods[];
913: PETSC_EXTERN PetscErrorCode VecTaggerCDFSetMethod(VecTagger, VecTaggerCDFMethod);
914: PETSC_EXTERN PetscErrorCode VecTaggerCDFGetMethod(VecTagger, VecTaggerCDFMethod *);
915: PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeSetTolerances(VecTagger, PetscInt, PetscReal, PetscReal);
916: PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeGetTolerances(VecTagger, PetscInt *, PetscReal *, PetscReal *);
918: PETSC_EXTERN PetscErrorCode VecTaggerOrSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
919: PETSC_EXTERN PetscErrorCode VecTaggerOrGetSubs(VecTagger, PetscInt *, VecTagger **);
921: PETSC_EXTERN PetscErrorCode VecTaggerAndSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
922: PETSC_EXTERN PetscErrorCode VecTaggerAndGetSubs(VecTagger, PetscInt *, VecTagger **);
924: PETSC_EXTERN PetscErrorCode VecTaggerInitializePackage(void);
925: PETSC_EXTERN PetscErrorCode VecTaggerFinalizePackage(void);
927: #if PetscDefined(USE_DEBUG)
928: /* This is an internal debug-only routine that should not be used by users */
929: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecValidValues_Internal(Vec, PetscInt, PetscBool);
930: #else
931: #define VecValidValues_Internal(...) PETSC_SUCCESS
932: #endif /* PETSC_USE_DEBUG */
934: #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)
936: #if PetscDefined(HAVE_CUDA)
937: #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
938: #else
939: #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) \
940: static inline __decl__ \
941: { \
942: __VA_ARGS__; \
943: VEC_CUPM_NOT_CONFIGURED(cuda); \
944: }
945: #endif /* PETSC_HAVE_CUDA */
947: /* extra underscore here to make it line up with the cuda versions */
948: #if PetscDefined(HAVE_HIP)
949: #define VEC_HIP__DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
950: #else
951: #define VEC_HIP__DECL_OR_STUB(__decl__, ...) \
952: static inline __decl__ \
953: { \
954: __VA_ARGS__; \
955: VEC_CUPM_NOT_CONFIGURED(hip); \
956: }
957: #endif /* PETSC_HAVE_HIP */
959: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateSeqCUDA(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)
960: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateSeqHIP(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)
962: 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)
963: 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)
965: 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)
966: 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)
968: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateMPICUDA(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)
969: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateMPIHIP(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)
971: 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)
972: 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)
974: 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)
975: 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)
977: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArray(Vec a, PetscScalar **b), (void)a, (void)b)
978: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArray(Vec a, PetscScalar **b), (void)a, (void)b)
980: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)
981: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)
983: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
984: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
986: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
987: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
989: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
990: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
992: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
993: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
995: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
996: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
998: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
999: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
1001: VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAResetArray(Vec a), (void)a)
1002: VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPResetArray(Vec a), (void)a)
1004: #undef VEC_CUPM_NOT_CONFIGURED
1005: #undef VEC_CUDA_DECL_OR_STUB
1006: #undef VEC_HIP__DECL_OR_STUB