Line data Source code
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
7 :
8 : #include <petscsys.h>
9 : #include <petscsftypes.h> /* for VecScatter, VecScatterType */
10 : #include <petscis.h>
11 : #include <petscdevicetypes.h>
12 : #include <petscviewer.h>
13 :
14 : /* SUBMANSEC = Vec */
15 :
16 : /*S
17 : Vec - Abstract PETSc vector object. Used for holding solutions and right-hand sides for linear systems, nonlinear systems, and time integrators
18 :
19 : Level: beginner
20 :
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
23 :
24 : .seealso: [](doc_vector), [](ch_vectors), `VecCreate()`, `VecType`, `VecSetType()`
25 : S*/
26 : typedef struct _p_Vec *Vec;
27 :
28 : /*E
29 : ScatterMode - Determines the direction of a scatter in `VecScatterBegin()` and `VecScatterEnd()`
30 :
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
37 :
38 : Level: beginner
39 :
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;
48 :
49 : /*MC
50 : SCATTER_FORWARD - Scatters the values as dictated by the `VecScatterCreate()` call during `VecScatterBegin()` and `VecScatterEnd()`
51 :
52 : Level: beginner
53 :
54 : .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_REVERSE`, `SCATTER_FORWARD_LOCAL`,
55 : `SCATTER_REVERSE_LOCAL`
56 : M*/
57 :
58 : /*MC
59 : SCATTER_REVERSE - Moves the values in the opposite direction then the directions indicated
60 : in the `VecScatterCreate()` during `VecScatterBegin()` and `VecScatterEnd()`
61 :
62 : Level: beginner
63 :
64 : .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_FORWARD`, `SCATTER_FORWARD_LOCAL`,
65 : `SCATTER_REVERSE_LOCAL`
66 : M*/
67 :
68 : /*MC
69 : SCATTER_FORWARD_LOCAL - Scatters the values as dictated by the `VecScatterCreate()` during `VecScatterBegin()` and `VecScatterEnd()` call except NO parallel communication
70 : is done. Any variables that have be moved between processes are ignored
71 :
72 : Level: developer
73 :
74 : .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_REVERSE`, `SCATTER_FORWARD`,
75 : `SCATTER_REVERSE_LOCAL`
76 : M*/
77 :
78 : /*MC
79 : SCATTER_REVERSE_LOCAL - Moves the values in the opposite direction then the directions indicated
80 : in the `VecScatterCreate()` during `VecScatterBegin()` and `VecScatterEnd()` except NO parallel communication
81 : is done. Any variables that have be moved between processes are ignored
82 :
83 : Level: developer
84 :
85 : .seealso: [](ch_vectors), `VecScatter`, `ScatterMode`, `VecScatterCreate()`, `VecScatterBegin()`, `VecScatterEnd()`, `SCATTER_FORWARD`, `SCATTER_FORWARD_LOCAL`,
86 : `SCATTER_REVERSE`
87 : M*/
88 :
89 : /*J
90 : VecType - String with the name of a PETSc vector, `Vec`, type
91 :
92 : Level: beginner
93 :
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 */
114 :
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 *);
121 :
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;
127 :
128 : PETSC_EXTERN PetscErrorCode VecInitializePackage(void);
129 : PETSC_EXTERN PetscErrorCode VecFinalizePackage(void);
130 :
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 *);
138 :
139 : PETSC_EXTERN PetscErrorCode VecSetFromOptions(Vec);
140 : PETSC_EXTERN PetscErrorCode VecViewFromOptions(Vec, PetscObject, const char[]);
141 :
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 *[]);
148 : PETSC_EXTERN PetscErrorCode VecGetState(Vec, PetscObjectState *);
149 :
150 : PETSC_EXTERN PetscErrorCode VecSetSizes(Vec, PetscInt, PetscInt);
151 :
152 : PETSC_EXTERN PetscErrorCode VecDotNorm2(Vec, Vec, PetscScalar *, PetscReal *);
153 : PETSC_EXTERN PetscErrorCode VecDot(Vec, Vec, PetscScalar *);
154 : PETSC_EXTERN PetscErrorCode VecDotRealPart(Vec, Vec, PetscReal *);
155 : PETSC_EXTERN PetscErrorCode VecTDot(Vec, Vec, PetscScalar *);
156 : PETSC_EXTERN PetscErrorCode VecMDot(Vec, PetscInt, const Vec[], PetscScalar[]);
157 : PETSC_EXTERN PetscErrorCode VecMTDot(Vec, PetscInt, const Vec[], PetscScalar[]);
158 : PETSC_EXTERN PetscErrorCode VecGetSubVector(Vec, IS, Vec *);
159 : PETSC_EXTERN PetscErrorCode VecRestoreSubVector(Vec, IS, Vec *);
160 : PETSC_EXTERN PetscErrorCode VecConcatenate(PetscInt, const Vec[], Vec *, IS *[]);
161 :
162 : /*E
163 : NormType - determines what type of norm to compute
164 :
165 : Values:
166 : + `NORM_1` - the one norm, $||v|| = \sum_i | v_i |$. $||A|| = \max_j || v_*j ||$, maximum column sum
167 : . `NORM_2` - the two norm, $||v|| = sqrt(\sum_i |v_i|^2)$ (vectors only)
168 : . `NORM_FROBENIUS` - $||A|| = sqrt(\sum_ij |A_ij|^2)$, same as `NORM_2` for vectors
169 : . `NORM_INFINITY` - $||v|| = \max_i |v_i|$. $||A|| = \max_i || v_i* ||$, maximum row sum
170 : - `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 :
172 : Level: beginner
173 :
174 : Note:
175 : The `v` above represents a `Vec` while the `A` represents a `Mat`
176 :
177 : .seealso: [](ch_vectors), `Vec`, `Mat`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `MatNorm()`, `NORM_1`,
178 : `NORM_2`, `NORM_FROBENIUS`, `NORM_INFINITY`, `NORM_1_AND_2`, `ReductionType`
179 : E*/
180 : typedef enum {
181 : NORM_1 = 0,
182 : NORM_2 = 1,
183 : NORM_FROBENIUS = 2,
184 : NORM_INFINITY = 3,
185 : NORM_1_AND_2 = 4
186 : } NormType;
187 : PETSC_EXTERN const char *const NormTypes[];
188 : #define NORM_MAX NORM_INFINITY
189 :
190 : /*MC
191 : NORM_1 - the one norm, $||v|| = \sum_i | v_i |$. $||A|| = \max_j || v_{*,j} ||$, maximum column sum
192 :
193 : Level: beginner
194 :
195 : .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_2`, `NORM_FROBENIUS`,
196 : `NORM_INFINITY`, `NORM_1_AND_2`
197 : M*/
198 :
199 : /*MC
200 : NORM_2 - the two norm, $||v|| = \sqrt{\sum_i |v_i|^2}$ (vectors only)
201 :
202 : Level: beginner
203 :
204 : .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_FROBENIUS`,
205 : `NORM_INFINITY`, `NORM_1_AND_2`
206 : M*/
207 :
208 : /*MC
209 : NORM_FROBENIUS - $||A|| = \sqrt{\sum_{i,j} |A_{i,j}|^2}$, same as `NORM_2` for vectors
210 :
211 : Level: beginner
212 :
213 : .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
214 : `NORM_INFINITY`, `NORM_1_AND_2`
215 : M*/
216 :
217 : /*MC
218 : NORM_INFINITY - $||v|| = \max_i |v_i|$. $||A|| = \max_i || v_{i,*} ||$, maximum row sum
219 :
220 : Level: beginner
221 :
222 : .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
223 : `NORM_FROBENIUS`, `NORM_1_AND_2`
224 : M*/
225 :
226 : /*MC
227 : 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 :
229 : Level: beginner
230 :
231 : .seealso: [](ch_vectors), `NormType`, `MatNorm()`, `VecNorm()`, `VecNormBegin()`, `VecNormEnd()`, `NORM_1`, `NORM_2`,
232 : `NORM_FROBENIUS`, `NORM_INFINITY`
233 : M*/
234 :
235 : /*MC
236 : NORM_MAX - see `NORM_INFINITY`
237 :
238 : Level: beginner
239 : M*/
240 :
241 : /*E
242 : ReductionType - determines what type of column reduction (one that is not a type of norm defined in `NormType`) to obtain with `MatGetColumnReductions()`
243 :
244 : Values:
245 : + `REDUCTION_SUM_REALPART` - sum of real part of each matrix column
246 : . `REDUCTION_SUM_IMAGINARYPART` - sum of imaginary part of each matrix column
247 : . `REDUCTION_MEAN_REALPART` - arithmetic mean of real part of each matrix column
248 : - `REDUCTION_MEAN_IMAGINARYPART` - arithmetic mean of imaginary part of each matrix column
249 :
250 : Level: beginner
251 :
252 : Developer Note:
253 : The constants defined in `ReductionType` MUST BE DISTINCT from those defined in `NormType`.
254 : This is because `MatGetColumnReductions()` is used to compute both norms and other types of reductions,
255 : and the constants defined in both `NormType` and `ReductionType` are used to designate the desired operation.
256 :
257 : .seealso: [](ch_vectors), `MatGetColumnReductions()`, `MatGetColumnNorms()`, `NormType`, `REDUCTION_SUM_REALPART`,
258 : `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
259 : E*/
260 : typedef enum {
261 : REDUCTION_SUM_REALPART = 10,
262 : REDUCTION_MEAN_REALPART = 11,
263 : REDUCTION_SUM_IMAGINARYPART = 12,
264 : REDUCTION_MEAN_IMAGINARYPART = 13
265 : } ReductionType;
266 :
267 : /*MC
268 : REDUCTION_SUM_REALPART - sum of real part of a matrix column to obtain with `MatGetColumnReductions()`
269 :
270 : Level: beginner
271 :
272 : .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_SUM_IMAGINARYPART`, `REDUCTION_MEAN_REALPART`
273 : M*/
274 :
275 : /*MC
276 : REDUCTION_SUM_IMAGINARYPART - sum of imaginary part of matrix column to obtain with `MatGetColumnReductions()`
277 :
278 : Level: beginner
279 :
280 : .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_SUM_REALPART`, `REDUCTION_MEAN_IMAGINARYPART`
281 : M*/
282 :
283 : /*MC
284 : REDUCTION_MEAN_REALPART - arithmetic mean of real part of matrix column to obtain with `MatGetColumnReductions()`
285 :
286 : Level: beginner
287 :
288 : .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_IMAGINARYPART`, `REDUCTION_SUM_REALPART`
289 : M*/
290 :
291 : /*MC
292 : REDUCTION_MEAN_IMAGINARYPART - arithmetic mean of imaginary part of matrix column to obtain with `MatGetColumnReductions()`
293 :
294 : Level: beginner
295 :
296 : .seealso: [](ch_vectors), `ReductionType`, `MatGetColumnReductions()`, `REDUCTION_MEAN_REALPART`, `REDUCTION_SUM_IMAGINARYPART`
297 : M*/
298 :
299 : PETSC_EXTERN PetscErrorCode VecNorm(Vec, NormType, PetscReal *);
300 : PETSC_EXTERN PetscErrorCode VecNormAvailable(Vec, NormType, PetscBool *, PetscReal *);
301 : PETSC_EXTERN PetscErrorCode VecFlag(Vec, PetscInt);
302 : PETSC_EXTERN PetscErrorCode VecNormalize(Vec, PetscReal *);
303 : PETSC_EXTERN PetscErrorCode VecSum(Vec, PetscScalar *);
304 : PETSC_EXTERN PetscErrorCode VecMean(Vec, PetscScalar *);
305 : PETSC_EXTERN PetscErrorCode VecMax(Vec, PetscInt *, PetscReal *);
306 : PETSC_EXTERN PetscErrorCode VecMin(Vec, PetscInt *, PetscReal *);
307 : PETSC_EXTERN PetscErrorCode VecScale(Vec, PetscScalar);
308 : PETSC_EXTERN PetscErrorCode VecCopy(Vec, Vec);
309 : PETSC_EXTERN PetscErrorCode VecSetRandom(Vec, PetscRandom);
310 : PETSC_EXTERN PetscErrorCode VecSet(Vec, PetscScalar);
311 : PETSC_DEPRECATED_FUNCTION(3, 22, 0, "VecFlag()", ) PetscErrorCode VecSetInf(Vec);
312 : PETSC_EXTERN PetscErrorCode VecSwap(Vec, Vec);
313 : PETSC_EXTERN PetscErrorCode VecAXPY(Vec, PetscScalar, Vec);
314 : PETSC_EXTERN PetscErrorCode VecAXPBY(Vec, PetscScalar, PetscScalar, Vec);
315 : PETSC_EXTERN PetscErrorCode VecMAXPY(Vec, PetscInt, const PetscScalar[], Vec[]);
316 : PETSC_EXTERN PetscErrorCode VecMAXPBY(Vec, PetscInt, const PetscScalar[], PetscScalar, Vec[]);
317 : PETSC_EXTERN PetscErrorCode VecAYPX(Vec, PetscScalar, Vec);
318 : PETSC_EXTERN PetscErrorCode VecWAXPY(Vec, PetscScalar, Vec, Vec);
319 : PETSC_EXTERN PetscErrorCode VecAXPBYPCZ(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec);
320 : PETSC_EXTERN PetscErrorCode VecPointwiseMax(Vec, Vec, Vec);
321 : PETSC_EXTERN PetscErrorCode VecPointwiseMaxAbs(Vec, Vec, Vec);
322 : PETSC_EXTERN PetscErrorCode VecPointwiseMin(Vec, Vec, Vec);
323 : PETSC_EXTERN PetscErrorCode VecPointwiseMult(Vec, Vec, Vec);
324 : PETSC_EXTERN PetscErrorCode VecPointwiseDivide(Vec, Vec, Vec);
325 : PETSC_EXTERN PetscErrorCode VecMaxPointwiseDivide(Vec, Vec, PetscReal *);
326 : PETSC_EXTERN PetscErrorCode VecShift(Vec, PetscScalar);
327 : PETSC_EXTERN PetscErrorCode VecReciprocal(Vec);
328 : PETSC_EXTERN PetscErrorCode VecPermute(Vec, IS, PetscBool);
329 : PETSC_EXTERN PetscErrorCode VecSqrtAbs(Vec);
330 : PETSC_EXTERN PetscErrorCode VecLog(Vec);
331 : PETSC_EXTERN PetscErrorCode VecExp(Vec);
332 : PETSC_EXTERN PetscErrorCode VecAbs(Vec);
333 : PETSC_EXTERN PetscErrorCode VecDuplicate(Vec, Vec *);
334 : PETSC_EXTERN PetscErrorCode VecDuplicateVecs(Vec, PetscInt, Vec *[]);
335 : PETSC_EXTERN PetscErrorCode VecDestroyVecs(PetscInt, Vec *[]);
336 : PETSC_EXTERN PetscErrorCode VecStrideNormAll(Vec, NormType, PetscReal[]);
337 : PETSC_EXTERN PetscErrorCode VecStrideMaxAll(Vec, PetscInt[], PetscReal[]);
338 : PETSC_EXTERN PetscErrorCode VecStrideMinAll(Vec, PetscInt[], PetscReal[]);
339 : PETSC_EXTERN PetscErrorCode VecStrideScaleAll(Vec, const PetscScalar[]);
340 : PETSC_EXTERN PetscErrorCode VecStrideSumAll(Vec, PetscScalar *);
341 : PETSC_EXTERN PetscErrorCode VecUniqueEntries(Vec, PetscInt *, PetscScalar **);
342 :
343 : PETSC_EXTERN PetscErrorCode VecStrideNorm(Vec, PetscInt, NormType, PetscReal *);
344 : PETSC_EXTERN PetscErrorCode VecStrideMax(Vec, PetscInt, PetscInt *, PetscReal *);
345 : PETSC_EXTERN PetscErrorCode VecStrideMin(Vec, PetscInt, PetscInt *, PetscReal *);
346 : PETSC_EXTERN PetscErrorCode VecStrideScale(Vec, PetscInt, PetscScalar);
347 : PETSC_EXTERN PetscErrorCode VecStrideSum(Vec, PetscInt, PetscScalar *);
348 : PETSC_EXTERN PetscErrorCode VecStrideSet(Vec, PetscInt, PetscScalar);
349 :
350 : PETSC_EXTERN PetscErrorCode VecStrideGather(Vec, PetscInt, Vec, InsertMode);
351 : PETSC_EXTERN PetscErrorCode VecStrideScatter(Vec, PetscInt, Vec, InsertMode);
352 : PETSC_EXTERN PetscErrorCode VecStrideGatherAll(Vec, Vec[], InsertMode);
353 : PETSC_EXTERN PetscErrorCode VecStrideScatterAll(Vec[], Vec, InsertMode);
354 :
355 : PETSC_EXTERN PetscErrorCode VecStrideSubSetScatter(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
356 : PETSC_EXTERN PetscErrorCode VecStrideSubSetGather(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
357 :
358 : PETSC_EXTERN PetscErrorCode VecSetValues(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
359 : PETSC_EXTERN PetscErrorCode VecGetValues(Vec, PetscInt, const PetscInt[], PetscScalar[]);
360 : PETSC_EXTERN PetscErrorCode VecAssemblyBegin(Vec);
361 : PETSC_EXTERN PetscErrorCode VecAssemblyEnd(Vec);
362 : PETSC_EXTERN PetscErrorCode VecStashSetInitialSize(Vec, PetscInt, PetscInt);
363 : PETSC_EXTERN PetscErrorCode VecStashView(Vec, PetscViewer);
364 : PETSC_EXTERN PetscErrorCode VecStashViewFromOptions(Vec, PetscObject, const char[]);
365 : PETSC_EXTERN PetscErrorCode VecStashGetInfo(Vec, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
366 :
367 : PETSC_EXTERN PetscErrorCode VecSetPreallocationCOO(Vec, PetscCount, const PetscInt[]);
368 : PETSC_EXTERN PetscErrorCode VecSetPreallocationCOOLocal(Vec, PetscCount, PetscInt[]);
369 : PETSC_EXTERN PetscErrorCode VecSetValuesCOO(Vec, const PetscScalar[], InsertMode);
370 :
371 : /*@C
372 : VecSetValue - Set a single entry into a PETSc vector, `Vec`.
373 :
374 : Not Collective
375 :
376 : Input Parameters:
377 : + v - the vector
378 : . row - the row location of the entry
379 : . value - the value to insert
380 : - mode - either `INSERT_VALUES` or `ADD_VALUES`
381 :
382 : Level: beginner
383 :
384 : Notes:
385 : For efficiency one should use `VecSetValues()` and set several or
386 : many values simultaneously if possible.
387 :
388 : These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
389 : MUST be called after all calls to `VecSetValue()` have been completed.
390 :
391 : `VecSetValue()` uses 0-based indices in Python, C, and Fortran
392 :
393 : .seealso: [](ch_vectors), `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValueLocal()`
394 : @*/
395 261678 : static inline PetscErrorCode VecSetValue(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
396 : {
397 261678 : return VecSetValues(v, 1, &i, &va, mode);
398 : }
399 :
400 : PETSC_EXTERN PetscErrorCode VecSetBlockSize(Vec, PetscInt);
401 : PETSC_EXTERN PetscErrorCode VecGetBlockSize(Vec, PetscInt *);
402 : PETSC_EXTERN PetscErrorCode VecSetValuesBlocked(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
403 :
404 : /* Dynamic creation and loading functions */
405 : PETSC_EXTERN PetscFunctionList VecList;
406 : PETSC_EXTERN PetscErrorCode VecSetType(Vec, VecType);
407 : PETSC_EXTERN PetscErrorCode VecGetType(Vec, VecType *);
408 : PETSC_EXTERN PetscErrorCode VecRegister(const char[], PetscErrorCode (*)(Vec));
409 : PETSC_EXTERN PetscErrorCode VecRegisterAll(void);
410 :
411 : PETSC_EXTERN PetscErrorCode VecScatterBegin(VecScatter, Vec, Vec, InsertMode, ScatterMode);
412 : PETSC_EXTERN PetscErrorCode VecScatterEnd(VecScatter, Vec, Vec, InsertMode, ScatterMode);
413 : PETSC_EXTERN PetscErrorCode VecScatterDestroy(VecScatter *);
414 : PETSC_EXTERN PetscErrorCode VecScatterSetUp(VecScatter);
415 : PETSC_EXTERN PetscErrorCode VecScatterCopy(VecScatter, VecScatter *);
416 : PETSC_EXTERN PetscErrorCode VecScatterView(VecScatter, PetscViewer);
417 : PETSC_EXTERN PetscErrorCode VecScatterViewFromOptions(VecScatter, PetscObject, const char[]);
418 : PETSC_EXTERN PetscErrorCode VecScatterRemap(VecScatter, PetscInt[], PetscInt[]);
419 : PETSC_EXTERN PetscErrorCode VecScatterGetMerged(VecScatter, PetscBool *);
420 :
421 : PETSC_EXTERN PetscErrorCode VecGetArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
422 : PETSC_EXTERN PetscErrorCode VecRestoreArray4d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
423 : PETSC_EXTERN PetscErrorCode VecGetArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
424 : PETSC_EXTERN PetscErrorCode VecRestoreArray3d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
425 : PETSC_EXTERN PetscErrorCode VecGetArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
426 : PETSC_EXTERN PetscErrorCode VecRestoreArray2d(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
427 : PETSC_EXTERN PetscErrorCode VecGetArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);
428 : PETSC_EXTERN PetscErrorCode VecRestoreArray1d(Vec, PetscInt, PetscInt, PetscScalar *[]);
429 :
430 : PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
431 : PETSC_EXTERN PetscErrorCode VecGetArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
432 : PETSC_EXTERN PetscErrorCode VecRestoreArray4dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
433 : PETSC_EXTERN PetscErrorCode VecGetArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
434 : PETSC_EXTERN PetscErrorCode VecRestoreArray3dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
435 : PETSC_EXTERN PetscErrorCode VecGetArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
436 : PETSC_EXTERN PetscErrorCode VecRestoreArray2dWrite(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
437 : PETSC_EXTERN PetscErrorCode VecGetArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);
438 : PETSC_EXTERN PetscErrorCode VecRestoreArray1dWrite(Vec, PetscInt, PetscInt, PetscScalar *[]);
439 :
440 : PETSC_EXTERN PetscErrorCode VecGetArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
441 : PETSC_EXTERN PetscErrorCode VecRestoreArray4dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ****[]);
442 : PETSC_EXTERN PetscErrorCode VecGetArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
443 : PETSC_EXTERN PetscErrorCode VecRestoreArray3dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar ***[]);
444 : PETSC_EXTERN PetscErrorCode VecGetArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
445 : PETSC_EXTERN PetscErrorCode VecRestoreArray2dRead(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar **[]);
446 : PETSC_EXTERN PetscErrorCode VecGetArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);
447 : PETSC_EXTERN PetscErrorCode VecRestoreArray1dRead(Vec, PetscInt, PetscInt, PetscScalar *[]);
448 :
449 : PETSC_EXTERN PetscErrorCode VecPlaceArray(Vec, const PetscScalar[]);
450 : PETSC_EXTERN PetscErrorCode VecResetArray(Vec);
451 : PETSC_EXTERN PetscErrorCode VecReplaceArray(Vec, const PetscScalar[]);
452 :
453 : PETSC_EXTERN PetscErrorCode VecGetArrays(const Vec[], PetscInt, PetscScalar **[]);
454 : PETSC_EXTERN PetscErrorCode VecRestoreArrays(const Vec[], PetscInt, PetscScalar **[]);
455 :
456 : PETSC_EXTERN PetscErrorCode VecView(Vec, PetscViewer);
457 : PETSC_EXTERN PetscErrorCode VecViewNative(Vec, PetscViewer);
458 : PETSC_EXTERN PetscErrorCode VecEqual(Vec, Vec, PetscBool *);
459 : PETSC_EXTERN PetscErrorCode VecLoad(Vec, PetscViewer);
460 :
461 : PETSC_EXTERN PetscErrorCode VecGetSize(Vec, PetscInt *);
462 : PETSC_EXTERN PetscErrorCode VecGetLocalSize(Vec, PetscInt *);
463 : PETSC_EXTERN PetscErrorCode VecGetOwnershipRange(Vec, PetscInt *, PetscInt *);
464 : PETSC_EXTERN PetscErrorCode VecGetOwnershipRanges(Vec, const PetscInt *[]);
465 :
466 : PETSC_EXTERN PetscErrorCode VecSetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping);
467 : PETSC_EXTERN PetscErrorCode VecSetValuesLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
468 :
469 : PETSC_EXTERN PetscErrorCode VecViennaCLGetCLContext(Vec, PETSC_UINTPTR_T *);
470 : PETSC_EXTERN PetscErrorCode VecViennaCLGetCLQueue(Vec, PETSC_UINTPTR_T *);
471 : PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemRead(Vec, PETSC_UINTPTR_T *);
472 : PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMemWrite(Vec, PETSC_UINTPTR_T *);
473 : PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMemWrite(Vec);
474 : PETSC_EXTERN PetscErrorCode VecViennaCLGetCLMem(Vec, PETSC_UINTPTR_T *);
475 : PETSC_EXTERN PetscErrorCode VecViennaCLRestoreCLMem(Vec);
476 :
477 : /*@C
478 : VecSetValueLocal - Set a single entry into a vector using the local numbering of the vector, see `VecSetValuesLocal()`
479 :
480 : Not Collective
481 :
482 : Input Parameters:
483 : + v - the vector
484 : . row - the local row location of the entry
485 : . value - the value to insert
486 : - mode - either `INSERT_VALUES` or `ADD_VALUES`
487 :
488 : Level: beginner
489 :
490 : Notes:
491 : For efficiency one should use `VecSetValuesLocal()` and set several or
492 : many values simultaneously if possible.
493 :
494 : These values may be cached, so `VecAssemblyBegin()` and `VecAssemblyEnd()`
495 : MUST be called after all calls to `VecSetValueLocal()` have been completed.
496 :
497 : See `VecSetLocalToGlobalMapping()` for how the local numbering is defined
498 :
499 : `VecSetValueLocal()` uses 0-based indices in Fortran as well as in C.
500 :
501 : .seealso: [](ch_vectors), `VecSetValuesLocal()`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`, `VecSetValue()`,
502 : `VecSetLocalToGlobalMapping()`
503 : @*/
504 : static inline PetscErrorCode VecSetValueLocal(Vec v, PetscInt i, PetscScalar va, InsertMode mode)
505 : {
506 : return VecSetValuesLocal(v, 1, &i, &va, mode);
507 : }
508 :
509 : /*MC
510 : VecCheckAssembled - checks if values have been changed in the vector, by `VecSetValues()` or related routines, but it has not been assembled
511 :
512 : Synopsis:
513 : #include <petscvec.h>
514 : VecCheckAssembled(Vec v);
515 :
516 : Not Collective
517 :
518 : Input Parameter:
519 : . v - the vector to check
520 :
521 : Level: developer
522 :
523 : Note:
524 : After calls to `VecSetValues()` and related routines one must call `VecAssemblyBegin()` and `VecAssemblyEnd()` before using the vector
525 :
526 : .seealso: [](ch_vectors), `Vec`, `VecSetValues()`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`
527 : M*/
528 : #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()?")
529 :
530 : PETSC_EXTERN PetscErrorCode VecSetValuesBlockedLocal(Vec, PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
531 : PETSC_EXTERN PetscErrorCode VecGetLocalToGlobalMapping(Vec, ISLocalToGlobalMapping *);
532 :
533 : PETSC_EXTERN PetscErrorCode VecDotBegin(Vec, Vec, PetscScalar *);
534 : PETSC_EXTERN PetscErrorCode VecDotEnd(Vec, Vec, PetscScalar *);
535 : PETSC_EXTERN PetscErrorCode VecTDotBegin(Vec, Vec, PetscScalar *);
536 : PETSC_EXTERN PetscErrorCode VecTDotEnd(Vec, Vec, PetscScalar *);
537 : PETSC_EXTERN PetscErrorCode VecNormBegin(Vec, NormType, PetscReal *);
538 : PETSC_EXTERN PetscErrorCode VecNormEnd(Vec, NormType, PetscReal *);
539 : PETSC_EXTERN PetscErrorCode VecErrorWeightedNorms(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *);
540 :
541 : PETSC_EXTERN PetscErrorCode VecMDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
542 : PETSC_EXTERN PetscErrorCode VecMDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
543 : PETSC_EXTERN PetscErrorCode VecMTDotBegin(Vec, PetscInt, const Vec[], PetscScalar[]);
544 : PETSC_EXTERN PetscErrorCode VecMTDotEnd(Vec, PetscInt, const Vec[], PetscScalar[]);
545 : PETSC_EXTERN PetscErrorCode PetscCommSplitReductionBegin(MPI_Comm);
546 :
547 : PETSC_EXTERN PetscErrorCode VecBindToCPU(Vec, PetscBool);
548 : PETSC_DEPRECATED_FUNCTION(3, 13, 0, "VecBindToCPU()", ) static inline PetscErrorCode VecPinToCPU(Vec v, PetscBool flg)
549 : {
550 : return VecBindToCPU(v, flg);
551 : }
552 : PETSC_EXTERN PetscErrorCode VecBoundToCPU(Vec, PetscBool *);
553 : PETSC_EXTERN PetscErrorCode VecSetBindingPropagates(Vec, PetscBool);
554 : PETSC_EXTERN PetscErrorCode VecGetBindingPropagates(Vec, PetscBool *);
555 : PETSC_EXTERN PetscErrorCode VecSetPinnedMemoryMin(Vec, size_t);
556 : PETSC_EXTERN PetscErrorCode VecGetPinnedMemoryMin(Vec, size_t *);
557 :
558 : PETSC_EXTERN PetscErrorCode VecGetOffloadMask(Vec, PetscOffloadMask *);
559 :
560 : typedef enum {
561 : VEC_IGNORE_OFF_PROC_ENTRIES,
562 : VEC_IGNORE_NEGATIVE_INDICES,
563 : VEC_SUBSET_OFF_PROC_ENTRIES
564 : } VecOption;
565 : PETSC_EXTERN PetscErrorCode VecSetOption(Vec, VecOption, PetscBool);
566 :
567 : PETSC_EXTERN PetscErrorCode VecGetArray(Vec, PetscScalar **);
568 : PETSC_EXTERN PetscErrorCode VecGetArrayWrite(Vec, PetscScalar **);
569 : PETSC_EXTERN PetscErrorCode VecGetArrayRead(Vec, const PetscScalar **);
570 : PETSC_EXTERN PetscErrorCode VecRestoreArray(Vec, PetscScalar **);
571 : PETSC_EXTERN PetscErrorCode VecRestoreArrayWrite(Vec, PetscScalar **);
572 : PETSC_EXTERN PetscErrorCode VecRestoreArrayRead(Vec, const PetscScalar **);
573 : PETSC_EXTERN PetscErrorCode VecCreateLocalVector(Vec, Vec *);
574 : PETSC_EXTERN PetscErrorCode VecGetLocalVector(Vec, Vec);
575 : PETSC_EXTERN PetscErrorCode VecRestoreLocalVector(Vec, Vec);
576 : PETSC_EXTERN PetscErrorCode VecGetLocalVectorRead(Vec, Vec);
577 : PETSC_EXTERN PetscErrorCode VecRestoreLocalVectorRead(Vec, Vec);
578 : PETSC_EXTERN PetscErrorCode VecGetArrayAndMemType(Vec, PetscScalar **, PetscMemType *);
579 : PETSC_EXTERN PetscErrorCode VecRestoreArrayAndMemType(Vec, PetscScalar **);
580 : PETSC_EXTERN PetscErrorCode VecGetArrayReadAndMemType(Vec, const PetscScalar **, PetscMemType *);
581 : PETSC_EXTERN PetscErrorCode VecRestoreArrayReadAndMemType(Vec, const PetscScalar **);
582 : PETSC_EXTERN PetscErrorCode VecGetArrayWriteAndMemType(Vec, PetscScalar **, PetscMemType *);
583 : PETSC_EXTERN PetscErrorCode VecRestoreArrayWriteAndMemType(Vec, PetscScalar **);
584 :
585 : /*@C
586 : 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
587 :
588 : Logically Collective; No Fortran Support
589 :
590 : Input Parameters:
591 : + x - the vector
592 : - y - the second vector
593 :
594 : Output Parameters:
595 : + xv - location to put pointer to the first array
596 : - yv - location to put pointer to the second array
597 :
598 : Level: developer
599 :
600 : .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecRestoreArrayPair()`
601 : @*/
602 : static inline PetscErrorCode VecGetArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
603 : {
604 : PetscFunctionBegin;
605 : PetscCall(VecGetArray(y, yv));
606 : if (x == y) *xv = *yv;
607 : else PetscCall(VecGetArrayRead(x, (const PetscScalar **)xv));
608 : PetscFunctionReturn(PETSC_SUCCESS);
609 : }
610 :
611 : /*@C
612 : VecRestoreArrayPair - Returns a pair of pointers for two vectors that may be common obtained with `VecGetArrayPair()`
613 :
614 : Logically Collective; No Fortran Support
615 :
616 : Input Parameters:
617 : + x - the vector
618 : . y - the second vector
619 : . xv - location to put pointer to the first array
620 : - yv - location to put pointer to the second array
621 :
622 : Level: developer
623 :
624 : .seealso: [](ch_vectors), `VecGetArray()`, `VecGetArrayRead()`, `VecGetArrayPair()`
625 : @*/
626 : static inline PetscErrorCode VecRestoreArrayPair(Vec x, Vec y, PetscScalar **xv, PetscScalar **yv)
627 : {
628 : PetscFunctionBegin;
629 : PetscCall(VecRestoreArray(y, yv));
630 : if (x != y) PetscCall(VecRestoreArrayRead(x, (const PetscScalar **)xv));
631 : PetscFunctionReturn(PETSC_SUCCESS);
632 : }
633 :
634 : PETSC_EXTERN PetscErrorCode VecLockReadPush(Vec);
635 : PETSC_EXTERN PetscErrorCode VecLockReadPop(Vec);
636 : PETSC_EXTERN PetscErrorCode VecLockWriteSet(Vec, PetscBool);
637 : PETSC_EXTERN PetscErrorCode VecLockGet(Vec, PetscInt *);
638 : PETSC_EXTERN PetscErrorCode VecLockGetLocation(Vec, const char *[], const char *[], int *);
639 :
640 230423 : static inline PetscErrorCode VecSetErrorIfLocked(Vec x, PetscInt arg)
641 : {
642 230423 : PetscInt state;
643 :
644 230423 : PetscFunctionBegin;
645 230423 : PetscCall(VecLockGet(x, &state));
646 230423 : if (PetscUnlikely(state != 0)) {
647 0 : const char *file, *func, *name;
648 0 : int line;
649 :
650 0 : PetscCall(VecLockGetLocation(x, &file, &func, &line));
651 0 : PetscCall(PetscObjectGetName((PetscObject)x, &name));
652 0 : 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);
653 : }
654 230423 : PetscFunctionReturn(PETSC_SUCCESS);
655 : }
656 :
657 : /* The three are deprecated */
658 : PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPush()", ) static inline PetscErrorCode VecLockPush(Vec v)
659 : {
660 : return VecLockReadPush(v);
661 : }
662 :
663 : PETSC_DEPRECATED_FUNCTION(3, 11, 0, "VecLockReadPop()", ) static inline PetscErrorCode VecLockPop(Vec v)
664 : {
665 : return VecLockReadPop(v);
666 : }
667 :
668 : #define VecLocked(x, arg) VecSetErrorIfLocked(x, arg) PETSC_DEPRECATED_MACRO(3, 11, 0, "VecSetErrorIfLocked()", )
669 :
670 : /*E
671 : VecOperation - Enumeration of overide-able methods in the `Vec` implementation function-table.
672 :
673 : Values:
674 : + `VECOP_DUPLICATE` - `VecDuplicate()`
675 : . `VECOP_SET` - `VecSet()`
676 : . `VECOP_VIEW` - `VecView()`
677 : . `VECOP_LOAD` - `VecLoad()`
678 : . `VECOP_VIEWNATIVE` - `VecViewNative()`
679 : - `VECOP_LOADNATIVE` - `VecLoadNative()`
680 :
681 : Level: advanced
682 :
683 : Notes:
684 : Some operations may serve as the implementation for other routines not listed above. For
685 : example `VECOP_SET` can be used to simultaneously overriding the implementation used in
686 : `VecSet()`, `VecSetInf()`, and `VecZeroEntries()`.
687 :
688 : Entries to `VecOperation` are added as needed so if you do not see the operation listed which
689 : you'd like to replace, please send mail to `petsc-maint@mcs.anl.gov`!
690 :
691 : .seealso: [](ch_vectors), `Vec`, `VecSetOperation()`
692 : E*/
693 : typedef enum {
694 : VECOP_DUPLICATE = 0,
695 : VECOP_SET = 10,
696 : VECOP_VIEW = 33,
697 : VECOP_LOAD = 41,
698 : VECOP_VIEWNATIVE = 69,
699 : VECOP_LOADNATIVE = 70
700 : } VecOperation;
701 : PETSC_EXTERN PetscErrorCode VecSetOperation(Vec, VecOperation, void (*)(void));
702 :
703 : /*
704 : Routines for dealing with ghosted vectors:
705 : vectors with ghost elements at the end of the array.
706 : */
707 : PETSC_EXTERN PetscErrorCode VecMPISetGhost(Vec, PetscInt, const PetscInt[]);
708 : PETSC_EXTERN PetscErrorCode VecCreateGhost(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
709 : PETSC_EXTERN PetscErrorCode VecCreateGhostWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
710 : PETSC_EXTERN PetscErrorCode VecCreateGhostBlock(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Vec *);
711 : PETSC_EXTERN PetscErrorCode VecCreateGhostBlockWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], Vec *);
712 : PETSC_EXTERN PetscErrorCode VecGhostGetGhostIS(Vec, IS *);
713 : PETSC_EXTERN PetscErrorCode VecGhostGetLocalForm(Vec, Vec *);
714 : PETSC_EXTERN PetscErrorCode VecGhostRestoreLocalForm(Vec, Vec *);
715 : PETSC_EXTERN PetscErrorCode VecGhostIsLocalForm(Vec, Vec, PetscBool *);
716 : PETSC_EXTERN PetscErrorCode VecGhostUpdateBegin(Vec, InsertMode, ScatterMode);
717 : PETSC_EXTERN PetscErrorCode VecGhostUpdateEnd(Vec, InsertMode, ScatterMode);
718 :
719 : PETSC_EXTERN PetscErrorCode VecConjugate(Vec);
720 : PETSC_EXTERN PetscErrorCode VecImaginaryPart(Vec);
721 : PETSC_EXTERN PetscErrorCode VecRealPart(Vec);
722 :
723 : PETSC_EXTERN PetscErrorCode VecScatterCreateToAll(Vec, VecScatter *, Vec *);
724 : PETSC_EXTERN PetscErrorCode VecScatterCreateToZero(Vec, VecScatter *, Vec *);
725 :
726 : /* These vector calls were included from TAO. They miss vtable entries and GPU implementation */
727 : PETSC_EXTERN PetscErrorCode ISComplementVec(IS, Vec, IS *);
728 : PETSC_EXTERN PetscErrorCode VecPow(Vec, PetscScalar);
729 : PETSC_EXTERN PetscErrorCode VecMedian(Vec, Vec, Vec, Vec);
730 : PETSC_EXTERN PetscErrorCode VecWhichInactive(Vec, Vec, Vec, Vec, PetscBool, IS *);
731 : PETSC_EXTERN PetscErrorCode VecWhichBetween(Vec, Vec, Vec, IS *);
732 : PETSC_EXTERN PetscErrorCode VecWhichBetweenOrEqual(Vec, Vec, Vec, IS *);
733 : PETSC_EXTERN PetscErrorCode VecWhichGreaterThan(Vec, Vec, IS *);
734 : PETSC_EXTERN PetscErrorCode VecWhichLessThan(Vec, Vec, IS *);
735 : PETSC_EXTERN PetscErrorCode VecWhichEqual(Vec, Vec, IS *);
736 : PETSC_EXTERN PetscErrorCode VecISAXPY(Vec, IS, PetscScalar, Vec);
737 : PETSC_EXTERN PetscErrorCode VecISCopy(Vec, IS, ScatterMode, Vec);
738 : PETSC_EXTERN PetscErrorCode VecISSet(Vec, IS, PetscScalar);
739 : PETSC_EXTERN PetscErrorCode VecISShift(Vec, IS, PetscScalar);
740 : PETSC_EXTERN PetscErrorCode VecBoundGradientProjection(Vec, Vec, Vec, Vec, Vec);
741 : PETSC_EXTERN PetscErrorCode VecStepBoundInfo(Vec, Vec, Vec, Vec, PetscReal *, PetscReal *, PetscReal *);
742 : PETSC_EXTERN PetscErrorCode VecStepMax(Vec, Vec, PetscReal *);
743 : PETSC_EXTERN PetscErrorCode VecStepMaxBounded(Vec, Vec, Vec, Vec, PetscReal *);
744 :
745 : PETSC_EXTERN PetscErrorCode PetscViewerMathematicaGetVector(PetscViewer, Vec);
746 : PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutVector(PetscViewer, Vec);
747 :
748 : /*S
749 : Vecs - Collection of vectors where the data for the vectors is stored in
750 : one contiguous memory
751 :
752 : Level: advanced
753 :
754 : Notes:
755 : Temporary construct for handling multiply right-hand side solves
756 :
757 : This is faked by storing a single vector that has enough array space for
758 : n vectors
759 :
760 : S*/
761 : struct _n_Vecs {
762 : PetscInt n;
763 : Vec v;
764 : };
765 : typedef struct _n_Vecs *Vecs;
766 : PETSC_EXTERN PetscErrorCode VecsDestroy(Vecs);
767 : PETSC_EXTERN PetscErrorCode VecsCreateSeq(MPI_Comm, PetscInt, PetscInt, Vecs *);
768 : PETSC_EXTERN PetscErrorCode VecsCreateSeqWithArray(MPI_Comm, PetscInt, PetscInt, PetscScalar *, Vecs *);
769 : PETSC_EXTERN PetscErrorCode VecsDuplicate(Vecs, Vecs *);
770 :
771 : #if PetscDefined(HAVE_VIENNACL)
772 : typedef struct _p_PetscViennaCLIndices *PetscViennaCLIndices;
773 : PETSC_EXTERN PetscErrorCode VecViennaCLCopyToGPUSome_Public(Vec, PetscViennaCLIndices);
774 : PETSC_EXTERN PetscErrorCode VecViennaCLCopyFromGPUSome_Public(Vec, PetscViennaCLIndices);
775 : PETSC_EXTERN PetscErrorCode VecCreateSeqViennaCL(MPI_Comm, PetscInt, Vec *);
776 : PETSC_EXTERN PetscErrorCode VecCreateMPIViennaCL(MPI_Comm, PetscInt, PetscInt, Vec *);
777 : #endif
778 : #if PetscDefined(HAVE_CUDA) || PetscDefined(HAVE_HIP)
779 : PETSC_EXTERN PetscErrorCode VecScatterInitializeForGPU(VecScatter, Vec);
780 : PETSC_EXTERN PetscErrorCode VecScatterFinalizeForGPU(VecScatter);
781 : #endif
782 : #if PetscDefined(HAVE_KOKKOS_KERNELS)
783 : PETSC_EXTERN PetscErrorCode VecCreateSeqKokkos(MPI_Comm, PetscInt, Vec *);
784 : PETSC_EXTERN PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm, PetscInt, PetscInt, const PetscScalar *, Vec *);
785 : PETSC_EXTERN PetscErrorCode VecCreateMPIKokkos(MPI_Comm, PetscInt, PetscInt, Vec *);
786 : PETSC_EXTERN PetscErrorCode VecCreateMPIKokkosWithArray(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar *, Vec *);
787 : #endif
788 :
789 : PETSC_EXTERN PetscErrorCode VecNestGetSubVecs(Vec, PetscInt *, Vec **);
790 : PETSC_EXTERN PetscErrorCode VecNestGetSubVec(Vec, PetscInt, Vec *);
791 : PETSC_EXTERN PetscErrorCode VecNestSetSubVecs(Vec, PetscInt, PetscInt *, Vec *);
792 : PETSC_EXTERN PetscErrorCode VecNestSetSubVec(Vec, PetscInt, Vec);
793 : PETSC_EXTERN PetscErrorCode VecCreateNest(MPI_Comm, PetscInt, IS *, Vec *, Vec *);
794 : PETSC_EXTERN PetscErrorCode VecNestGetSize(Vec, PetscInt *);
795 :
796 : PETSC_EXTERN PetscErrorCode PetscOptionsGetVec(PetscOptions, const char[], const char[], Vec, PetscBool *);
797 : PETSC_EXTERN PetscErrorCode VecFilter(Vec, PetscReal);
798 : PETSC_DEPRECATED_FUNCTION(3, 20, 0, "VecFilter()", ) static inline PetscErrorCode VecChop(Vec v, PetscReal tol)
799 : {
800 : return VecFilter(v, tol);
801 : }
802 :
803 : PETSC_EXTERN PetscErrorCode VecGetLayout(Vec, PetscLayout *);
804 : PETSC_EXTERN PetscErrorCode VecSetLayout(Vec, PetscLayout);
805 :
806 : PETSC_EXTERN PetscErrorCode PetscSectionVecView(PetscSection, Vec, PetscViewer);
807 : PETSC_EXTERN PetscErrorCode VecGetValuesSection(Vec, PetscSection, PetscInt, PetscScalar **);
808 : PETSC_EXTERN PetscErrorCode VecSetValuesSection(Vec, PetscSection, PetscInt, const PetscScalar[], InsertMode);
809 : PETSC_EXTERN PetscErrorCode PetscSectionVecNorm(PetscSection, PetscSection, Vec, NormType, PetscReal[]);
810 :
811 : /*S
812 : VecTagger - Object used to manage the tagging of a subset of indices based on the values of a vector. The
813 : motivating application is the selection of cells for refinement or coarsening based on vector containing
814 : the values in an error indicator metric.
815 :
816 : Values:
817 : + `VECTAGGERABSOLUTE` - "absolute" values are in a interval (box for complex values) of explicitly defined values
818 : . `VECTAGGERRELATIVE` - "relative" values are in a interval (box for complex values) of values relative to the set of all values in the vector
819 : . `VECTAGGERCDF` - "cdf" values are in a relative range of the *cumulative distribution* of values in the vector
820 : . `VECTAGGEROR` - "or" values are in the union of other tags
821 : - `VECTAGGERAND` - "and" values are in the intersection of other tags
822 :
823 : Level: advanced
824 :
825 : Developer Note:
826 : Why not use a `DMLabel` or similar object
827 :
828 : .seealso: [](ch_vectors), `Vec`, `VecTaggerType`, `VecTaggerCreate()`
829 : S*/
830 : typedef struct _p_VecTagger *VecTagger;
831 :
832 : /*J
833 : VecTaggerType - String with the name of a `VecTagger` type
834 :
835 : Level: advanced
836 :
837 : .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerCreate()`
838 : J*/
839 : typedef const char *VecTaggerType;
840 : #define VECTAGGERABSOLUTE "absolute"
841 : #define VECTAGGERRELATIVE "relative"
842 : #define VECTAGGERCDF "cdf"
843 : #define VECTAGGEROR "or"
844 : #define VECTAGGERAND "and"
845 :
846 : PETSC_EXTERN PetscClassId VEC_TAGGER_CLASSID;
847 : PETSC_EXTERN PetscFunctionList VecTaggerList;
848 : PETSC_EXTERN PetscErrorCode VecTaggerRegister(const char[], PetscErrorCode (*)(VecTagger));
849 : PETSC_EXTERN PetscErrorCode VecTaggerRegisterAll(void);
850 :
851 : PETSC_EXTERN PetscErrorCode VecTaggerCreate(MPI_Comm, VecTagger *);
852 : PETSC_EXTERN PetscErrorCode VecTaggerSetBlockSize(VecTagger, PetscInt);
853 : PETSC_EXTERN PetscErrorCode VecTaggerGetBlockSize(VecTagger, PetscInt *);
854 : PETSC_EXTERN PetscErrorCode VecTaggerSetType(VecTagger, VecTaggerType);
855 : PETSC_EXTERN PetscErrorCode VecTaggerGetType(VecTagger, VecTaggerType *);
856 : PETSC_EXTERN PetscErrorCode VecTaggerSetInvert(VecTagger, PetscBool);
857 : PETSC_EXTERN PetscErrorCode VecTaggerGetInvert(VecTagger, PetscBool *);
858 : PETSC_EXTERN PetscErrorCode VecTaggerSetFromOptions(VecTagger);
859 : PETSC_EXTERN PetscErrorCode VecTaggerSetUp(VecTagger);
860 : PETSC_EXTERN PetscErrorCode VecTaggerView(VecTagger, PetscViewer);
861 : PETSC_EXTERN PetscErrorCode VecTaggerComputeIS(VecTagger, Vec, IS *, PetscBool *);
862 : PETSC_EXTERN PetscErrorCode VecTaggerDestroy(VecTagger *);
863 :
864 : /*S
865 : VecTaggerBox - A interval (box for complex numbers) range used to tag values. For real scalars, this is just a closed interval; for complex scalars,
866 : 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.
867 :
868 : Level: beginner
869 :
870 : .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerComputeIntervals()`
871 : S*/
872 : typedef struct {
873 : PetscScalar min;
874 : PetscScalar max;
875 : } VecTaggerBox;
876 : PETSC_EXTERN PetscErrorCode VecTaggerComputeBoxes(VecTagger, Vec, PetscInt *, VecTaggerBox **, PetscBool *);
877 :
878 : PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteSetBox(VecTagger, VecTaggerBox *);
879 : PETSC_EXTERN PetscErrorCode VecTaggerAbsoluteGetBox(VecTagger, const VecTaggerBox **);
880 :
881 : PETSC_EXTERN PetscErrorCode VecTaggerRelativeSetBox(VecTagger, VecTaggerBox *);
882 : PETSC_EXTERN PetscErrorCode VecTaggerRelativeGetBox(VecTagger, const VecTaggerBox **);
883 :
884 : PETSC_EXTERN PetscErrorCode VecTaggerCDFSetBox(VecTagger, VecTaggerBox *);
885 : PETSC_EXTERN PetscErrorCode VecTaggerCDFGetBox(VecTagger, const VecTaggerBox **);
886 :
887 : /*E
888 : 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).
889 :
890 : Values:
891 : + `VECTAGGER_CDF_GATHER` - gather the data to MPI rank 0, perform the computation and broadcast the result
892 : - `VECTAGGER_CDF_ITERATIVE` - compute the results on all ranks iteratively using `MPI_Allreduce()`
893 :
894 : Level: advanced
895 :
896 : Note:
897 : Relevant only in parallel: in serial it is directly computed.
898 :
899 : .seealso: [](ch_vectors), `Vec`, `VecTagger`, `VecTaggerType`, `VecTaggerCreate()`, `VecTaggerCDFSetMethod()`, `VecTaggerCDFMethods`
900 : E*/
901 : typedef enum {
902 : VECTAGGER_CDF_GATHER,
903 : VECTAGGER_CDF_ITERATIVE,
904 : VECTAGGER_CDF_NUM_METHODS
905 : } VecTaggerCDFMethod;
906 : PETSC_EXTERN const char *const VecTaggerCDFMethods[];
907 :
908 : PETSC_EXTERN PetscErrorCode VecTaggerCDFSetMethod(VecTagger, VecTaggerCDFMethod);
909 : PETSC_EXTERN PetscErrorCode VecTaggerCDFGetMethod(VecTagger, VecTaggerCDFMethod *);
910 : PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeSetTolerances(VecTagger, PetscInt, PetscReal, PetscReal);
911 : PETSC_EXTERN PetscErrorCode VecTaggerCDFIterativeGetTolerances(VecTagger, PetscInt *, PetscReal *, PetscReal *);
912 :
913 : PETSC_EXTERN PetscErrorCode VecTaggerOrSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
914 : PETSC_EXTERN PetscErrorCode VecTaggerOrGetSubs(VecTagger, PetscInt *, VecTagger **);
915 :
916 : PETSC_EXTERN PetscErrorCode VecTaggerAndSetSubs(VecTagger, PetscInt, VecTagger *, PetscCopyMode);
917 : PETSC_EXTERN PetscErrorCode VecTaggerAndGetSubs(VecTagger, PetscInt *, VecTagger **);
918 :
919 : PETSC_EXTERN PetscErrorCode VecTaggerInitializePackage(void);
920 : PETSC_EXTERN PetscErrorCode VecTaggerFinalizePackage(void);
921 :
922 : #if PetscDefined(USE_DEBUG)
923 : /* This is an internal debug-only routine that should not be used by users */
924 : PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecValidValues_Internal(Vec, PetscInt, PetscBool);
925 : #else
926 : #define VecValidValues_Internal(...) PETSC_SUCCESS
927 : #endif /* PETSC_USE_DEBUG */
928 :
929 : #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)
930 :
931 : #if PetscDefined(HAVE_CUDA)
932 : #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
933 : #else
934 : #define VEC_CUDA_DECL_OR_STUB(__decl__, ...) \
935 : static inline __decl__ \
936 : { \
937 : __VA_ARGS__; \
938 : VEC_CUPM_NOT_CONFIGURED(cuda); \
939 : }
940 : #endif /* PETSC_HAVE_CUDA */
941 :
942 : /* extra underscore here to make it line up with the cuda versions */
943 : #if PetscDefined(HAVE_HIP)
944 : #define VEC_HIP__DECL_OR_STUB(__decl__, ...) PETSC_EXTERN __decl__;
945 : #else
946 : #define VEC_HIP__DECL_OR_STUB(__decl__, ...) \
947 : static inline __decl__ \
948 : { \
949 : __VA_ARGS__; \
950 : VEC_CUPM_NOT_CONFIGURED(hip); \
951 : }
952 : #endif /* PETSC_HAVE_HIP */
953 :
954 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateSeqCUDA(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)
955 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateSeqHIP(MPI_Comm a, PetscInt b, Vec *c), (void)a, (void)b, (void)c)
956 :
957 : 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)
958 : 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)
959 :
960 : 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)
961 : 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)
962 :
963 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCreateMPICUDA(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)
964 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecCreateMPIHIP(MPI_Comm a, PetscInt b, PetscInt c, Vec *d), (void)a, (void)b, (void)c, (void)d)
965 :
966 : 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)
967 : 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)
968 :
969 : 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)
970 : 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)
971 :
972 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArray(Vec a, PetscScalar **b), (void)a, (void)b)
973 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArray(Vec a, PetscScalar **b), (void)a, (void)b)
974 :
975 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)
976 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArray(Vec a, PetscScalar **b), (void)a, (void)b)
977 :
978 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
979 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
980 :
981 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
982 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayRead(Vec a, const PetscScalar **b), (void)a, (void)b)
983 :
984 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
985 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPGetArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
986 :
987 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDARestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
988 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPRestoreArrayWrite(Vec a, PetscScalar **b), (void)a, (void)b)
989 :
990 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
991 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPPlaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
992 :
993 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
994 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPReplaceArray(Vec a, const PetscScalar b[]), (void)a, (void)b)
995 :
996 : VEC_CUDA_DECL_OR_STUB(PetscErrorCode VecCUDAResetArray(Vec a), (void)a)
997 : VEC_HIP__DECL_OR_STUB(PetscErrorCode VecHIPResetArray(Vec a), (void)a)
998 :
999 : #undef VEC_CUPM_NOT_CONFIGURED
1000 : #undef VEC_CUDA_DECL_OR_STUB
1001 : #undef VEC_HIP__DECL_OR_STUB
|