Actual source code: rvector.c

  1: /*
  2:      Provides the interface functions for vector operations that have PetscScalar/PetscReal in the signature
  3:    These are the vector functions the user calls.
  4: */
  5: #include <petsc/private/vecimpl.h>

  7: PetscInt VecGetSubVectorSavedStateId = -1;

  9: #if PetscDefined(USE_DEBUG)
 10: // this is a no-op '0' macro in optimized builds
 11: PetscErrorCode VecValidValues_Internal(Vec vec, PetscInt argnum, PetscBool begin)
 12: {
 13:   PetscFunctionBegin;
 14:   if (vec->petscnative || vec->ops->getarray) {
 15:     PetscInt           n;
 16:     const PetscScalar *x;
 17:     PetscOffloadMask   mask;

 19:     PetscCall(VecGetOffloadMask(vec, &mask));
 20:     if (!PetscOffloadHost(mask)) PetscFunctionReturn(PETSC_SUCCESS);
 21:     PetscCall(VecGetLocalSize(vec, &n));
 22:     PetscCall(VecGetArrayRead(vec, &x));
 23:     for (PetscInt i = 0; i < n; i++) {
 24:       if (begin) {
 25:         PetscCheck(!PetscIsInfOrNanScalar(x[i]), PETSC_COMM_SELF, PETSC_ERR_FP, "Vec entry at local location %" PetscInt_FMT " is not-a-number or infinite at beginning of function: Parameter number %" PetscInt_FMT, i, argnum);
 26:       } else {
 27:         PetscCheck(!PetscIsInfOrNanScalar(x[i]), PETSC_COMM_SELF, PETSC_ERR_FP, "Vec entry at local location %" PetscInt_FMT " is not-a-number or infinite at end of function: Parameter number %" PetscInt_FMT, i, argnum);
 28:       }
 29:     }
 30:     PetscCall(VecRestoreArrayRead(vec, &x));
 31:   }
 32:   PetscFunctionReturn(PETSC_SUCCESS);
 33: }
 34: #endif

 36: /*@
 37:   VecMaxPointwiseDivide - Computes the maximum of the componentwise division `max = max_i abs(x[i]/y[i])`.

 39:   Logically Collective

 41:   Input Parameters:
 42: + x - the numerators
 43: - y - the denominators

 45:   Output Parameter:
 46: . max - the result

 48:   Level: advanced

 50:   Notes:
 51:   `x` and `y` may be the same vector

 53:   if a particular `y[i]` is zero, it is treated as 1 in the above formula

 55: .seealso: [](ch_vectors), `Vec`, `VecPointwiseDivide()`, `VecPointwiseMult()`, `VecPointwiseMax()`, `VecPointwiseMin()`, `VecPointwiseMaxAbs()`
 56: @*/
 57: PetscErrorCode VecMaxPointwiseDivide(Vec x, Vec y, PetscReal *max)
 58: {
 59:   PetscFunctionBegin;
 62:   PetscAssertPointer(max, 3);
 65:   PetscCheckSameTypeAndComm(x, 1, y, 2);
 66:   VecCheckSameSize(x, 1, y, 2);
 67:   VecCheckAssembled(x);
 68:   VecCheckAssembled(y);
 69:   PetscCall(VecLockReadPush(x));
 70:   PetscCall(VecLockReadPush(y));
 71:   PetscUseTypeMethod(x, maxpointwisedivide, y, max);
 72:   PetscCall(VecLockReadPop(x));
 73:   PetscCall(VecLockReadPop(y));
 74:   PetscFunctionReturn(PETSC_SUCCESS);
 75: }

 77: /*@
 78:   VecDot - Computes the vector dot product.

 80:   Collective

 82:   Input Parameters:
 83: + x - first vector
 84: - y - second vector

 86:   Output Parameter:
 87: . val - the dot product

 89:   Level: intermediate

 91:   Notes for Users of Complex Numbers:
 92:   For complex vectors, `VecDot()` computes
 93: $     val = (x,y) = y^H x,
 94:   where y^H denotes the conjugate transpose of y. Note that this corresponds to the usual "mathematicians" complex
 95:   inner product where the SECOND argument gets the complex conjugate. Since the `BLASdot()` complex conjugates the first
 96:   first argument we call the `BLASdot()` with the arguments reversed.

 98:   Use `VecTDot()` for the indefinite form
 99: $     val = (x,y) = y^T x,
100:   where y^T denotes the transpose of y.

102: .seealso: [](ch_vectors), `Vec`, `VecMDot()`, `VecTDot()`, `VecNorm()`, `VecDotBegin()`, `VecDotEnd()`, `VecDotRealPart()`
103: @*/
104: PetscErrorCode VecDot(Vec x, Vec y, PetscScalar *val)
105: {
106:   PetscFunctionBegin;
109:   PetscAssertPointer(val, 3);
112:   PetscCheckSameTypeAndComm(x, 1, y, 2);
113:   VecCheckSameSize(x, 1, y, 2);
114:   VecCheckAssembled(x);
115:   VecCheckAssembled(y);

117:   PetscCall(VecLockReadPush(x));
118:   PetscCall(VecLockReadPush(y));
119:   PetscCall(PetscLogEventBegin(VEC_Dot, x, y, 0, 0));
120:   PetscUseTypeMethod(x, dot, y, val);
121:   PetscCall(PetscLogEventEnd(VEC_Dot, x, y, 0, 0));
122:   PetscCall(VecLockReadPop(x));
123:   PetscCall(VecLockReadPop(y));
124:   PetscFunctionReturn(PETSC_SUCCESS);
125: }

127: /*@
128:   VecDotRealPart - Computes the real part of the vector dot product.

130:   Collective

132:   Input Parameters:
133: + x - first vector
134: - y - second vector

136:   Output Parameter:
137: . val - the real part of the dot product;

139:   Level: intermediate

141:   Notes for Users of Complex Numbers:
142:   See `VecDot()` for more details on the definition of the dot product for complex numbers

144:   For real numbers this returns the same value as `VecDot()`

146:   For complex numbers in C^n (that is a vector of n components with a complex number for each component) this is equal to the usual real dot product on the
147:   the space R^{2n} (that is a vector of 2n components with the real or imaginary part of the complex numbers for components)

149:   Developer Notes:
150:   This is not currently optimized to compute only the real part of the dot product.

152: .seealso: [](ch_vectors), `Vec`, `VecMDot()`, `VecTDot()`, `VecNorm()`, `VecDotBegin()`, `VecDotEnd()`, `VecDot()`, `VecDotNorm2()`
153: @*/
154: PetscErrorCode VecDotRealPart(Vec x, Vec y, PetscReal *val)
155: {
156:   PetscScalar fdot;

158:   PetscFunctionBegin;
159:   PetscCall(VecDot(x, y, &fdot));
160:   *val = PetscRealPart(fdot);
161:   PetscFunctionReturn(PETSC_SUCCESS);
162: }

164: /*@
165:   VecNorm  - Computes the vector norm.

167:   Collective

169:   Input Parameters:
170: + x    - the vector
171: - type - the type of the norm requested

173:   Output Parameter:
174: . val - the norm

176:   Level: intermediate

178:   Notes:
179:   See `NormType` for descriptions of each norm.

181:   For complex numbers `NORM_1` will return the traditional 1 norm of the 2 norm of the complex
182:   numbers; that is the 1 norm of the absolute values of the complex entries. In PETSc 3.6 and
183:   earlier releases it returned the 1 norm of the 1 norm of the complex entries (what is
184:   returned by the BLAS routine `asum()`). Both are valid norms but most people expect the former.

186:   This routine stashes the computed norm value, repeated calls before the vector entries are
187:   changed are then rapid since the precomputed value is immediately available. Certain vector
188:   operations such as `VecSet()` store the norms so the value is immediately available and does
189:   not need to be explicitly computed. `VecScale()` updates any stashed norm values, thus calls
190:   after `VecScale()` do not need to explicitly recompute the norm.

192: .seealso: [](ch_vectors), `Vec`, `NormType`, `VecDot()`, `VecTDot()`, `VecDotBegin()`, `VecDotEnd()`, `VecNormAvailable()`,
193:           `VecNormBegin()`, `VecNormEnd()`, `NormType()`
194: @*/
195: PetscErrorCode VecNorm(Vec x, NormType type, PetscReal *val)
196: {
197:   PetscBool flg = PETSC_TRUE;

199:   PetscFunctionBegin;
202:   VecCheckAssembled(x);
204:   PetscAssertPointer(val, 3);

206:   PetscCall(VecNormAvailable(x, type, &flg, val));
207:   // check that all MPI processes call this routine together and have same availability
208:   if (PetscDefined(USE_DEBUG)) {
209:     PetscMPIInt b0 = (PetscMPIInt)flg, b1[2], b2[2];
210:     b1[0]          = -b0;
211:     b1[1]          = b0;
212:     PetscCall(MPIU_Allreduce(b1, b2, 2, MPI_INT, MPI_MAX, PetscObjectComm((PetscObject)x)));
213:     PetscCheck(-b2[0] == b2[1], PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_WRONGSTATE, "Some MPI processes have cached %s norm, others do not. This may happen when some MPI processes call VecGetArray() and some others do not.", NormTypes[type]);
214:     if (flg) {
215:       PetscReal b1[2], b2[2];
216:       b1[0] = -(*val);
217:       b1[1] = *val;
218:       PetscCall(MPIU_Allreduce(b1, b2, 2, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)x)));
219:       PetscCheck(-b2[0] == b2[1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Difference in cached %s norms: local %g", NormTypes[type], (double)*val);
220:     }
221:   }
222:   if (flg) PetscFunctionReturn(PETSC_SUCCESS);

224:   PetscCall(VecLockReadPush(x));
225:   PetscCall(PetscLogEventBegin(VEC_Norm, x, 0, 0, 0));
226:   PetscUseTypeMethod(x, norm, type, val);
227:   PetscCall(PetscLogEventEnd(VEC_Norm, x, 0, 0, 0));
228:   PetscCall(VecLockReadPop(x));

230:   if (type != NORM_1_AND_2) PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[type], *val));
231:   PetscFunctionReturn(PETSC_SUCCESS);
232: }

234: /*@
235:   VecNormAvailable  - Returns the vector norm if it is already known. That is, it has been previously computed and cached in the vector

237:   Not Collective

239:   Input Parameters:
240: + x    - the vector
241: - type - one of `NORM_1` (sum_i |x[i]|), `NORM_2` sqrt(sum_i (x[i])^2), `NORM_INFINITY` max_i |x[i]|.  Also available
242:           `NORM_1_AND_2`, which computes both norms and stores them
243:           in a two element array.

245:   Output Parameters:
246: + available - `PETSC_TRUE` if the val returned is valid
247: - val       - the norm

249:   Level: intermediate

251:   Developer Notes:
252:   `PETSC_HAVE_SLOW_BLAS_NORM2` will cause a C (loop unrolled) version of the norm to be used, rather
253:   than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines
254:   (as opposed to vendor provided) because the FORTRAN BLAS `NRM2()` routine is very slow.

256: .seealso: [](ch_vectors), `Vec`, `VecDot()`, `VecTDot()`, `VecNorm()`, `VecDotBegin()`, `VecDotEnd()`,
257:           `VecNormBegin()`, `VecNormEnd()`
258: @*/
259: PetscErrorCode VecNormAvailable(Vec x, NormType type, PetscBool *available, PetscReal *val)
260: {
261:   PetscFunctionBegin;
264:   PetscAssertPointer(available, 3);
265:   PetscAssertPointer(val, 4);

267:   if (type == NORM_1_AND_2) {
268:     *available = PETSC_FALSE;
269:   } else {
270:     PetscCall(PetscObjectComposedDataGetReal((PetscObject)x, NormIds[type], *val, *available));
271:   }
272:   PetscFunctionReturn(PETSC_SUCCESS);
273: }

275: /*@
276:   VecNormalize - Normalizes a vector by its 2-norm.

278:   Collective

280:   Input Parameter:
281: . x - the vector

283:   Output Parameter:
284: . val - the vector norm before normalization. May be `NULL` if the value is not needed.

286:   Level: intermediate

288: .seealso: [](ch_vectors), `Vec`, `VecNorm()`, `NORM_2`, `NormType`
289: @*/
290: PetscErrorCode VecNormalize(Vec x, PetscReal *val)
291: {
292:   PetscReal norm;

294:   PetscFunctionBegin;
297:   PetscCall(VecSetErrorIfLocked(x, 1));
298:   if (val) PetscAssertPointer(val, 2);
299:   PetscCall(PetscLogEventBegin(VEC_Normalize, x, 0, 0, 0));
300:   PetscCall(VecNorm(x, NORM_2, &norm));
301:   if (norm == 0.0) PetscCall(PetscInfo(x, "Vector of zero norm can not be normalized; Returning only the zero norm\n"));
302:   else if (PetscIsInfOrNanReal(norm)) PetscCall(PetscInfo(x, "Vector with Inf or Nan norm can not be normalized; Returning only the norm\n"));
303:   else {
304:     PetscScalar s = 1.0 / norm;
305:     PetscCall(VecScale(x, s));
306:   }
307:   PetscCall(PetscLogEventEnd(VEC_Normalize, x, 0, 0, 0));
308:   if (val) *val = norm;
309:   PetscFunctionReturn(PETSC_SUCCESS);
310: }

312: /*@C
313:   VecMax - Determines the vector component with maximum real part and its location.

315:   Collective

317:   Input Parameter:
318: . x - the vector

320:   Output Parameters:
321: + p   - the index of `val` (pass `NULL` if you don't want this) in the vector
322: - val - the maximum component

324:   Level: intermediate

326:   Notes:
327:   Returns the value `PETSC_MIN_REAL` and negative `p` if the vector is of length 0.

329:   Returns the smallest index with the maximum value

331: .seealso: [](ch_vectors), `Vec`, `VecNorm()`, `VecMin()`
332: @*/
333: PetscErrorCode VecMax(Vec x, PetscInt *p, PetscReal *val)
334: {
335:   PetscFunctionBegin;
338:   VecCheckAssembled(x);
339:   if (p) PetscAssertPointer(p, 2);
340:   PetscAssertPointer(val, 3);
341:   PetscCall(VecLockReadPush(x));
342:   PetscCall(PetscLogEventBegin(VEC_Max, x, 0, 0, 0));
343:   PetscUseTypeMethod(x, max, p, val);
344:   PetscCall(PetscLogEventEnd(VEC_Max, x, 0, 0, 0));
345:   PetscCall(VecLockReadPop(x));
346:   PetscFunctionReturn(PETSC_SUCCESS);
347: }

349: /*@C
350:   VecMin - Determines the vector component with minimum real part and its location.

352:   Collective

354:   Input Parameter:
355: . x - the vector

357:   Output Parameters:
358: + p   - the index of `val` (pass `NULL` if you don't want this location) in the vector
359: - val - the minimum component

361:   Level: intermediate

363:   Notes:
364:   Returns the value `PETSC_MAX_REAL` and negative `p` if the vector is of length 0.

366:   This returns the smallest index with the minimum value

368: .seealso: [](ch_vectors), `Vec`, `VecMax()`
369: @*/
370: PetscErrorCode VecMin(Vec x, PetscInt *p, PetscReal *val)
371: {
372:   PetscFunctionBegin;
375:   VecCheckAssembled(x);
376:   if (p) PetscAssertPointer(p, 2);
377:   PetscAssertPointer(val, 3);
378:   PetscCall(VecLockReadPush(x));
379:   PetscCall(PetscLogEventBegin(VEC_Min, x, 0, 0, 0));
380:   PetscUseTypeMethod(x, min, p, val);
381:   PetscCall(PetscLogEventEnd(VEC_Min, x, 0, 0, 0));
382:   PetscCall(VecLockReadPop(x));
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

386: /*@
387:   VecTDot - Computes an indefinite vector dot product. That is, this
388:   routine does NOT use the complex conjugate.

390:   Collective

392:   Input Parameters:
393: + x - first vector
394: - y - second vector

396:   Output Parameter:
397: . val - the dot product

399:   Level: intermediate

401:   Notes for Users of Complex Numbers:
402:   For complex vectors, `VecTDot()` computes the indefinite form
403: $     val = (x,y) = y^T x,
404:   where y^T denotes the transpose of y.

406:   Use `VecDot()` for the inner product
407: $     val = (x,y) = y^H x,
408:   where y^H denotes the conjugate transpose of y.

410: .seealso: [](ch_vectors), `Vec`, `VecDot()`, `VecMTDot()`
411: @*/
412: PetscErrorCode VecTDot(Vec x, Vec y, PetscScalar *val)
413: {
414:   PetscFunctionBegin;
417:   PetscAssertPointer(val, 3);
420:   PetscCheckSameTypeAndComm(x, 1, y, 2);
421:   VecCheckSameSize(x, 1, y, 2);
422:   VecCheckAssembled(x);
423:   VecCheckAssembled(y);

425:   PetscCall(VecLockReadPush(x));
426:   PetscCall(VecLockReadPush(y));
427:   PetscCall(PetscLogEventBegin(VEC_TDot, x, y, 0, 0));
428:   PetscUseTypeMethod(x, tdot, y, val);
429:   PetscCall(PetscLogEventEnd(VEC_TDot, x, y, 0, 0));
430:   PetscCall(VecLockReadPop(x));
431:   PetscCall(VecLockReadPop(y));
432:   PetscFunctionReturn(PETSC_SUCCESS);
433: }

435: PetscErrorCode VecScaleAsync_Private(Vec x, PetscScalar alpha, PetscDeviceContext dctx)
436: {
437:   PetscReal   norms[4];
438:   PetscBool   flgs[4];
439:   PetscScalar one = 1.0;

441:   PetscFunctionBegin;
444:   VecCheckAssembled(x);
445:   PetscCall(VecSetErrorIfLocked(x, 1));
447:   if (alpha == one) PetscFunctionReturn(PETSC_SUCCESS);

449:   /* get current stashed norms */
450:   for (PetscInt i = 0; i < 4; i++) PetscCall(PetscObjectComposedDataGetReal((PetscObject)x, NormIds[i], norms[i], flgs[i]));

452:   PetscCall(PetscLogEventBegin(VEC_Scale, x, 0, 0, 0));
453:   VecMethodDispatch(x, dctx, VecAsyncFnName(Scale), scale, (Vec, PetscScalar, PetscDeviceContext), alpha);
454:   PetscCall(PetscLogEventEnd(VEC_Scale, x, 0, 0, 0));

456:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
457:   /* put the scaled stashed norms back into the Vec */
458:   for (PetscInt i = 0; i < 4; i++) {
459:     PetscReal ar = PetscAbsScalar(alpha);
460:     if (flgs[i]) PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[i], ar * norms[i]));
461:   }
462:   PetscFunctionReturn(PETSC_SUCCESS);
463: }

465: /*@
466:   VecScale - Scales a vector.

468:   Logically Collective

470:   Input Parameters:
471: + x     - the vector
472: - alpha - the scalar

474:   Level: intermediate

476:   Note:
477:   For a vector with n components, `VecScale()` computes  x[i] = alpha * x[i], for i=1,...,n.

479: .seealso: [](ch_vectors), `Vec`, `VecSet()`
480: @*/
481: PetscErrorCode VecScale(Vec x, PetscScalar alpha)
482: {
483:   PetscFunctionBegin;
484:   PetscCall(VecScaleAsync_Private(x, alpha, NULL));
485:   PetscFunctionReturn(PETSC_SUCCESS);
486: }

488: PetscErrorCode VecSetAsync_Private(Vec x, PetscScalar alpha, PetscDeviceContext dctx)
489: {
490:   PetscFunctionBegin;
493:   VecCheckAssembled(x);
495:   PetscCall(VecSetErrorIfLocked(x, 1));

497:   if (alpha == 0) {
498:     PetscReal norm;
499:     PetscBool set;

501:     PetscCall(VecNormAvailable(x, NORM_2, &set, &norm));
502:     if (set == PETSC_TRUE && norm == 0) PetscFunctionReturn(PETSC_SUCCESS);
503:   }
504:   PetscCall(PetscLogEventBegin(VEC_Set, x, 0, 0, 0));
505:   VecMethodDispatch(x, dctx, VecAsyncFnName(Set), set, (Vec, PetscScalar, PetscDeviceContext), alpha);
506:   PetscCall(PetscLogEventEnd(VEC_Set, x, 0, 0, 0));
507:   PetscCall(PetscObjectStateIncrease((PetscObject)x));

509:   /*  norms can be simply set (if |alpha|*N not too large) */
510:   {
511:     PetscReal      val = PetscAbsScalar(alpha);
512:     const PetscInt N   = x->map->N;

514:     if (N == 0) {
515:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_1], 0.0));
516:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_INFINITY], 0.0));
517:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_2], 0.0));
518:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_FROBENIUS], 0.0));
519:     } else if (val > PETSC_MAX_REAL / N) {
520:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_INFINITY], val));
521:     } else {
522:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_1], N * val));
523:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_INFINITY], val));
524:       val *= PetscSqrtReal((PetscReal)N);
525:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_2], val));
526:       PetscCall(PetscObjectComposedDataSetReal((PetscObject)x, NormIds[NORM_FROBENIUS], val));
527:     }
528:   }
529:   PetscFunctionReturn(PETSC_SUCCESS);
530: }

532: /*@
533:   VecSet - Sets all components of a vector to a single scalar value.

535:   Logically Collective

537:   Input Parameters:
538: + x     - the vector
539: - alpha - the scalar

541:   Level: beginner

543:   Notes:
544:   For a vector of dimension n, `VecSet()` sets x[i] = alpha, for i=1,...,n,
545:   so that all vector entries then equal the identical
546:   scalar value, `alpha`.  Use the more general routine
547:   `VecSetValues()` to set different vector entries.

549:   You CANNOT call this after you have called `VecSetValues()` but before you call
550:   `VecAssemblyBegin()`

552:   If `alpha` is zero and the norm of the vector is known to be zero then this skips the unneeded zeroing process

554: .seealso: [](ch_vectors), `Vec`, `VecSetValues()`, `VecSetValuesBlocked()`, `VecSetRandom()`
555: @*/
556: PetscErrorCode VecSet(Vec x, PetscScalar alpha)
557: {
558:   PetscFunctionBegin;
559:   PetscCall(VecSetAsync_Private(x, alpha, NULL));
560:   PetscFunctionReturn(PETSC_SUCCESS);
561: }

563: PetscErrorCode VecAXPYAsync_Private(Vec y, PetscScalar alpha, Vec x, PetscDeviceContext dctx)
564: {
565:   PetscFunctionBegin;
570:   PetscCheckSameTypeAndComm(x, 3, y, 1);
571:   VecCheckSameSize(x, 3, y, 1);
572:   VecCheckAssembled(x);
573:   VecCheckAssembled(y);
575:   if (alpha == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
576:   PetscCall(VecSetErrorIfLocked(y, 1));
577:   if (x == y) {
578:     PetscCall(VecScale(y, alpha + 1.0));
579:     PetscFunctionReturn(PETSC_SUCCESS);
580:   }
581:   PetscCall(VecLockReadPush(x));
582:   PetscCall(PetscLogEventBegin(VEC_AXPY, x, y, 0, 0));
583:   VecMethodDispatch(y, dctx, VecAsyncFnName(AXPY), axpy, (Vec, PetscScalar, Vec, PetscDeviceContext), alpha, x);
584:   PetscCall(PetscLogEventEnd(VEC_AXPY, x, y, 0, 0));
585:   PetscCall(VecLockReadPop(x));
586:   PetscCall(PetscObjectStateIncrease((PetscObject)y));
587:   PetscFunctionReturn(PETSC_SUCCESS);
588: }
589: /*@
590:   VecAXPY - Computes `y = alpha x + y`.

592:   Logically Collective

594:   Input Parameters:
595: + alpha - the scalar
596: . x     - vector scale by `alpha`
597: - y     - vector accumulated into

599:   Output Parameter:
600: . y - output vector

602:   Level: intermediate

604:   Notes:
605:   This routine is optimized for alpha of 0.0, otherwise it calls the BLAS routine
606: .vb
607:     VecAXPY(y,alpha,x)                   y = alpha x           +      y
608:     VecAYPX(y,beta,x)                    y =       x           + beta y
609:     VecAXPBY(y,alpha,beta,x)             y = alpha x           + beta y
610:     VecWAXPY(w,alpha,x,y)                w = alpha x           +      y
611:     VecAXPBYPCZ(z,alpha,beta,gamma,x,y)  z = alpha x           + beta y + gamma z
612:     VecMAXPY(y,nv,alpha[],x[])           y = sum alpha[i] x[i] +      y
613: .ve

615: .seealso: [](ch_vectors), `Vec`, `VecAYPX()`, `VecMAXPY()`, `VecWAXPY()`, `VecAXPBYPCZ()`, `VecAXPBY()`
616: @*/
617: PetscErrorCode VecAXPY(Vec y, PetscScalar alpha, Vec x)
618: {
619:   PetscFunctionBegin;
620:   PetscCall(VecAXPYAsync_Private(y, alpha, x, NULL));
621:   PetscFunctionReturn(PETSC_SUCCESS);
622: }

624: PetscErrorCode VecAYPXAsync_Private(Vec y, PetscScalar beta, Vec x, PetscDeviceContext dctx)
625: {
626:   PetscFunctionBegin;
631:   PetscCheckSameTypeAndComm(x, 3, y, 1);
632:   VecCheckSameSize(x, 1, y, 3);
633:   VecCheckAssembled(x);
634:   VecCheckAssembled(y);
636:   PetscCall(VecSetErrorIfLocked(y, 1));
637:   if (x == y) {
638:     PetscCall(VecScale(y, beta + 1.0));
639:     PetscFunctionReturn(PETSC_SUCCESS);
640:   }
641:   PetscCall(VecLockReadPush(x));
642:   if (beta == (PetscScalar)0.0) {
643:     PetscCall(VecCopy(x, y));
644:   } else {
645:     PetscCall(PetscLogEventBegin(VEC_AYPX, x, y, 0, 0));
646:     VecMethodDispatch(y, dctx, VecAsyncFnName(AYPX), aypx, (Vec, PetscScalar, Vec, PetscDeviceContext), beta, x);
647:     PetscCall(PetscLogEventEnd(VEC_AYPX, x, y, 0, 0));
648:     PetscCall(PetscObjectStateIncrease((PetscObject)y));
649:   }
650:   PetscCall(VecLockReadPop(x));
651:   PetscFunctionReturn(PETSC_SUCCESS);
652: }

654: /*@
655:   VecAYPX - Computes `y = x + beta y`.

657:   Logically Collective

659:   Input Parameters:
660: + beta - the scalar
661: . x    - the unscaled vector
662: - y    - the vector to be scaled

664:   Output Parameter:
665: . y - output vector

667:   Level: intermediate

669:   Developer Notes:
670:   The implementation is optimized for `beta` of -1.0, 0.0, and 1.0

672: .seealso: [](ch_vectors), `Vec`, `VecMAXPY()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBYPCZ()`, `VecAXPBY()`
673: @*/
674: PetscErrorCode VecAYPX(Vec y, PetscScalar beta, Vec x)
675: {
676:   PetscFunctionBegin;
677:   PetscCall(VecAYPXAsync_Private(y, beta, x, NULL));
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: PetscErrorCode VecAXPBYAsync_Private(Vec y, PetscScalar alpha, PetscScalar beta, Vec x, PetscDeviceContext dctx)
682: {
683:   PetscFunctionBegin;
688:   PetscCheckSameTypeAndComm(x, 4, y, 1);
689:   VecCheckSameSize(y, 1, x, 4);
690:   VecCheckAssembled(x);
691:   VecCheckAssembled(y);
694:   if (alpha == (PetscScalar)0.0 && beta == (PetscScalar)1.0) PetscFunctionReturn(PETSC_SUCCESS);
695:   if (x == y) {
696:     PetscCall(VecScale(y, alpha + beta));
697:     PetscFunctionReturn(PETSC_SUCCESS);
698:   }

700:   PetscCall(VecSetErrorIfLocked(y, 1));
701:   PetscCall(VecLockReadPush(x));
702:   PetscCall(PetscLogEventBegin(VEC_AXPY, y, x, 0, 0));
703:   VecMethodDispatch(y, dctx, VecAsyncFnName(AXPBY), axpby, (Vec, PetscScalar, PetscScalar, Vec, PetscDeviceContext), alpha, beta, x);
704:   PetscCall(PetscLogEventEnd(VEC_AXPY, y, x, 0, 0));
705:   PetscCall(PetscObjectStateIncrease((PetscObject)y));
706:   PetscCall(VecLockReadPop(x));
707:   PetscFunctionReturn(PETSC_SUCCESS);
708: }

710: /*@
711:   VecAXPBY - Computes `y = alpha x + beta y`.

713:   Logically Collective

715:   Input Parameters:
716: + alpha - first scalar
717: . beta  - second scalar
718: . x     - the first scaled vector
719: - y     - the second scaled vector

721:   Output Parameter:
722: . y - output vector

724:   Level: intermediate

726:   Developer Notes:
727:   The implementation is optimized for `alpha` and/or `beta` values of 0.0 and 1.0

729: .seealso: [](ch_vectors), `Vec`, `VecAYPX()`, `VecMAXPY()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBYPCZ()`
730: @*/
731: PetscErrorCode VecAXPBY(Vec y, PetscScalar alpha, PetscScalar beta, Vec x)
732: {
733:   PetscFunctionBegin;
734:   PetscCall(VecAXPBYAsync_Private(y, alpha, beta, x, NULL));
735:   PetscFunctionReturn(PETSC_SUCCESS);
736: }

738: PetscErrorCode VecAXPBYPCZAsync_Private(Vec z, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec x, Vec y, PetscDeviceContext dctx)
739: {
740:   PetscFunctionBegin;
747:   PetscCheckSameTypeAndComm(x, 5, y, 6);
748:   PetscCheckSameTypeAndComm(x, 5, z, 1);
749:   VecCheckSameSize(x, 5, y, 6);
750:   VecCheckSameSize(x, 5, z, 1);
751:   PetscCheck(x != y && x != z, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_IDN, "x, y, and z must be different vectors");
752:   PetscCheck(y != z, PetscObjectComm((PetscObject)y), PETSC_ERR_ARG_IDN, "x, y, and z must be different vectors");
753:   VecCheckAssembled(x);
754:   VecCheckAssembled(y);
755:   VecCheckAssembled(z);
759:   if (alpha == (PetscScalar)0.0 && beta == (PetscScalar)0.0 && gamma == (PetscScalar)1.0) PetscFunctionReturn(PETSC_SUCCESS);

761:   PetscCall(VecSetErrorIfLocked(z, 1));
762:   PetscCall(VecLockReadPush(x));
763:   PetscCall(VecLockReadPush(y));
764:   PetscCall(PetscLogEventBegin(VEC_AXPBYPCZ, x, y, z, 0));
765:   VecMethodDispatch(z, dctx, VecAsyncFnName(AXPBYPCZ), axpbypcz, (Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec, PetscDeviceContext), alpha, beta, gamma, x, y);
766:   PetscCall(PetscLogEventEnd(VEC_AXPBYPCZ, x, y, z, 0));
767:   PetscCall(PetscObjectStateIncrease((PetscObject)z));
768:   PetscCall(VecLockReadPop(x));
769:   PetscCall(VecLockReadPop(y));
770:   PetscFunctionReturn(PETSC_SUCCESS);
771: }
772: /*@
773:   VecAXPBYPCZ - Computes `z = alpha x + beta y + gamma z`

775:   Logically Collective

777:   Input Parameters:
778: + alpha - first scalar
779: . beta  - second scalar
780: . gamma - third scalar
781: . x     - first vector
782: . y     - second vector
783: - z     - third vector

785:   Output Parameter:
786: . z - output vector

788:   Level: intermediate

790:   Note:
791:   `x`, `y` and `z` must be different vectors

793:   Developer Notes:
794:   The implementation is optimized for `alpha` of 1.0 and `gamma` of 1.0 or 0.0

796: .seealso: [](ch_vectors), `Vec`, `VecAYPX()`, `VecMAXPY()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBY()`
797: @*/
798: PetscErrorCode VecAXPBYPCZ(Vec z, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec x, Vec y)
799: {
800:   PetscFunctionBegin;
801:   PetscCall(VecAXPBYPCZAsync_Private(z, alpha, beta, gamma, x, y, NULL));
802:   PetscFunctionReturn(PETSC_SUCCESS);
803: }

805: PetscErrorCode VecWAXPYAsync_Private(Vec w, PetscScalar alpha, Vec x, Vec y, PetscDeviceContext dctx)
806: {
807:   PetscFunctionBegin;
814:   PetscCheckSameTypeAndComm(x, 3, y, 4);
815:   PetscCheckSameTypeAndComm(y, 4, w, 1);
816:   VecCheckSameSize(x, 3, y, 4);
817:   VecCheckSameSize(x, 3, w, 1);
818:   PetscCheck(w != y, PETSC_COMM_SELF, PETSC_ERR_SUP, "Result vector w cannot be same as input vector y, suggest VecAXPY()");
819:   PetscCheck(w != x, PETSC_COMM_SELF, PETSC_ERR_SUP, "Result vector w cannot be same as input vector x, suggest VecAYPX()");
820:   VecCheckAssembled(x);
821:   VecCheckAssembled(y);
823:   PetscCall(VecSetErrorIfLocked(w, 1));

825:   PetscCall(VecLockReadPush(x));
826:   PetscCall(VecLockReadPush(y));
827:   if (alpha == (PetscScalar)0.0) {
828:     PetscCall(VecCopyAsync_Private(y, w, dctx));
829:   } else {
830:     PetscCall(PetscLogEventBegin(VEC_WAXPY, x, y, w, 0));
831:     VecMethodDispatch(w, dctx, VecAsyncFnName(WAXPY), waxpy, (Vec, PetscScalar, Vec, Vec, PetscDeviceContext), alpha, x, y);
832:     PetscCall(PetscLogEventEnd(VEC_WAXPY, x, y, w, 0));
833:     PetscCall(PetscObjectStateIncrease((PetscObject)w));
834:   }
835:   PetscCall(VecLockReadPop(x));
836:   PetscCall(VecLockReadPop(y));
837:   PetscFunctionReturn(PETSC_SUCCESS);
838: }

840: /*@
841:   VecWAXPY - Computes `w = alpha x + y`.

843:   Logically Collective

845:   Input Parameters:
846: + alpha - the scalar
847: . x     - first vector, multiplied by `alpha`
848: - y     - second vector

850:   Output Parameter:
851: . w - the result

853:   Level: intermediate

855:   Note:
856:   `w` cannot be either `x` or `y`, but `x` and `y` can be the same

858:   Developer Notes:
859:   The implementation is optimized for alpha of -1.0, 0.0, and 1.0

861: .seealso: [](ch_vectors), `Vec`, `VecAXPY()`, `VecAYPX()`, `VecAXPBY()`, `VecMAXPY()`, `VecAXPBYPCZ()`
862: @*/
863: PetscErrorCode VecWAXPY(Vec w, PetscScalar alpha, Vec x, Vec y)
864: {
865:   PetscFunctionBegin;
866:   PetscCall(VecWAXPYAsync_Private(w, alpha, x, y, NULL));
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: /*@C
871:   VecSetValues - Inserts or adds values into certain locations of a vector.

873:   Not Collective

875:   Input Parameters:
876: + x    - vector to insert in
877: . ni   - number of elements to add
878: . ix   - indices where to add
879: . y    - array of values
880: - iora - either `INSERT_VALUES` to replace the current values or `ADD_VALUES` to add values to any existing entries

882:   Level: beginner

884:   Notes:
885: .vb
886:    `VecSetValues()` sets x[ix[i]] = y[i], for i=0,...,ni-1.
887: .ve

889:   Calls to `VecSetValues()` with the `INSERT_VALUES` and `ADD_VALUES`
890:   options cannot be mixed without intervening calls to the assembly
891:   routines.

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

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

898:   If you call `VecSetOption`(x, `VEC_IGNORE_NEGATIVE_INDICES`,`PETSC_TRUE`),
899:   negative indices may be passed in ix. These rows are
900:   simply ignored. This allows easily inserting element load matrices
901:   with homogeneous Dirichlet boundary conditions that you don't want represented
902:   in the vector.

904: .seealso: [](ch_vectors), `Vec`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesLocal()`,
905:           `VecSetValue()`, `VecSetValuesBlocked()`, `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `VecGetValues()`
906: @*/
907: PetscErrorCode VecSetValues(Vec x, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode iora)
908: {
909:   PetscFunctionBeginHot;
911:   if (!ni) PetscFunctionReturn(PETSC_SUCCESS);
912:   PetscAssertPointer(ix, 3);
913:   PetscAssertPointer(y, 4);

916:   PetscCall(PetscLogEventBegin(VEC_SetValues, x, 0, 0, 0));
917:   PetscUseTypeMethod(x, setvalues, ni, ix, y, iora);
918:   PetscCall(PetscLogEventEnd(VEC_SetValues, x, 0, 0, 0));
919:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
920:   PetscFunctionReturn(PETSC_SUCCESS);
921: }

923: /*@C
924:   VecGetValues - Gets values from certain locations of a vector. Currently
925:   can only get values on the same processor on which they are owned

927:   Not Collective

929:   Input Parameters:
930: + x  - vector to get values from
931: . ni - number of elements to get
932: - ix - indices where to get them from (in global 1d numbering)

934:   Output Parameter:
935: . y - array of values

937:   Level: beginner

939:   Notes:
940:   The user provides the allocated array y; it is NOT allocated in this routine

942:   `VecGetValues()` gets y[i] = x[ix[i]], for i=0,...,ni-1.

944:   `VecAssemblyBegin()` and `VecAssemblyEnd()`  MUST be called before calling this if `VecSetValues()` or related routine has been called

946:   VecGetValues() uses 0-based indices in Fortran as well as in C.

948:   If you call `VecSetOption`(x, `VEC_IGNORE_NEGATIVE_INDICES`,`PETSC_TRUE`),
949:   negative indices may be passed in ix. These rows are
950:   simply ignored.

952: .seealso: [](ch_vectors), `Vec`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValues()`
953: @*/
954: PetscErrorCode VecGetValues(Vec x, PetscInt ni, const PetscInt ix[], PetscScalar y[])
955: {
956:   PetscFunctionBegin;
958:   if (!ni) PetscFunctionReturn(PETSC_SUCCESS);
959:   PetscAssertPointer(ix, 3);
960:   PetscAssertPointer(y, 4);
962:   VecCheckAssembled(x);
963:   PetscUseTypeMethod(x, getvalues, ni, ix, y);
964:   PetscFunctionReturn(PETSC_SUCCESS);
965: }

967: /*@C
968:   VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector.

970:   Not Collective

972:   Input Parameters:
973: + x    - vector to insert in
974: . ni   - number of blocks to add
975: . ix   - indices where to add in block count, rather than element count
976: . y    - array of values
977: - iora - either `INSERT_VALUES` replaces existing entries with new values, `ADD_VALUES`, adds values to any existing entries

979:   Level: intermediate

981:   Notes:
982:   `VecSetValuesBlocked()` sets x[bs*ix[i]+j] = y[bs*i+j],
983:   for j=0,...,bs-1, for i=0,...,ni-1. where bs was set with VecSetBlockSize().

985:   Calls to `VecSetValuesBlocked()` with the `INSERT_VALUES` and `ADD_VALUES`
986:   options cannot be mixed without intervening calls to the assembly
987:   routines.

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

992:   `VecSetValuesBlocked()` uses 0-based indices in Fortran as well as in C.

994:   Negative indices may be passed in ix, these rows are
995:   simply ignored. This allows easily inserting element load matrices
996:   with homogeneous Dirichlet boundary conditions that you don't want represented
997:   in the vector.

999: .seealso: [](ch_vectors), `Vec`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValuesBlockedLocal()`,
1000:           `VecSetValues()`
1001: @*/
1002: PetscErrorCode VecSetValuesBlocked(Vec x, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode iora)
1003: {
1004:   PetscFunctionBeginHot;
1006:   if (!ni) PetscFunctionReturn(PETSC_SUCCESS);
1007:   PetscAssertPointer(ix, 3);
1008:   PetscAssertPointer(y, 4);

1011:   PetscCall(PetscLogEventBegin(VEC_SetValues, x, 0, 0, 0));
1012:   PetscUseTypeMethod(x, setvaluesblocked, ni, ix, y, iora);
1013:   PetscCall(PetscLogEventEnd(VEC_SetValues, x, 0, 0, 0));
1014:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
1015:   PetscFunctionReturn(PETSC_SUCCESS);
1016: }

1018: /*@C
1019:   VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
1020:   using a local ordering of the nodes.

1022:   Not Collective

1024:   Input Parameters:
1025: + x    - vector to insert in
1026: . ni   - number of elements to add
1027: . ix   - indices where to add
1028: . y    - array of values
1029: - iora - either `INSERT_VALUES` replaces existing entries with new values, `ADD_VALUES` adds values to any existing entries

1031:   Level: intermediate

1033:   Notes:
1034:   `VecSetValuesLocal()` sets x[ix[i]] = y[i], for i=0,...,ni-1.

1036:   Calls to `VecSetValues()` with the `INSERT_VALUES` and `ADD_VALUES`
1037:   options cannot be mixed without intervening calls to the assembly
1038:   routines.

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

1043:   `VecSetValuesLocal()` uses 0-based indices in Fortran as well as in C.

1045: .seealso: [](ch_vectors), `Vec`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValues()`, `VecSetLocalToGlobalMapping()`,
1046:           `VecSetValuesBlockedLocal()`
1047: @*/
1048: PetscErrorCode VecSetValuesLocal(Vec x, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode iora)
1049: {
1050:   PetscInt lixp[128], *lix = lixp;

1052:   PetscFunctionBeginHot;
1054:   if (!ni) PetscFunctionReturn(PETSC_SUCCESS);
1055:   PetscAssertPointer(ix, 3);
1056:   PetscAssertPointer(y, 4);

1059:   PetscCall(PetscLogEventBegin(VEC_SetValues, x, 0, 0, 0));
1060:   if (!x->ops->setvalueslocal) {
1061:     if (x->map->mapping) {
1062:       if (ni > 128) PetscCall(PetscMalloc1(ni, &lix));
1063:       PetscCall(ISLocalToGlobalMappingApply(x->map->mapping, ni, (PetscInt *)ix, lix));
1064:       PetscUseTypeMethod(x, setvalues, ni, lix, y, iora);
1065:       if (ni > 128) PetscCall(PetscFree(lix));
1066:     } else PetscUseTypeMethod(x, setvalues, ni, ix, y, iora);
1067:   } else PetscUseTypeMethod(x, setvalueslocal, ni, ix, y, iora);
1068:   PetscCall(PetscLogEventEnd(VEC_SetValues, x, 0, 0, 0));
1069:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
1070:   PetscFunctionReturn(PETSC_SUCCESS);
1071: }

1073: /*@
1074:   VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1075:   using a local ordering of the nodes.

1077:   Not Collective

1079:   Input Parameters:
1080: + x    - vector to insert in
1081: . ni   - number of blocks to add
1082: . ix   - indices where to add in block count, not element count
1083: . y    - array of values
1084: - iora - either `INSERT_VALUES` replaces existing entries with new values, `ADD_VALUES` adds values to any existing entries

1086:   Level: intermediate

1088:   Notes:
1089:   `VecSetValuesBlockedLocal()` sets x[bs*ix[i]+j] = y[bs*i+j],
1090:   for j=0,..bs-1, for i=0,...,ni-1, where bs has been set with `VecSetBlockSize()`.

1092:   Calls to `VecSetValuesBlockedLocal()` with the `INSERT_VALUES` and `ADD_VALUES`
1093:   options cannot be mixed without intervening calls to the assembly
1094:   routines.

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

1099:   `VecSetValuesBlockedLocal()` uses 0-based indices in Fortran as well as in C.

1101: .seealso: [](ch_vectors), `Vec`, `VecAssemblyBegin()`, `VecAssemblyEnd()`, `VecSetValues()`, `VecSetValuesBlocked()`,
1102:           `VecSetLocalToGlobalMapping()`
1103: @*/
1104: PetscErrorCode VecSetValuesBlockedLocal(Vec x, PetscInt ni, const PetscInt ix[], const PetscScalar y[], InsertMode iora)
1105: {
1106:   PetscInt lixp[128], *lix = lixp;

1108:   PetscFunctionBeginHot;
1110:   if (!ni) PetscFunctionReturn(PETSC_SUCCESS);
1111:   PetscAssertPointer(ix, 3);
1112:   PetscAssertPointer(y, 4);
1114:   PetscCall(PetscLogEventBegin(VEC_SetValues, x, 0, 0, 0));
1115:   if (x->map->mapping) {
1116:     if (ni > 128) PetscCall(PetscMalloc1(ni, &lix));
1117:     PetscCall(ISLocalToGlobalMappingApplyBlock(x->map->mapping, ni, (PetscInt *)ix, lix));
1118:     PetscUseTypeMethod(x, setvaluesblocked, ni, lix, y, iora);
1119:     if (ni > 128) PetscCall(PetscFree(lix));
1120:   } else {
1121:     PetscUseTypeMethod(x, setvaluesblocked, ni, ix, y, iora);
1122:   }
1123:   PetscCall(PetscLogEventEnd(VEC_SetValues, x, 0, 0, 0));
1124:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
1125:   PetscFunctionReturn(PETSC_SUCCESS);
1126: }

1128: static PetscErrorCode VecMXDot_Private(Vec x, PetscInt nv, const Vec y[], PetscScalar result[], PetscErrorCode (*mxdot)(Vec, PetscInt, const Vec[], PetscScalar[]), PetscLogEvent event)
1129: {
1130:   PetscFunctionBegin;
1133:   VecCheckAssembled(x);
1135:   if (!nv) PetscFunctionReturn(PETSC_SUCCESS);
1136:   PetscAssertPointer(y, 3);
1137:   for (PetscInt i = 0; i < nv; ++i) {
1140:     PetscCheckSameTypeAndComm(x, 1, y[i], 3);
1141:     VecCheckSameSize(x, 1, y[i], 3);
1142:     VecCheckAssembled(y[i]);
1143:     PetscCall(VecLockReadPush(y[i]));
1144:   }
1145:   PetscAssertPointer(result, 4);

1148:   PetscCall(VecLockReadPush(x));
1149:   PetscCall(PetscLogEventBegin(event, x, *y, 0, 0));
1150:   PetscCall((*mxdot)(x, nv, y, result));
1151:   PetscCall(PetscLogEventEnd(event, x, *y, 0, 0));
1152:   PetscCall(VecLockReadPop(x));
1153:   for (PetscInt i = 0; i < nv; ++i) PetscCall(VecLockReadPop(y[i]));
1154:   PetscFunctionReturn(PETSC_SUCCESS);
1155: }

1157: /*@
1158:   VecMTDot - Computes indefinite vector multiple dot products.
1159:   That is, it does NOT use the complex conjugate.

1161:   Collective

1163:   Input Parameters:
1164: + x  - one vector
1165: . nv - number of vectors
1166: - y  - array of vectors.  Note that vectors are pointers

1168:   Output Parameter:
1169: . val - array of the dot products

1171:   Level: intermediate

1173:   Notes for Users of Complex Numbers:
1174:   For complex vectors, `VecMTDot()` computes the indefinite form
1175: $      val = (x,y) = y^T x,
1176:   where y^T denotes the transpose of y.

1178:   Use `VecMDot()` for the inner product
1179: $      val = (x,y) = y^H x,
1180:   where y^H denotes the conjugate transpose of y.

1182: .seealso: [](ch_vectors), `Vec`, `VecMDot()`, `VecTDot()`
1183: @*/
1184: PetscErrorCode VecMTDot(Vec x, PetscInt nv, const Vec y[], PetscScalar val[])
1185: {
1186:   PetscFunctionBegin;
1188:   PetscCall(VecMXDot_Private(x, nv, y, val, x->ops->mtdot, VEC_MTDot));
1189:   PetscFunctionReturn(PETSC_SUCCESS);
1190: }

1192: /*@
1193:   VecMDot - Computes multiple vector dot products.

1195:   Collective

1197:   Input Parameters:
1198: + x  - one vector
1199: . nv - number of vectors
1200: - y  - array of vectors.

1202:   Output Parameter:
1203: . val - array of the dot products (does not allocate the array)

1205:   Level: intermediate

1207:   Notes for Users of Complex Numbers:
1208:   For complex vectors, `VecMDot()` computes
1209: $     val = (x,y) = y^H x,
1210:   where y^H denotes the conjugate transpose of y.

1212:   Use `VecMTDot()` for the indefinite form
1213: $     val = (x,y) = y^T x,
1214:   where y^T denotes the transpose of y.

1216: .seealso: [](ch_vectors), `Vec`, `VecMTDot()`, `VecDot()`
1217: @*/
1218: PetscErrorCode VecMDot(Vec x, PetscInt nv, const Vec y[], PetscScalar val[])
1219: {
1220:   PetscFunctionBegin;
1222:   PetscCall(VecMXDot_Private(x, nv, y, val, x->ops->mdot, VEC_MDot));
1223:   PetscFunctionReturn(PETSC_SUCCESS);
1224: }

1226: PetscErrorCode VecMAXPYAsync_Private(Vec y, PetscInt nv, const PetscScalar alpha[], Vec x[], PetscDeviceContext dctx)
1227: {
1228:   PetscFunctionBegin;
1230:   VecCheckAssembled(y);
1232:   PetscCall(VecSetErrorIfLocked(y, 1));
1233:   PetscCheck(nv >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vectors (given %" PetscInt_FMT ") cannot be negative", nv);
1234:   if (nv) {
1235:     PetscInt zeros = 0;

1237:     PetscAssertPointer(alpha, 3);
1238:     PetscAssertPointer(x, 4);
1239:     for (PetscInt i = 0; i < nv; ++i) {
1243:       PetscCheckSameTypeAndComm(y, 1, x[i], 4);
1244:       VecCheckSameSize(y, 1, x[i], 4);
1245:       PetscCheck(y != x[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Array of vectors 'x' cannot contain y, found x[%" PetscInt_FMT "] == y", i);
1246:       VecCheckAssembled(x[i]);
1247:       PetscCall(VecLockReadPush(x[i]));
1248:       zeros += alpha[i] == (PetscScalar)0.0;
1249:     }

1251:     if (zeros < nv) {
1252:       PetscCall(PetscLogEventBegin(VEC_MAXPY, y, *x, 0, 0));
1253:       VecMethodDispatch(y, dctx, VecAsyncFnName(MAXPY), maxpy, (Vec, PetscInt, const PetscScalar[], Vec[], PetscDeviceContext), nv, alpha, x);
1254:       PetscCall(PetscLogEventEnd(VEC_MAXPY, y, *x, 0, 0));
1255:       PetscCall(PetscObjectStateIncrease((PetscObject)y));
1256:     }

1258:     for (PetscInt i = 0; i < nv; ++i) PetscCall(VecLockReadPop(x[i]));
1259:   }
1260:   PetscFunctionReturn(PETSC_SUCCESS);
1261: }

1263: /*@
1264:   VecMAXPY - Computes `y = y + sum alpha[i] x[i]`

1266:   Logically Collective

1268:   Input Parameters:
1269: + nv    - number of scalars and x-vectors
1270: . alpha - array of scalars
1271: . y     - one vector
1272: - x     - array of vectors

1274:   Level: intermediate

1276:   Note:
1277:   `y` cannot be any of the `x` vectors

1279: .seealso: [](ch_vectors), `Vec`, `VecMAXPBY()`,`VecAYPX()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBYPCZ()`, `VecAXPBY()`
1280: @*/
1281: PetscErrorCode VecMAXPY(Vec y, PetscInt nv, const PetscScalar alpha[], Vec x[])
1282: {
1283:   PetscFunctionBegin;
1284:   PetscCall(VecMAXPYAsync_Private(y, nv, alpha, x, NULL));
1285:   PetscFunctionReturn(PETSC_SUCCESS);
1286: }

1288: /*@
1289:   VecMAXPBY - Computes `y = beta y + sum alpha[i] x[i]`

1291:   Logically Collective

1293:   Input Parameters:
1294: + nv    - number of scalars and x-vectors
1295: . alpha - array of scalars
1296: . beta  - scalar
1297: . y     - one vector
1298: - x     - array of vectors

1300:   Level: intermediate

1302:   Note:
1303:   `y` cannot be any of the `x` vectors.

1305:   Developer Notes:
1306:   This is a convenience routine, but implementations might be able to optimize it, for example, when `beta` is zero.

1308: .seealso: [](ch_vectors), `Vec`, `VecMAXPY()`, `VecAYPX()`, `VecWAXPY()`, `VecAXPY()`, `VecAXPBYPCZ()`, `VecAXPBY()`
1309: @*/
1310: PetscErrorCode VecMAXPBY(Vec y, PetscInt nv, const PetscScalar alpha[], PetscScalar beta, Vec x[])
1311: {
1312:   PetscFunctionBegin;
1314:   VecCheckAssembled(y);
1316:   PetscCall(VecSetErrorIfLocked(y, 1));
1317:   PetscCheck(nv >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vectors (given %" PetscInt_FMT ") cannot be negative", nv);

1320:   if (y->ops->maxpby) {
1321:     PetscInt zeros = 0;

1323:     if (nv) {
1324:       PetscAssertPointer(alpha, 3);
1325:       PetscAssertPointer(x, 5);
1326:     }

1328:     for (PetscInt i = 0; i < nv; ++i) { // scan all alpha[]
1332:       PetscCheckSameTypeAndComm(y, 1, x[i], 5);
1333:       VecCheckSameSize(y, 1, x[i], 5);
1334:       PetscCheck(y != x[i], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Array of vectors 'x' cannot contain y, found x[%" PetscInt_FMT "] == y", i);
1335:       VecCheckAssembled(x[i]);
1336:       PetscCall(VecLockReadPush(x[i]));
1337:       zeros += alpha[i] == (PetscScalar)0.0;
1338:     }

1340:     if (zeros < nv) { // has nonzero alpha
1341:       PetscCall(PetscLogEventBegin(VEC_MAXPY, y, *x, 0, 0));
1342:       PetscUseTypeMethod(y, maxpby, nv, alpha, beta, x);
1343:       PetscCall(PetscLogEventEnd(VEC_MAXPY, y, *x, 0, 0));
1344:       PetscCall(PetscObjectStateIncrease((PetscObject)y));
1345:     } else {
1346:       PetscCall(VecScale(y, beta));
1347:     }

1349:     for (PetscInt i = 0; i < nv; ++i) PetscCall(VecLockReadPop(x[i]));
1350:   } else { // no maxpby
1351:     if (beta == 0.0) PetscCall(VecSet(y, 0.0));
1352:     else PetscCall(VecScale(y, beta));
1353:     PetscCall(VecMAXPY(y, nv, alpha, x));
1354:   }
1355:   PetscFunctionReturn(PETSC_SUCCESS);
1356: }

1358: /*@
1359:   VecConcatenate - Creates a new vector that is a vertical concatenation of all the given array of vectors
1360:   in the order they appear in the array. The concatenated vector resides on the same
1361:   communicator and is the same type as the source vectors.

1363:   Collective

1365:   Input Parameters:
1366: + nx - number of vectors to be concatenated
1367: - X  - array containing the vectors to be concatenated in the order of concatenation

1369:   Output Parameters:
1370: + Y    - concatenated vector
1371: - x_is - array of index sets corresponding to the concatenated components of `Y` (pass `NULL` if not needed)

1373:   Level: advanced

1375:   Notes:
1376:   Concatenation is similar to the functionality of a `VECNEST` object; they both represent combination of
1377:   different vector spaces. However, concatenated vectors do not store any information about their
1378:   sub-vectors and own their own data. Consequently, this function provides index sets to enable the
1379:   manipulation of data in the concatenated vector that corresponds to the original components at creation.

1381:   This is a useful tool for outer loop algorithms, particularly constrained optimizers, where the solver
1382:   has to operate on combined vector spaces and cannot utilize `VECNEST` objects due to incompatibility with
1383:   bound projections.

1385: .seealso: [](ch_vectors), `Vec`, `VECNEST`, `VECSCATTER`, `VecScatterCreate()`
1386: @*/
1387: PetscErrorCode VecConcatenate(PetscInt nx, const Vec X[], Vec *Y, IS *x_is[])
1388: {
1389:   MPI_Comm comm;
1390:   VecType  vec_type;
1391:   Vec      Ytmp, Xtmp;
1392:   IS      *is_tmp;
1393:   PetscInt i, shift = 0, Xnl, Xng, Xbegin;

1395:   PetscFunctionBegin;
1399:   PetscAssertPointer(Y, 3);

1401:   if ((*X)->ops->concatenate) {
1402:     /* use the dedicated concatenation function if available */
1403:     PetscCall((*(*X)->ops->concatenate)(nx, X, Y, x_is));
1404:   } else {
1405:     /* loop over vectors and start creating IS */
1406:     comm = PetscObjectComm((PetscObject)*X);
1407:     PetscCall(VecGetType(*X, &vec_type));
1408:     PetscCall(PetscMalloc1(nx, &is_tmp));
1409:     for (i = 0; i < nx; i++) {
1410:       PetscCall(VecGetSize(X[i], &Xng));
1411:       PetscCall(VecGetLocalSize(X[i], &Xnl));
1412:       PetscCall(VecGetOwnershipRange(X[i], &Xbegin, NULL));
1413:       PetscCall(ISCreateStride(comm, Xnl, shift + Xbegin, 1, &is_tmp[i]));
1414:       shift += Xng;
1415:     }
1416:     /* create the concatenated vector */
1417:     PetscCall(VecCreate(comm, &Ytmp));
1418:     PetscCall(VecSetType(Ytmp, vec_type));
1419:     PetscCall(VecSetSizes(Ytmp, PETSC_DECIDE, shift));
1420:     PetscCall(VecSetUp(Ytmp));
1421:     /* copy data from X array to Y and return */
1422:     for (i = 0; i < nx; i++) {
1423:       PetscCall(VecGetSubVector(Ytmp, is_tmp[i], &Xtmp));
1424:       PetscCall(VecCopy(X[i], Xtmp));
1425:       PetscCall(VecRestoreSubVector(Ytmp, is_tmp[i], &Xtmp));
1426:     }
1427:     *Y = Ytmp;
1428:     if (x_is) {
1429:       *x_is = is_tmp;
1430:     } else {
1431:       for (i = 0; i < nx; i++) PetscCall(ISDestroy(&is_tmp[i]));
1432:       PetscCall(PetscFree(is_tmp));
1433:     }
1434:   }
1435:   PetscFunctionReturn(PETSC_SUCCESS);
1436: }

1438: /* A helper function for VecGetSubVector to check if we can implement it with no-copy (i.e. the subvector shares
1439:    memory with the original vector), and the block size of the subvector.

1441:     Input Parameters:
1442: +   X - the original vector
1443: -   is - the index set of the subvector

1445:     Output Parameters:
1446: +   contig - PETSC_TRUE if the index set refers to contiguous entries on this process, else PETSC_FALSE
1447: .   start  - start of contiguous block, as an offset from the start of the ownership range of the original vector
1448: -   blocksize - the block size of the subvector

1450: */
1451: PetscErrorCode VecGetSubVectorContiguityAndBS_Private(Vec X, IS is, PetscBool *contig, PetscInt *start, PetscInt *blocksize)
1452: {
1453:   PetscInt  gstart, gend, lstart;
1454:   PetscBool red[2] = {PETSC_TRUE /*contiguous*/, PETSC_TRUE /*validVBS*/};
1455:   PetscInt  n, N, ibs, vbs, bs = -1;

1457:   PetscFunctionBegin;
1458:   PetscCall(ISGetLocalSize(is, &n));
1459:   PetscCall(ISGetSize(is, &N));
1460:   PetscCall(ISGetBlockSize(is, &ibs));
1461:   PetscCall(VecGetBlockSize(X, &vbs));
1462:   PetscCall(VecGetOwnershipRange(X, &gstart, &gend));
1463:   PetscCall(ISContiguousLocal(is, gstart, gend, &lstart, &red[0]));
1464:   /* block size is given by IS if ibs > 1; otherwise, check the vector */
1465:   if (ibs > 1) {
1466:     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, red, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)is)));
1467:     bs = ibs;
1468:   } else {
1469:     if (n % vbs || vbs == 1) red[1] = PETSC_FALSE; /* this process invalidate the collectiveness of block size */
1470:     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, red, 2, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)is)));
1471:     if (red[0] && red[1]) bs = vbs; /* all processes have a valid block size and the access will be contiguous */
1472:   }

1474:   *contig    = red[0];
1475:   *start     = lstart;
1476:   *blocksize = bs;
1477:   PetscFunctionReturn(PETSC_SUCCESS);
1478: }

1480: /* A helper function for VecGetSubVector, to be used when we have to build a standalone subvector through VecScatter

1482:     Input Parameters:
1483: +   X - the original vector
1484: .   is - the index set of the subvector
1485: -   bs - the block size of the subvector, gotten from VecGetSubVectorContiguityAndBS_Private()

1487:     Output Parameter:
1488: .   Z  - the subvector, which will compose the VecScatter context on output
1489: */
1490: PetscErrorCode VecGetSubVectorThroughVecScatter_Private(Vec X, IS is, PetscInt bs, Vec *Z)
1491: {
1492:   PetscInt   n, N;
1493:   VecScatter vscat;
1494:   Vec        Y;

1496:   PetscFunctionBegin;
1497:   PetscCall(ISGetLocalSize(is, &n));
1498:   PetscCall(ISGetSize(is, &N));
1499:   PetscCall(VecCreate(PetscObjectComm((PetscObject)is), &Y));
1500:   PetscCall(VecSetSizes(Y, n, N));
1501:   PetscCall(VecSetBlockSize(Y, bs));
1502:   PetscCall(VecSetType(Y, ((PetscObject)X)->type_name));
1503:   PetscCall(VecScatterCreate(X, is, Y, NULL, &vscat));
1504:   PetscCall(VecScatterBegin(vscat, X, Y, INSERT_VALUES, SCATTER_FORWARD));
1505:   PetscCall(VecScatterEnd(vscat, X, Y, INSERT_VALUES, SCATTER_FORWARD));
1506:   PetscCall(PetscObjectCompose((PetscObject)Y, "VecGetSubVector_Scatter", (PetscObject)vscat));
1507:   PetscCall(VecScatterDestroy(&vscat));
1508:   *Z = Y;
1509:   PetscFunctionReturn(PETSC_SUCCESS);
1510: }

1512: /*@
1513:   VecGetSubVector - Gets a vector representing part of another vector

1515:   Collective

1517:   Input Parameters:
1518: + X  - vector from which to extract a subvector
1519: - is - index set representing portion of `X` to extract

1521:   Output Parameter:
1522: . Y - subvector corresponding to `is`

1524:   Level: advanced

1526:   Notes:
1527:   The subvector `Y` should be returned with `VecRestoreSubVector()`.
1528:   `X` and must be defined on the same communicator

1530:   This function may return a subvector without making a copy, therefore it is not safe to use the original vector while
1531:   modifying the subvector.  Other non-overlapping subvectors can still be obtained from X using this function.

1533:   The resulting subvector inherits the block size from `is` if greater than one. Otherwise, the block size is guessed from the block size of the original `X`.

1535: .seealso: [](ch_vectors), `Vec`, `IS`, `VECNEST`, `MatCreateSubMatrix()`
1536: @*/
1537: PetscErrorCode VecGetSubVector(Vec X, IS is, Vec *Y)
1538: {
1539:   Vec Z;

1541:   PetscFunctionBegin;
1544:   PetscCheckSameComm(X, 1, is, 2);
1545:   PetscAssertPointer(Y, 3);
1546:   if (X->ops->getsubvector) {
1547:     PetscUseTypeMethod(X, getsubvector, is, &Z);
1548:   } else { /* Default implementation currently does no caching */
1549:     PetscBool contig;
1550:     PetscInt  n, N, start, bs;

1552:     PetscCall(ISGetLocalSize(is, &n));
1553:     PetscCall(ISGetSize(is, &N));
1554:     PetscCall(VecGetSubVectorContiguityAndBS_Private(X, is, &contig, &start, &bs));
1555:     if (contig) { /* We can do a no-copy implementation */
1556:       const PetscScalar *x;
1557:       PetscInt           state = 0;
1558:       PetscBool          isstd, iscuda, iship;

1560:       PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &isstd, VECSEQ, VECMPI, VECSTANDARD, ""));
1561:       PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &iscuda, VECSEQCUDA, VECMPICUDA, ""));
1562:       PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &iship, VECSEQHIP, VECMPIHIP, ""));
1563:       if (iscuda) {
1564: #if defined(PETSC_HAVE_CUDA)
1565:         const PetscScalar *x_d;
1566:         PetscMPIInt        size;
1567:         PetscOffloadMask   flg;

1569:         PetscCall(VecCUDAGetArrays_Private(X, &x, &x_d, &flg));
1570:         PetscCheck(flg != PETSC_OFFLOAD_UNALLOCATED, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for PETSC_OFFLOAD_UNALLOCATED");
1571:         PetscCheck(!n || x || x_d, PETSC_COMM_SELF, PETSC_ERR_SUP, "Missing vector data");
1572:         if (x) x += start;
1573:         if (x_d) x_d += start;
1574:         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)X), &size));
1575:         if (size == 1) {
1576:           PetscCall(VecCreateSeqCUDAWithArrays(PetscObjectComm((PetscObject)X), bs, n, x, x_d, &Z));
1577:         } else {
1578:           PetscCall(VecCreateMPICUDAWithArrays(PetscObjectComm((PetscObject)X), bs, n, N, x, x_d, &Z));
1579:         }
1580:         Z->offloadmask = flg;
1581: #endif
1582:       } else if (iship) {
1583: #if defined(PETSC_HAVE_HIP)
1584:         const PetscScalar *x_d;
1585:         PetscMPIInt        size;
1586:         PetscOffloadMask   flg;

1588:         PetscCall(VecHIPGetArrays_Private(X, &x, &x_d, &flg));
1589:         PetscCheck(flg != PETSC_OFFLOAD_UNALLOCATED, PETSC_COMM_SELF, PETSC_ERR_SUP, "Not for PETSC_OFFLOAD_UNALLOCATED");
1590:         PetscCheck(!n || x || x_d, PETSC_COMM_SELF, PETSC_ERR_SUP, "Missing vector data");
1591:         if (x) x += start;
1592:         if (x_d) x_d += start;
1593:         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)X), &size));
1594:         if (size == 1) {
1595:           PetscCall(VecCreateSeqHIPWithArrays(PetscObjectComm((PetscObject)X), bs, n, x, x_d, &Z));
1596:         } else {
1597:           PetscCall(VecCreateMPIHIPWithArrays(PetscObjectComm((PetscObject)X), bs, n, N, x, x_d, &Z));
1598:         }
1599:         Z->offloadmask = flg;
1600: #endif
1601:       } else if (isstd) {
1602:         PetscMPIInt size;

1604:         PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)X), &size));
1605:         PetscCall(VecGetArrayRead(X, &x));
1606:         if (x) x += start;
1607:         if (size == 1) {
1608:           PetscCall(VecCreateSeqWithArray(PetscObjectComm((PetscObject)X), bs, n, x, &Z));
1609:         } else {
1610:           PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)X), bs, n, N, x, &Z));
1611:         }
1612:         PetscCall(VecRestoreArrayRead(X, &x));
1613:       } else { /* default implementation: use place array */
1614:         PetscCall(VecGetArrayRead(X, &x));
1615:         PetscCall(VecCreate(PetscObjectComm((PetscObject)X), &Z));
1616:         PetscCall(VecSetType(Z, ((PetscObject)X)->type_name));
1617:         PetscCall(VecSetSizes(Z, n, N));
1618:         PetscCall(VecSetBlockSize(Z, bs));
1619:         PetscCall(VecPlaceArray(Z, PetscSafePointerPlusOffset(x, start)));
1620:         PetscCall(VecRestoreArrayRead(X, &x));
1621:       }

1623:       /* this is relevant only in debug mode */
1624:       PetscCall(VecLockGet(X, &state));
1625:       if (state) PetscCall(VecLockReadPush(Z));
1626:       Z->ops->placearray   = NULL;
1627:       Z->ops->replacearray = NULL;
1628:     } else { /* Have to create a scatter and do a copy */
1629:       PetscCall(VecGetSubVectorThroughVecScatter_Private(X, is, bs, &Z));
1630:     }
1631:   }
1632:   /* Record the state when the subvector was gotten so we know whether its values need to be put back */
1633:   if (VecGetSubVectorSavedStateId < 0) PetscCall(PetscObjectComposedDataRegister(&VecGetSubVectorSavedStateId));
1634:   PetscCall(PetscObjectComposedDataSetInt((PetscObject)Z, VecGetSubVectorSavedStateId, 1));
1635:   *Y = Z;
1636:   PetscFunctionReturn(PETSC_SUCCESS);
1637: }

1639: /*@
1640:   VecRestoreSubVector - Restores a subvector extracted using `VecGetSubVector()`

1642:   Collective

1644:   Input Parameters:
1645: + X  - vector from which subvector was obtained
1646: . is - index set representing the subset of `X`
1647: - Y  - subvector being restored

1649:   Level: advanced

1651: .seealso: [](ch_vectors), `Vec`, `IS`, `VecGetSubVector()`
1652: @*/
1653: PetscErrorCode VecRestoreSubVector(Vec X, IS is, Vec *Y)
1654: {
1655:   PETSC_UNUSED PetscObjectState dummystate = 0;
1656:   PetscBool                     unchanged;

1658:   PetscFunctionBegin;
1661:   PetscCheckSameComm(X, 1, is, 2);
1662:   PetscAssertPointer(Y, 3);

1665:   if (X->ops->restoresubvector) PetscUseTypeMethod(X, restoresubvector, is, Y);
1666:   else {
1667:     PetscCall(PetscObjectComposedDataGetInt((PetscObject)*Y, VecGetSubVectorSavedStateId, dummystate, unchanged));
1668:     if (!unchanged) { /* If Y's state has not changed since VecGetSubVector(), we only need to destroy Y */
1669:       VecScatter scatter;
1670:       PetscInt   state;

1672:       PetscCall(VecLockGet(X, &state));
1673:       PetscCheck(state == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vec X is locked for read-only or read/write access");

1675:       PetscCall(PetscObjectQuery((PetscObject)*Y, "VecGetSubVector_Scatter", (PetscObject *)&scatter));
1676:       if (scatter) {
1677:         PetscCall(VecScatterBegin(scatter, *Y, X, INSERT_VALUES, SCATTER_REVERSE));
1678:         PetscCall(VecScatterEnd(scatter, *Y, X, INSERT_VALUES, SCATTER_REVERSE));
1679:       } else {
1680:         PetscBool iscuda, iship;
1681:         PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &iscuda, VECSEQCUDA, VECMPICUDA, ""));
1682:         PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &iship, VECSEQHIP, VECMPIHIP, ""));

1684:         if (iscuda) {
1685: #if defined(PETSC_HAVE_CUDA)
1686:           PetscOffloadMask ymask = (*Y)->offloadmask;

1688:           /* The offloadmask of X dictates where to move memory
1689:               If X GPU data is valid, then move Y data on GPU if needed
1690:               Otherwise, move back to the CPU */
1691:           switch (X->offloadmask) {
1692:           case PETSC_OFFLOAD_BOTH:
1693:             if (ymask == PETSC_OFFLOAD_CPU) {
1694:               PetscCall(VecCUDAResetArray(*Y));
1695:             } else if (ymask == PETSC_OFFLOAD_GPU) {
1696:               X->offloadmask = PETSC_OFFLOAD_GPU;
1697:             }
1698:             break;
1699:           case PETSC_OFFLOAD_GPU:
1700:             if (ymask == PETSC_OFFLOAD_CPU) PetscCall(VecCUDAResetArray(*Y));
1701:             break;
1702:           case PETSC_OFFLOAD_CPU:
1703:             if (ymask == PETSC_OFFLOAD_GPU) PetscCall(VecResetArray(*Y));
1704:             break;
1705:           case PETSC_OFFLOAD_UNALLOCATED:
1706:           case PETSC_OFFLOAD_KOKKOS:
1707:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen");
1708:           }
1709: #endif
1710:         } else if (iship) {
1711: #if defined(PETSC_HAVE_HIP)
1712:           PetscOffloadMask ymask = (*Y)->offloadmask;

1714:           /* The offloadmask of X dictates where to move memory
1715:               If X GPU data is valid, then move Y data on GPU if needed
1716:               Otherwise, move back to the CPU */
1717:           switch (X->offloadmask) {
1718:           case PETSC_OFFLOAD_BOTH:
1719:             if (ymask == PETSC_OFFLOAD_CPU) {
1720:               PetscCall(VecHIPResetArray(*Y));
1721:             } else if (ymask == PETSC_OFFLOAD_GPU) {
1722:               X->offloadmask = PETSC_OFFLOAD_GPU;
1723:             }
1724:             break;
1725:           case PETSC_OFFLOAD_GPU:
1726:             if (ymask == PETSC_OFFLOAD_CPU) PetscCall(VecHIPResetArray(*Y));
1727:             break;
1728:           case PETSC_OFFLOAD_CPU:
1729:             if (ymask == PETSC_OFFLOAD_GPU) PetscCall(VecResetArray(*Y));
1730:             break;
1731:           case PETSC_OFFLOAD_UNALLOCATED:
1732:           case PETSC_OFFLOAD_KOKKOS:
1733:             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "This should not happen");
1734:           }
1735: #endif
1736:         } else {
1737:           /* If OpenCL vecs updated the device memory, this triggers a copy on the CPU */
1738:           PetscCall(VecResetArray(*Y));
1739:         }
1740:         PetscCall(PetscObjectStateIncrease((PetscObject)X));
1741:       }
1742:     }
1743:   }
1744:   PetscCall(VecDestroy(Y));
1745:   PetscFunctionReturn(PETSC_SUCCESS);
1746: }

1748: /*@
1749:   VecCreateLocalVector - Creates a vector object suitable for use with `VecGetLocalVector()` and friends. You must call `VecDestroy()` when the
1750:   vector is no longer needed.

1752:   Not Collective.

1754:   Input Parameter:
1755: . v - The vector for which the local vector is desired.

1757:   Output Parameter:
1758: . w - Upon exit this contains the local vector.

1760:   Level: beginner

1762: .seealso: [](ch_vectors), `Vec`, `VecGetLocalVectorRead()`, `VecRestoreLocalVectorRead()`, `VecGetLocalVector()`, `VecRestoreLocalVector()`
1763: @*/
1764: PetscErrorCode VecCreateLocalVector(Vec v, Vec *w)
1765: {
1766:   PetscMPIInt size;

1768:   PetscFunctionBegin;
1770:   PetscAssertPointer(w, 2);
1771:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
1772:   if (size == 1) PetscCall(VecDuplicate(v, w));
1773:   else if (v->ops->createlocalvector) PetscUseTypeMethod(v, createlocalvector, w);
1774:   else {
1775:     VecType  type;
1776:     PetscInt n;

1778:     PetscCall(VecCreate(PETSC_COMM_SELF, w));
1779:     PetscCall(VecGetLocalSize(v, &n));
1780:     PetscCall(VecSetSizes(*w, n, n));
1781:     PetscCall(VecGetBlockSize(v, &n));
1782:     PetscCall(VecSetBlockSize(*w, n));
1783:     PetscCall(VecGetType(v, &type));
1784:     PetscCall(VecSetType(*w, type));
1785:   }
1786:   PetscFunctionReturn(PETSC_SUCCESS);
1787: }

1789: /*@
1790:   VecGetLocalVectorRead - Maps the local portion of a vector into a
1791:   vector.

1793:   Not Collective.

1795:   Input Parameter:
1796: . v - The vector for which the local vector is desired.

1798:   Output Parameter:
1799: . w - Upon exit this contains the local vector.

1801:   Level: beginner

1803:   Notes:
1804:   You must call `VecRestoreLocalVectorRead()` when the local
1805:   vector is no longer needed.

1807:   This function is similar to `VecGetArrayRead()` which maps the local
1808:   portion into a raw pointer.  `VecGetLocalVectorRead()` is usually
1809:   almost as efficient as `VecGetArrayRead()` but in certain circumstances
1810:   `VecGetLocalVectorRead()` can be much more efficient than
1811:   `VecGetArrayRead()`.  This is because the construction of a contiguous
1812:   array representing the vector data required by `VecGetArrayRead()` can
1813:   be an expensive operation for certain vector types.  For example, for
1814:   GPU vectors `VecGetArrayRead()` requires that the data between device
1815:   and host is synchronized.

1817:   Unlike `VecGetLocalVector()`, this routine is not collective and
1818:   preserves cached information.

1820: .seealso: [](ch_vectors), `Vec`, `VecCreateLocalVector()`, `VecRestoreLocalVectorRead()`, `VecGetLocalVector()`, `VecGetArrayRead()`, `VecGetArray()`
1821: @*/
1822: PetscErrorCode VecGetLocalVectorRead(Vec v, Vec w)
1823: {
1824:   PetscFunctionBegin;
1827:   VecCheckSameLocalSize(v, 1, w, 2);
1828:   if (v->ops->getlocalvectorread) {
1829:     PetscUseTypeMethod(v, getlocalvectorread, w);
1830:   } else {
1831:     PetscScalar *a;

1833:     PetscCall(VecGetArrayRead(v, (const PetscScalar **)&a));
1834:     PetscCall(VecPlaceArray(w, a));
1835:   }
1836:   PetscCall(PetscObjectStateIncrease((PetscObject)w));
1837:   PetscCall(VecLockReadPush(v));
1838:   PetscCall(VecLockReadPush(w));
1839:   PetscFunctionReturn(PETSC_SUCCESS);
1840: }

1842: /*@
1843:   VecRestoreLocalVectorRead - Unmaps the local portion of a vector
1844:   previously mapped into a vector using `VecGetLocalVectorRead()`.

1846:   Not Collective.

1848:   Input Parameters:
1849: + v - The local portion of this vector was previously mapped into `w` using `VecGetLocalVectorRead()`.
1850: - w - The vector into which the local portion of `v` was mapped.

1852:   Level: beginner

1854: .seealso: [](ch_vectors), `Vec`, `VecCreateLocalVector()`, `VecGetLocalVectorRead()`, `VecGetLocalVector()`, `VecGetArrayRead()`, `VecGetArray()`
1855: @*/
1856: PetscErrorCode VecRestoreLocalVectorRead(Vec v, Vec w)
1857: {
1858:   PetscFunctionBegin;
1861:   if (v->ops->restorelocalvectorread) {
1862:     PetscUseTypeMethod(v, restorelocalvectorread, w);
1863:   } else {
1864:     const PetscScalar *a;

1866:     PetscCall(VecGetArrayRead(w, &a));
1867:     PetscCall(VecRestoreArrayRead(v, &a));
1868:     PetscCall(VecResetArray(w));
1869:   }
1870:   PetscCall(VecLockReadPop(v));
1871:   PetscCall(VecLockReadPop(w));
1872:   PetscCall(PetscObjectStateIncrease((PetscObject)w));
1873:   PetscFunctionReturn(PETSC_SUCCESS);
1874: }

1876: /*@
1877:   VecGetLocalVector - Maps the local portion of a vector into a
1878:   vector.

1880:   Collective

1882:   Input Parameter:
1883: . v - The vector for which the local vector is desired.

1885:   Output Parameter:
1886: . w - Upon exit this contains the local vector.

1888:   Level: beginner

1890:   Notes:
1891:   You must call `VecRestoreLocalVector()` when the local
1892:   vector is no longer needed.

1894:   This function is similar to `VecGetArray()` which maps the local
1895:   portion into a raw pointer.  `VecGetLocalVector()` is usually about as
1896:   efficient as `VecGetArray()` but in certain circumstances
1897:   `VecGetLocalVector()` can be much more efficient than `VecGetArray()`.
1898:   This is because the construction of a contiguous array representing
1899:   the vector data required by `VecGetArray()` can be an expensive
1900:   operation for certain vector types.  For example, for GPU vectors
1901:   `VecGetArray()` requires that the data between device and host is
1902:   synchronized.

1904: .seealso: [](ch_vectors), `Vec`, `VecCreateLocalVector()`, `VecRestoreLocalVector()`, `VecGetLocalVectorRead()`, `VecGetArrayRead()`, `VecGetArray()`
1905: @*/
1906: PetscErrorCode VecGetLocalVector(Vec v, Vec w)
1907: {
1908:   PetscFunctionBegin;
1911:   VecCheckSameLocalSize(v, 1, w, 2);
1912:   if (v->ops->getlocalvector) {
1913:     PetscUseTypeMethod(v, getlocalvector, w);
1914:   } else {
1915:     PetscScalar *a;

1917:     PetscCall(VecGetArray(v, &a));
1918:     PetscCall(VecPlaceArray(w, a));
1919:   }
1920:   PetscCall(PetscObjectStateIncrease((PetscObject)w));
1921:   PetscFunctionReturn(PETSC_SUCCESS);
1922: }

1924: /*@
1925:   VecRestoreLocalVector - Unmaps the local portion of a vector
1926:   previously mapped into a vector using `VecGetLocalVector()`.

1928:   Logically Collective.

1930:   Input Parameters:
1931: + v - The local portion of this vector was previously mapped into `w` using `VecGetLocalVector()`.
1932: - w - The vector into which the local portion of `v` was mapped.

1934:   Level: beginner

1936: .seealso: [](ch_vectors), `Vec`, `VecCreateLocalVector()`, `VecGetLocalVector()`, `VecGetLocalVectorRead()`, `VecRestoreLocalVectorRead()`, `LocalVectorRead()`, `VecGetArrayRead()`, `VecGetArray()`
1937: @*/
1938: PetscErrorCode VecRestoreLocalVector(Vec v, Vec w)
1939: {
1940:   PetscFunctionBegin;
1943:   if (v->ops->restorelocalvector) {
1944:     PetscUseTypeMethod(v, restorelocalvector, w);
1945:   } else {
1946:     PetscScalar *a;
1947:     PetscCall(VecGetArray(w, &a));
1948:     PetscCall(VecRestoreArray(v, &a));
1949:     PetscCall(VecResetArray(w));
1950:   }
1951:   PetscCall(PetscObjectStateIncrease((PetscObject)w));
1952:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
1953:   PetscFunctionReturn(PETSC_SUCCESS);
1954: }

1956: /*@C
1957:   VecGetArray - Returns a pointer to a contiguous array that contains this
1958:   MPI processes's portion of the vector data

1960:   Logically Collective

1962:   Input Parameter:
1963: . x - the vector

1965:   Output Parameter:
1966: . a - location to put pointer to the array

1968:   Level: beginner

1970:   Notes:
1971:   For the standard PETSc vectors, `VecGetArray()` returns a pointer to the local data array and
1972:   does not use any copies. If the underlying vector data is not stored in a contiguous array
1973:   this routine will copy the data to a contiguous array and return a pointer to that. You MUST
1974:   call `VecRestoreArray()` when you no longer need access to the array.

1976:   Fortran Notes:
1977:   `VecGetArray()` Fortran binding is deprecated (since PETSc 3.19), use `VecGetArrayF90()`

1979: .seealso: [](ch_vectors), `Vec`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecGetArrayReadF90()`, `VecPlaceArray()`, `VecGetArray2d()`,
1980:           `VecGetArrayPair()`, `VecRestoreArrayPair()`, `VecGetArrayWrite()`, `VecRestoreArrayWrite()`
1981: @*/
1982: PetscErrorCode VecGetArray(Vec x, PetscScalar **a)
1983: {
1984:   PetscFunctionBegin;
1986:   PetscCall(VecSetErrorIfLocked(x, 1));
1987:   if (x->ops->getarray) { /* The if-else order matters! VECNEST, VECCUDA etc should have ops->getarray while VECCUDA etc are petscnative */
1988:     PetscUseTypeMethod(x, getarray, a);
1989:   } else if (x->petscnative) { /* VECSTANDARD */
1990:     *a = *((PetscScalar **)x->data);
1991:   } else SETERRQ(PetscObjectComm((PetscObject)x), PETSC_ERR_SUP, "Cannot get array for vector type \"%s\"", ((PetscObject)x)->type_name);
1992:   PetscFunctionReturn(PETSC_SUCCESS);
1993: }

1995: /*@C
1996:   VecRestoreArray - Restores a vector after `VecGetArray()` has been called and the array is no longer needed

1998:   Logically Collective

2000:   Input Parameters:
2001: + x - the vector
2002: - a - location of pointer to array obtained from `VecGetArray()`

2004:   Level: beginner

2006:   Fortran Notes:
2007:   `VecRestoreArray()` Fortran binding is deprecated (since PETSc 3.19), use `VecRestoreArrayF90()`

2009: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArrayRead()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecRestoreArrayReadF90()`, `VecPlaceArray()`, `VecRestoreArray2d()`,
2010:           `VecGetArrayPair()`, `VecRestoreArrayPair()`
2011: @*/
2012: PetscErrorCode VecRestoreArray(Vec x, PetscScalar **a)
2013: {
2014:   PetscFunctionBegin;
2016:   if (a) PetscAssertPointer(a, 2);
2017:   if (x->ops->restorearray) {
2018:     PetscUseTypeMethod(x, restorearray, a);
2019:   } else PetscCheck(x->petscnative, PetscObjectComm((PetscObject)x), PETSC_ERR_SUP, "Cannot restore array for vector type \"%s\"", ((PetscObject)x)->type_name);
2020:   if (a) *a = NULL;
2021:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
2022:   PetscFunctionReturn(PETSC_SUCCESS);
2023: }
2024: /*@C
2025:   VecGetArrayRead - Get read-only pointer to contiguous array containing this processor's portion of the vector data.

2027:   Not Collective

2029:   Input Parameter:
2030: . x - the vector

2032:   Output Parameter:
2033: . a - the array

2035:   Level: beginner

2037:   Notes:
2038:   The array must be returned using a matching call to `VecRestoreArrayRead()`.

2040:   Unlike `VecGetArray()`, preserves cached information like vector norms.

2042:   Standard PETSc vectors use contiguous storage so that this routine does not perform a copy.  Other vector
2043:   implementations may require a copy, but such implementations should cache the contiguous representation so that
2044:   only one copy is performed when this routine is called multiple times in sequence.

2046:   Fortran Notes:
2047:   `VecGetArrayRead()` Fortran binding is deprecated (since PETSc 3.19), use `VecGetArrayReadF90()`

2049: .seealso: [](ch_vectors), `Vec`, `VecGetArrayReadF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2050: @*/
2051: PetscErrorCode VecGetArrayRead(Vec x, const PetscScalar **a)
2052: {
2053:   PetscFunctionBegin;
2055:   PetscAssertPointer(a, 2);
2056:   if (x->ops->getarrayread) {
2057:     PetscUseTypeMethod(x, getarrayread, a);
2058:   } else if (x->ops->getarray) {
2059:     PetscObjectState state;

2061:     /* VECNEST, VECCUDA, VECKOKKOS etc */
2062:     // x->ops->getarray may bump the object state, but since we know this is a read-only get
2063:     // we can just undo that
2064:     PetscCall(PetscObjectStateGet((PetscObject)x, &state));
2065:     PetscUseTypeMethod(x, getarray, (PetscScalar **)a);
2066:     PetscCall(PetscObjectStateSet((PetscObject)x, state));
2067:   } else if (x->petscnative) {
2068:     /* VECSTANDARD */
2069:     *a = *((PetscScalar **)x->data);
2070:   } else SETERRQ(PetscObjectComm((PetscObject)x), PETSC_ERR_SUP, "Cannot get array read for vector type \"%s\"", ((PetscObject)x)->type_name);
2071:   PetscFunctionReturn(PETSC_SUCCESS);
2072: }

2074: /*@C
2075:   VecRestoreArrayRead - Restore array obtained with `VecGetArrayRead()`

2077:   Not Collective

2079:   Input Parameters:
2080: + x - the vector
2081: - a - the array

2083:   Level: beginner

2085:   Fortran Notes:
2086:   `VecRestoreArrayRead()` Fortran binding is deprecated (since PETSc 3.19), use `VecRestoreArrayReadF90()`

2088: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayReadF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2089: @*/
2090: PetscErrorCode VecRestoreArrayRead(Vec x, const PetscScalar **a)
2091: {
2092:   PetscFunctionBegin;
2094:   if (a) PetscAssertPointer(a, 2);
2095:   if (x->petscnative) { /* VECSTANDARD, VECCUDA, VECKOKKOS etc */
2096:     /* nothing */
2097:   } else if (x->ops->restorearrayread) { /* VECNEST */
2098:     PetscUseTypeMethod(x, restorearrayread, a);
2099:   } else { /* No one? */
2100:     PetscObjectState state;

2102:     // x->ops->restorearray may bump the object state, but since we know this is a read-restore
2103:     // we can just undo that
2104:     PetscCall(PetscObjectStateGet((PetscObject)x, &state));
2105:     PetscUseTypeMethod(x, restorearray, (PetscScalar **)a);
2106:     PetscCall(PetscObjectStateSet((PetscObject)x, state));
2107:   }
2108:   if (a) *a = NULL;
2109:   PetscFunctionReturn(PETSC_SUCCESS);
2110: }

2112: /*@C
2113:   VecGetArrayWrite - Returns a pointer to a contiguous array that WILL contain this
2114:   MPI processes's portion of the vector data.

2116:   Logically Collective

2118:   Input Parameter:
2119: . x - the vector

2121:   Output Parameter:
2122: . a - location to put pointer to the array

2124:   Level: intermediate

2126:   Note:
2127:   The values in this array are NOT valid, the caller of this routine is responsible for putting
2128:   values into the array; any values it does not set will be invalid.

2130:   The array must be returned using a matching call to `VecRestoreArrayRead()`.

2132:   For vectors associated with GPUs, the host and device vectors are not synchronized before
2133:   giving access. If you need correct values in the array use `VecGetArray()`

2135:   Fortran Notes:
2136:   `VecGetArrayWrite()` Fortran binding is deprecated (since PETSc 3.19), use `VecGetArrayWriteF90()`

2138: .seealso: [](ch_vectors), `Vec`, `VecGetArrayWriteF90()`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecGetArrayReadF90()`, `VecPlaceArray()`, `VecGetArray2d()`,
2139:           `VecGetArrayPair()`, `VecRestoreArrayPair()`, `VecGetArray()`, `VecRestoreArrayWrite()`
2140: @*/
2141: PetscErrorCode VecGetArrayWrite(Vec x, PetscScalar **a)
2142: {
2143:   PetscFunctionBegin;
2145:   PetscAssertPointer(a, 2);
2146:   PetscCall(VecSetErrorIfLocked(x, 1));
2147:   if (x->ops->getarraywrite) {
2148:     PetscUseTypeMethod(x, getarraywrite, a);
2149:   } else {
2150:     PetscCall(VecGetArray(x, a));
2151:   }
2152:   PetscFunctionReturn(PETSC_SUCCESS);
2153: }

2155: /*@C
2156:   VecRestoreArrayWrite - Restores a vector after `VecGetArrayWrite()` has been called.

2158:   Logically Collective

2160:   Input Parameters:
2161: + x - the vector
2162: - a - location of pointer to array obtained from `VecGetArray()`

2164:   Level: beginner

2166:   Fortran Notes:
2167:   `VecRestoreArrayWrite()` Fortran binding is deprecated (since PETSc 3.19), use `VecRestoreArrayWriteF90()`

2169: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArrayRead()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecRestoreArrayReadF90()`, `VecPlaceArray()`, `VecRestoreArray2d()`,
2170:           `VecGetArrayPair()`, `VecRestoreArrayPair()`, `VecGetArrayWrite()`
2171: @*/
2172: PetscErrorCode VecRestoreArrayWrite(Vec x, PetscScalar **a)
2173: {
2174:   PetscFunctionBegin;
2176:   if (a) PetscAssertPointer(a, 2);
2177:   if (x->ops->restorearraywrite) {
2178:     PetscUseTypeMethod(x, restorearraywrite, a);
2179:   } else if (x->ops->restorearray) {
2180:     PetscUseTypeMethod(x, restorearray, a);
2181:   }
2182:   if (a) *a = NULL;
2183:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
2184:   PetscFunctionReturn(PETSC_SUCCESS);
2185: }

2187: /*@C
2188:   VecGetArrays - Returns a pointer to the arrays in a set of vectors
2189:   that were created by a call to `VecDuplicateVecs()`.

2191:   Logically Collective; No Fortran Support

2193:   Input Parameters:
2194: + x - the vectors
2195: - n - the number of vectors

2197:   Output Parameter:
2198: . a - location to put pointer to the array

2200:   Level: intermediate

2202:   Note:
2203:   You MUST call `VecRestoreArrays()` when you no longer need access to the arrays.

2205: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArrays()`
2206: @*/
2207: PetscErrorCode VecGetArrays(const Vec x[], PetscInt n, PetscScalar **a[])
2208: {
2209:   PetscInt      i;
2210:   PetscScalar **q;

2212:   PetscFunctionBegin;
2213:   PetscAssertPointer(x, 1);
2215:   PetscAssertPointer(a, 3);
2216:   PetscCheck(n > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must get at least one array n = %" PetscInt_FMT, n);
2217:   PetscCall(PetscMalloc1(n, &q));
2218:   for (i = 0; i < n; ++i) PetscCall(VecGetArray(x[i], &q[i]));
2219:   *a = q;
2220:   PetscFunctionReturn(PETSC_SUCCESS);
2221: }

2223: /*@C
2224:   VecRestoreArrays - Restores a group of vectors after `VecGetArrays()`
2225:   has been called.

2227:   Logically Collective; No Fortran Support

2229:   Input Parameters:
2230: + x - the vector
2231: . n - the number of vectors
2232: - a - location of pointer to arrays obtained from `VecGetArrays()`

2234:   Notes:
2235:   For regular PETSc vectors this routine does not involve any copies. For
2236:   any special vectors that do not store local vector data in a contiguous
2237:   array, this routine will copy the data back into the underlying
2238:   vector data structure from the arrays obtained with `VecGetArrays()`.

2240:   Level: intermediate

2242: .seealso: [](ch_vectors), `Vec`, `VecGetArrays()`, `VecRestoreArray()`
2243: @*/
2244: PetscErrorCode VecRestoreArrays(const Vec x[], PetscInt n, PetscScalar **a[])
2245: {
2246:   PetscInt      i;
2247:   PetscScalar **q = *a;

2249:   PetscFunctionBegin;
2250:   PetscAssertPointer(x, 1);
2252:   PetscAssertPointer(a, 3);

2254:   for (i = 0; i < n; ++i) PetscCall(VecRestoreArray(x[i], &q[i]));
2255:   PetscCall(PetscFree(q));
2256:   PetscFunctionReturn(PETSC_SUCCESS);
2257: }

2259: /*@C
2260:   VecGetArrayAndMemType - Like `VecGetArray()`, but if this is a standard device vector (e.g.,
2261:   `VECCUDA`), the returned pointer will be a device pointer to the device memory that contains
2262:   this MPI processes's portion of the vector data.

2264:   Logically Collective; No Fortran Support

2266:   Input Parameter:
2267: . x - the vector

2269:   Output Parameters:
2270: + a     - location to put pointer to the array
2271: - mtype - memory type of the array

2273:   Level: beginner

2275:   Note:
2276:   Device data is guaranteed to have the latest value. Otherwise, when this is a host vector
2277:   (e.g., `VECMPI`), this routine functions the same as `VecGetArray()` and returns a host
2278:   pointer.

2280:   For `VECKOKKOS`, if Kokkos is configured without device (e.g., use serial or openmp), per
2281:   this function, the vector works like `VECSEQ`/`VECMPI`; otherwise, it works like `VECCUDA` or
2282:   `VECHIP` etc.

2284:   Use `VecRestoreArrayAndMemType()` when the array access is no longer needed.

2286: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayAndMemType()`, `VecGetArrayReadAndMemType()`, `VecGetArrayWriteAndMemType()`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecGetArrayReadF90()`,
2287:           `VecPlaceArray()`, `VecGetArray2d()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`, `VecGetArrayWrite()`, `VecRestoreArrayWrite()`
2288: @*/
2289: PetscErrorCode VecGetArrayAndMemType(Vec x, PetscScalar **a, PetscMemType *mtype)
2290: {
2291:   PetscFunctionBegin;
2294:   PetscAssertPointer(a, 2);
2295:   if (mtype) PetscAssertPointer(mtype, 3);
2296:   PetscCall(VecSetErrorIfLocked(x, 1));
2297:   if (x->ops->getarrayandmemtype) {
2298:     /* VECCUDA, VECKOKKOS etc */
2299:     PetscUseTypeMethod(x, getarrayandmemtype, a, mtype);
2300:   } else {
2301:     /* VECSTANDARD, VECNEST, VECVIENNACL */
2302:     PetscCall(VecGetArray(x, a));
2303:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
2304:   }
2305:   PetscFunctionReturn(PETSC_SUCCESS);
2306: }

2308: /*@C
2309:   VecRestoreArrayAndMemType - Restores a vector after `VecGetArrayAndMemType()` has been called.

2311:   Logically Collective; No Fortran Support

2313:   Input Parameters:
2314: + x - the vector
2315: - a - location of pointer to array obtained from `VecGetArrayAndMemType()`

2317:   Level: beginner

2319: .seealso: [](ch_vectors), `Vec`, `VecGetArrayAndMemType()`, `VecGetArray()`, `VecRestoreArrayRead()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecRestoreArrayReadF90()`,
2320:           `VecPlaceArray()`, `VecRestoreArray2d()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2321: @*/
2322: PetscErrorCode VecRestoreArrayAndMemType(Vec x, PetscScalar **a)
2323: {
2324:   PetscFunctionBegin;
2327:   if (a) PetscAssertPointer(a, 2);
2328:   if (x->ops->restorearrayandmemtype) {
2329:     /* VECCUDA, VECKOKKOS etc */
2330:     PetscUseTypeMethod(x, restorearrayandmemtype, a);
2331:   } else {
2332:     /* VECNEST, VECVIENNACL */
2333:     PetscCall(VecRestoreArray(x, a));
2334:   } /* VECSTANDARD does nothing */
2335:   if (a) *a = NULL;
2336:   PetscCall(PetscObjectStateIncrease((PetscObject)x));
2337:   PetscFunctionReturn(PETSC_SUCCESS);
2338: }

2340: /*@C
2341:   VecGetArrayReadAndMemType - Like `VecGetArrayRead()`, but if the input vector is a device vector, it will return a read-only device pointer.
2342:   The returned pointer is guaranteed to point to up-to-date data. For host vectors, it functions as `VecGetArrayRead()`.

2344:   Not Collective; No Fortran Support

2346:   Input Parameter:
2347: . x - the vector

2349:   Output Parameters:
2350: + a     - the array
2351: - mtype - memory type of the array

2353:   Level: beginner

2355:   Notes:
2356:   The array must be returned using a matching call to `VecRestoreArrayReadAndMemType()`.

2358: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayReadAndMemType()`, `VecGetArrayAndMemType()`, `VecGetArrayWriteAndMemType()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2359: @*/
2360: PetscErrorCode VecGetArrayReadAndMemType(Vec x, const PetscScalar **a, PetscMemType *mtype)
2361: {
2362:   PetscFunctionBegin;
2365:   PetscAssertPointer(a, 2);
2366:   if (mtype) PetscAssertPointer(mtype, 3);
2367:   if (x->ops->getarrayreadandmemtype) {
2368:     /* VECCUDA/VECHIP though they are also petscnative */
2369:     PetscUseTypeMethod(x, getarrayreadandmemtype, a, mtype);
2370:   } else if (x->ops->getarrayandmemtype) {
2371:     /* VECKOKKOS */
2372:     PetscObjectState state;

2374:     // see VecGetArrayRead() for why
2375:     PetscCall(PetscObjectStateGet((PetscObject)x, &state));
2376:     PetscUseTypeMethod(x, getarrayandmemtype, (PetscScalar **)a, mtype);
2377:     PetscCall(PetscObjectStateSet((PetscObject)x, state));
2378:   } else {
2379:     PetscCall(VecGetArrayRead(x, a));
2380:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
2381:   }
2382:   PetscFunctionReturn(PETSC_SUCCESS);
2383: }

2385: /*@C
2386:   VecRestoreArrayReadAndMemType - Restore array obtained with `VecGetArrayReadAndMemType()`

2388:   Not Collective; No Fortran Support

2390:   Input Parameters:
2391: + x - the vector
2392: - a - the array

2394:   Level: beginner

2396: .seealso: [](ch_vectors), `Vec`, `VecGetArrayReadAndMemType()`, `VecRestoreArrayAndMemType()`, `VecRestoreArrayWriteAndMemType()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2397: @*/
2398: PetscErrorCode VecRestoreArrayReadAndMemType(Vec x, const PetscScalar **a)
2399: {
2400:   PetscFunctionBegin;
2403:   if (a) PetscAssertPointer(a, 2);
2404:   if (x->ops->restorearrayreadandmemtype) {
2405:     /* VECCUDA/VECHIP */
2406:     PetscUseTypeMethod(x, restorearrayreadandmemtype, a);
2407:   } else if (!x->petscnative) {
2408:     /* VECNEST */
2409:     PetscCall(VecRestoreArrayRead(x, a));
2410:   }
2411:   if (a) *a = NULL;
2412:   PetscFunctionReturn(PETSC_SUCCESS);
2413: }

2415: /*@C
2416:   VecGetArrayWriteAndMemType - Like `VecGetArrayWrite()`, but if this is a device vector it will always return
2417:   a device pointer to the device memory that contains this processor's portion of the vector data.

2419:   Logically Collective; No Fortran Support

2421:   Input Parameter:
2422: . x - the vector

2424:   Output Parameters:
2425: + a     - the array
2426: - mtype - memory type of the array

2428:   Level: beginner

2430:   Note:
2431:   The array must be returned using a matching call to `VecRestoreArrayWriteAndMemType()`, where it will label the device memory as most recent.

2433: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayWriteAndMemType()`, `VecGetArrayReadAndMemType()`, `VecGetArrayAndMemType()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`,
2434: @*/
2435: PetscErrorCode VecGetArrayWriteAndMemType(Vec x, PetscScalar **a, PetscMemType *mtype)
2436: {
2437:   PetscFunctionBegin;
2440:   PetscCall(VecSetErrorIfLocked(x, 1));
2441:   PetscAssertPointer(a, 2);
2442:   if (mtype) PetscAssertPointer(mtype, 3);
2443:   if (x->ops->getarraywriteandmemtype) {
2444:     /* VECCUDA, VECHIP, VECKOKKOS etc, though they are also petscnative */
2445:     PetscUseTypeMethod(x, getarraywriteandmemtype, a, mtype);
2446:   } else if (x->ops->getarrayandmemtype) {
2447:     PetscCall(VecGetArrayAndMemType(x, a, mtype));
2448:   } else {
2449:     /* VECNEST, VECVIENNACL */
2450:     PetscCall(VecGetArrayWrite(x, a));
2451:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
2452:   }
2453:   PetscFunctionReturn(PETSC_SUCCESS);
2454: }

2456: /*@C
2457:   VecRestoreArrayWriteAndMemType - Restore array obtained with `VecGetArrayWriteAndMemType()`

2459:   Logically Collective; No Fortran Support

2461:   Input Parameters:
2462: + x - the vector
2463: - a - the array

2465:   Level: beginner

2467: .seealso: [](ch_vectors), `Vec`, `VecGetArrayWriteAndMemType()`, `VecRestoreArrayAndMemType()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayPair()`, `VecRestoreArrayPair()`
2468: @*/
2469: PetscErrorCode VecRestoreArrayWriteAndMemType(Vec x, PetscScalar **a)
2470: {
2471:   PetscFunctionBegin;
2474:   PetscCall(VecSetErrorIfLocked(x, 1));
2475:   if (a) PetscAssertPointer(a, 2);
2476:   if (x->ops->restorearraywriteandmemtype) {
2477:     /* VECCUDA/VECHIP */
2478:     PetscMemType PETSC_UNUSED mtype; // since this function doesn't accept a memtype?
2479:     PetscUseTypeMethod(x, restorearraywriteandmemtype, a, &mtype);
2480:   } else if (x->ops->restorearrayandmemtype) {
2481:     PetscCall(VecRestoreArrayAndMemType(x, a));
2482:   } else {
2483:     PetscCall(VecRestoreArray(x, a));
2484:   }
2485:   if (a) *a = NULL;
2486:   PetscFunctionReturn(PETSC_SUCCESS);
2487: }

2489: /*@
2490:   VecPlaceArray - Allows one to replace the array in a vector with an
2491:   array provided by the user. This is useful to avoid copying an array
2492:   into a vector.

2494:   Logically Collective; No Fortran Support

2496:   Input Parameters:
2497: + vec   - the vector
2498: - array - the array

2500:   Level: developer

2502:   Notes:
2503:   Use `VecReplaceArray()` instead to permanently replace the array

2505:   You can return to the original array with a call to `VecResetArray()`. `vec` does not take
2506:   ownership of `array` in any way.

2508:   The user must free `array` themselves but be careful not to
2509:   do so before the vector has either been destroyed, had its original array restored with
2510:   `VecResetArray()` or permanently replaced with `VecReplaceArray()`.

2512: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecReplaceArray()`, `VecResetArray()`
2513: @*/
2514: PetscErrorCode VecPlaceArray(Vec vec, const PetscScalar array[])
2515: {
2516:   PetscFunctionBegin;
2519:   if (array) PetscAssertPointer(array, 2);
2520:   PetscUseTypeMethod(vec, placearray, array);
2521:   PetscCall(PetscObjectStateIncrease((PetscObject)vec));
2522:   PetscFunctionReturn(PETSC_SUCCESS);
2523: }

2525: /*@C
2526:   VecReplaceArray - Allows one to replace the array in a vector with an
2527:   array provided by the user. This is useful to avoid copying an array
2528:   into a vector.

2530:   Logically Collective; No Fortran Support

2532:   Input Parameters:
2533: + vec   - the vector
2534: - array - the array

2536:   Level: developer

2538:   Notes:
2539:   This permanently replaces the array and frees the memory associated
2540:   with the old array. Use `VecPlaceArray()` to temporarily replace the array.

2542:   The memory passed in MUST be obtained with `PetscMalloc()` and CANNOT be
2543:   freed by the user. It will be freed when the vector is destroyed.

2545: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecPlaceArray()`, `VecResetArray()`
2546: @*/
2547: PetscErrorCode VecReplaceArray(Vec vec, const PetscScalar array[])
2548: {
2549:   PetscFunctionBegin;
2552:   PetscUseTypeMethod(vec, replacearray, array);
2553:   PetscCall(PetscObjectStateIncrease((PetscObject)vec));
2554:   PetscFunctionReturn(PETSC_SUCCESS);
2555: }

2557: /*MC
2558:     VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
2559:     and makes them accessible via a Fortran pointer.

2561:     Synopsis:
2562:     VecDuplicateVecsF90(Vec x,PetscInt n,{Vec, pointer :: y(:)},integer ierr)

2564:     Collective

2566:     Input Parameters:
2567: +   x - a vector to mimic
2568: -   n - the number of vectors to obtain

2570:     Output Parameters:
2571: +   y - Fortran pointer to the array of vectors
2572: -   ierr - error code

2574:     Example of Usage:
2575: .vb
2576: #include <petsc/finclude/petscvec.h>
2577:     use petscvec

2579:     Vec x
2580:     Vec, pointer :: y(:)
2581:     ....
2582:     call VecDuplicateVecsF90(x,2,y,ierr)
2583:     call VecSet(y(2),alpha,ierr)
2584:     call VecSet(y(2),alpha,ierr)
2585:     ....
2586:     call VecDestroyVecsF90(2,y,ierr)
2587: .ve

2589:     Level: beginner

2591:     Note:
2592:     Use `VecDestroyVecsF90()` to free the space.

2594: .seealso: [](ch_vectors), `Vec`, `VecDestroyVecsF90()`, `VecDuplicateVecs()`
2595: M*/

2597: /*MC
2598:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
2599:     `VecGetArrayF90()`.

2601:     Synopsis:
2602:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2604:     Logically Collective

2606:     Input Parameters:
2607: +   x - vector
2608: -   xx_v - the Fortran pointer to the array

2610:     Output Parameter:
2611: .   ierr - error code

2613:     Example of Usage:
2614: .vb
2615: #include <petsc/finclude/petscvec.h>
2616:     use petscvec

2618:     PetscScalar, pointer :: xx_v(:)
2619:     ....
2620:     call VecGetArrayF90(x,xx_v,ierr)
2621:     xx_v(3) = a
2622:     call VecRestoreArrayF90(x,xx_v,ierr)
2623: .ve

2625:     Level: beginner

2627: .seealso: [](ch_vectors), `Vec`, `VecGetArrayF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrayReadF90()`
2628: M*/

2630: /*MC
2631:     VecDestroyVecsF90 - Frees a block of vectors obtained with `VecDuplicateVecsF90()`.

2633:     Synopsis:
2634:     VecDestroyVecsF90(PetscInt n,{Vec, pointer :: x(:)},PetscErrorCode ierr)

2636:     Collective

2638:     Input Parameters:
2639: +   n - the number of vectors previously obtained
2640: -   x - pointer to array of vector pointers

2642:     Output Parameter:
2643: .   ierr - error code

2645:     Level: beginner

2647: .seealso: [](ch_vectors), `Vec`, `VecDestroyVecs()`, `VecDuplicateVecsF90()`
2648: M*/

2650: /*MC
2651:     VecGetArrayF90 - Accesses a vector array from Fortran. For default PETSc
2652:     vectors, `VecGetArrayF90()` returns a pointer to the local data array. Otherwise,
2653:     this routine is implementation dependent. You MUST call `VecRestoreArrayF90()`
2654:     when you no longer need access to the array.

2656:     Synopsis:
2657:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2659:     Logically Collective

2661:     Input Parameter:
2662: .   x - vector

2664:     Output Parameters:
2665: +   xx_v - the Fortran pointer to the array
2666: -   ierr - error code

2668:     Example of Usage:
2669: .vb
2670: #include <petsc/finclude/petscvec.h>
2671:     use petscvec

2673:     PetscScalar, pointer :: xx_v(:)
2674:     ....
2675:     call VecGetArrayF90(x,xx_v,ierr)
2676:     xx_v(3) = a
2677:     call VecRestoreArrayF90(x,xx_v,ierr)
2678: .ve

2680:      Level: beginner

2682:     Note:
2683:     If you ONLY intend to read entries from the array and not change any entries you should use `VecGetArrayReadF90()`.

2685: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayReadF90()`
2686: M*/

2688: /*MC
2689:     VecGetArrayReadF90 - Accesses a read only array from Fortran. For default PETSc
2690:     vectors, `VecGetArrayF90()` returns a pointer to the local data array. Otherwise,
2691:     this routine is implementation dependent. You MUST call `VecRestoreArrayReadF90()`
2692:     when you no longer need access to the array.

2694:     Synopsis:
2695:     VecGetArrayReadF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2697:     Logically Collective

2699:     Input Parameter:
2700: .   x - vector

2702:     Output Parameters:
2703: +   xx_v - the Fortran pointer to the array
2704: -   ierr - error code

2706:     Example of Usage:
2707: .vb
2708: #include <petsc/finclude/petscvec.h>
2709:     use petscvec

2711:     PetscScalar, pointer :: xx_v(:)
2712:     ....
2713:     call VecGetArrayReadF90(x,xx_v,ierr)
2714:     a = xx_v(3)
2715:     call VecRestoreArrayReadF90(x,xx_v,ierr)
2716: .ve

2718:     Level: beginner

2720:     Note:
2721:     If you intend to write entries into the array you must use `VecGetArrayF90()`.

2723: .seealso: [](ch_vectors), `Vec`, `VecRestoreArrayReadF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecRestoreArrayRead()`, `VecGetArrayF90()`
2724: M*/

2726: /*MC
2727:     VecRestoreArrayReadF90 - Restores a readonly vector to a usable state after a call to
2728:     `VecGetArrayReadF90()`.

2730:     Synopsis:
2731:     VecRestoreArrayReadF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2733:     Logically Collective

2735:     Input Parameters:
2736: +   x - vector
2737: -   xx_v - the Fortran pointer to the array

2739:     Output Parameter:
2740: .   ierr - error code

2742:     Example of Usage:
2743: .vb
2744: #include <petsc/finclude/petscvec.h>
2745:     use petscvec

2747:     PetscScalar, pointer :: xx_v(:)
2748:     ....
2749:     call VecGetArrayReadF90(x,xx_v,ierr)
2750:     a = xx_v(3)
2751:     call VecRestoreArrayReadF90(x,xx_v,ierr)
2752: .ve

2754:     Level: beginner

2756: .seealso: [](ch_vectors), `Vec`, `VecGetArrayReadF90()`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecRestoreArrayRead()`, `VecRestoreArrayF90()`
2757: M*/

2759: /*@C
2760:   VecGetArray2d - Returns a pointer to a 2d contiguous array that contains this
2761:   processor's portion of the vector data.  You MUST call `VecRestoreArray2d()`
2762:   when you no longer need access to the array.

2764:   Logically Collective

2766:   Input Parameters:
2767: + x      - the vector
2768: . m      - first dimension of two dimensional array
2769: . n      - second dimension of two dimensional array
2770: . mstart - first index you will use in first coordinate direction (often 0)
2771: - nstart - first index in the second coordinate direction (often 0)

2773:   Output Parameter:
2774: . a - location to put pointer to the array

2776:   Level: developer

2778:   Notes:
2779:   For a vector obtained from `DMCreateLocalVector()` `mstart` and `nstart` are likely
2780:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
2781:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
2782:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray2d()`.

2784:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2786: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
2787:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
2788:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
2789: @*/
2790: PetscErrorCode VecGetArray2d(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
2791: {
2792:   PetscInt     i, N;
2793:   PetscScalar *aa;

2795:   PetscFunctionBegin;
2797:   PetscAssertPointer(a, 6);
2799:   PetscCall(VecGetLocalSize(x, &N));
2800:   PetscCheck(m * n == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 2d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n);
2801:   PetscCall(VecGetArray(x, &aa));

2803:   PetscCall(PetscMalloc1(m, a));
2804:   for (i = 0; i < m; i++) (*a)[i] = aa + i * n - nstart;
2805:   *a -= mstart;
2806:   PetscFunctionReturn(PETSC_SUCCESS);
2807: }

2809: /*@C
2810:   VecGetArray2dWrite - Returns a pointer to a 2d contiguous array that will contain this
2811:   processor's portion of the vector data.  You MUST call `VecRestoreArray2dWrite()`
2812:   when you no longer need access to the array.

2814:   Logically Collective

2816:   Input Parameters:
2817: + x      - the vector
2818: . m      - first dimension of two dimensional array
2819: . n      - second dimension of two dimensional array
2820: . mstart - first index you will use in first coordinate direction (often 0)
2821: - nstart - first index in the second coordinate direction (often 0)

2823:   Output Parameter:
2824: . a - location to put pointer to the array

2826:   Level: developer

2828:   Notes:
2829:   For a vector obtained from `DMCreateLocalVector()` `mstart` and `nstart` are likely
2830:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
2831:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
2832:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray2d()`.

2834:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2836: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
2837:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
2838:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
2839: @*/
2840: PetscErrorCode VecGetArray2dWrite(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
2841: {
2842:   PetscInt     i, N;
2843:   PetscScalar *aa;

2845:   PetscFunctionBegin;
2847:   PetscAssertPointer(a, 6);
2849:   PetscCall(VecGetLocalSize(x, &N));
2850:   PetscCheck(m * n == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 2d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n);
2851:   PetscCall(VecGetArrayWrite(x, &aa));

2853:   PetscCall(PetscMalloc1(m, a));
2854:   for (i = 0; i < m; i++) (*a)[i] = aa + i * n - nstart;
2855:   *a -= mstart;
2856:   PetscFunctionReturn(PETSC_SUCCESS);
2857: }

2859: /*@C
2860:   VecRestoreArray2d - Restores a vector after `VecGetArray2d()` has been called.

2862:   Logically Collective

2864:   Input Parameters:
2865: + x      - the vector
2866: . m      - first dimension of two dimensional array
2867: . n      - second dimension of the two dimensional array
2868: . mstart - first index you will use in first coordinate direction (often 0)
2869: . nstart - first index in the second coordinate direction (often 0)
2870: - a      - location of pointer to array obtained from `VecGetArray2d()`

2872:   Level: developer

2874:   Notes:
2875:   For regular PETSc vectors this routine does not involve any copies. For
2876:   any special vectors that do not store local vector data in a contiguous
2877:   array, this routine will copy the data back into the underlying
2878:   vector data structure from the array obtained with `VecGetArray()`.

2880:   This routine actually zeros out the `a` pointer.

2882: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
2883:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
2884:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
2885: @*/
2886: PetscErrorCode VecRestoreArray2d(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
2887: {
2888:   void *dummy;

2890:   PetscFunctionBegin;
2892:   PetscAssertPointer(a, 6);
2894:   dummy = (void *)(*a + mstart);
2895:   PetscCall(PetscFree(dummy));
2896:   PetscCall(VecRestoreArray(x, NULL));
2897:   *a = NULL;
2898:   PetscFunctionReturn(PETSC_SUCCESS);
2899: }

2901: /*@C
2902:   VecRestoreArray2dWrite - Restores a vector after `VecGetArray2dWrite()` has been called.

2904:   Logically Collective

2906:   Input Parameters:
2907: + x      - the vector
2908: . m      - first dimension of two dimensional array
2909: . n      - second dimension of the two dimensional array
2910: . mstart - first index you will use in first coordinate direction (often 0)
2911: . nstart - first index in the second coordinate direction (often 0)
2912: - a      - location of pointer to array obtained from `VecGetArray2d()`

2914:   Level: developer

2916:   Notes:
2917:   For regular PETSc vectors this routine does not involve any copies. For
2918:   any special vectors that do not store local vector data in a contiguous
2919:   array, this routine will copy the data back into the underlying
2920:   vector data structure from the array obtained with `VecGetArray()`.

2922:   This routine actually zeros out the `a` pointer.

2924: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
2925:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
2926:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
2927: @*/
2928: PetscErrorCode VecRestoreArray2dWrite(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
2929: {
2930:   void *dummy;

2932:   PetscFunctionBegin;
2934:   PetscAssertPointer(a, 6);
2936:   dummy = (void *)(*a + mstart);
2937:   PetscCall(PetscFree(dummy));
2938:   PetscCall(VecRestoreArrayWrite(x, NULL));
2939:   PetscFunctionReturn(PETSC_SUCCESS);
2940: }

2942: /*@C
2943:   VecGetArray1d - Returns a pointer to a 1d contiguous array that contains this
2944:   processor's portion of the vector data.  You MUST call `VecRestoreArray1d()`
2945:   when you no longer need access to the array.

2947:   Logically Collective

2949:   Input Parameters:
2950: + x      - the vector
2951: . m      - first dimension of two dimensional array
2952: - mstart - first index you will use in first coordinate direction (often 0)

2954:   Output Parameter:
2955: . a - location to put pointer to the array

2957:   Level: developer

2959:   Notes:
2960:   For a vector obtained from `DMCreateLocalVector()` `mstart` is likely
2961:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
2962:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`.

2964:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

2966: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
2967:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
2968:           `VecGetArray2d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
2969: @*/
2970: PetscErrorCode VecGetArray1d(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
2971: {
2972:   PetscInt N;

2974:   PetscFunctionBegin;
2976:   PetscAssertPointer(a, 4);
2978:   PetscCall(VecGetLocalSize(x, &N));
2979:   PetscCheck(m == N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local array size %" PetscInt_FMT " does not match 1d array dimensions %" PetscInt_FMT, N, m);
2980:   PetscCall(VecGetArray(x, a));
2981:   *a -= mstart;
2982:   PetscFunctionReturn(PETSC_SUCCESS);
2983: }

2985: /*@C
2986:   VecGetArray1dWrite - Returns a pointer to a 1d contiguous array that will contain this
2987:   processor's portion of the vector data.  You MUST call `VecRestoreArray1dWrite()`
2988:   when you no longer need access to the array.

2990:   Logically Collective

2992:   Input Parameters:
2993: + x      - the vector
2994: . m      - first dimension of two dimensional array
2995: - mstart - first index you will use in first coordinate direction (often 0)

2997:   Output Parameter:
2998: . a - location to put pointer to the array

3000:   Level: developer

3002:   Notes:
3003:   For a vector obtained from `DMCreateLocalVector()` `mstart` is likely
3004:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3005:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`.

3007:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3009: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3010:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3011:           `VecGetArray2d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3012: @*/
3013: PetscErrorCode VecGetArray1dWrite(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
3014: {
3015:   PetscInt N;

3017:   PetscFunctionBegin;
3019:   PetscAssertPointer(a, 4);
3021:   PetscCall(VecGetLocalSize(x, &N));
3022:   PetscCheck(m == N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local array size %" PetscInt_FMT " does not match 1d array dimensions %" PetscInt_FMT, N, m);
3023:   PetscCall(VecGetArrayWrite(x, a));
3024:   *a -= mstart;
3025:   PetscFunctionReturn(PETSC_SUCCESS);
3026: }

3028: /*@C
3029:   VecRestoreArray1d - Restores a vector after `VecGetArray1d()` has been called.

3031:   Logically Collective

3033:   Input Parameters:
3034: + x      - the vector
3035: . m      - first dimension of two dimensional array
3036: . mstart - first index you will use in first coordinate direction (often 0)
3037: - a      - location of pointer to array obtained from `VecGetArray1d()`

3039:   Level: developer

3041:   Notes:
3042:   For regular PETSc vectors this routine does not involve any copies. For
3043:   any special vectors that do not store local vector data in a contiguous
3044:   array, this routine will copy the data back into the underlying
3045:   vector data structure from the array obtained with `VecGetArray1d()`.

3047:   This routine actually zeros out the `a` pointer.

3049: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3050:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3051:           `VecGetArray1d()`, `VecRestoreArray2d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3052: @*/
3053: PetscErrorCode VecRestoreArray1d(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
3054: {
3055:   PetscFunctionBegin;
3058:   PetscCall(VecRestoreArray(x, NULL));
3059:   *a = NULL;
3060:   PetscFunctionReturn(PETSC_SUCCESS);
3061: }

3063: /*@C
3064:   VecRestoreArray1dWrite - Restores a vector after `VecGetArray1dWrite()` has been called.

3066:   Logically Collective

3068:   Input Parameters:
3069: + x      - the vector
3070: . m      - first dimension of two dimensional array
3071: . mstart - first index you will use in first coordinate direction (often 0)
3072: - a      - location of pointer to array obtained from `VecGetArray1d()`

3074:   Level: developer

3076:   Notes:
3077:   For regular PETSc vectors this routine does not involve any copies. For
3078:   any special vectors that do not store local vector data in a contiguous
3079:   array, this routine will copy the data back into the underlying
3080:   vector data structure from the array obtained with `VecGetArray1d()`.

3082:   This routine actually zeros out the `a` pointer.

3084: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3085:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3086:           `VecGetArray1d()`, `VecRestoreArray2d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3087: @*/
3088: PetscErrorCode VecRestoreArray1dWrite(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
3089: {
3090:   PetscFunctionBegin;
3093:   PetscCall(VecRestoreArrayWrite(x, NULL));
3094:   *a = NULL;
3095:   PetscFunctionReturn(PETSC_SUCCESS);
3096: }

3098: /*@C
3099:   VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this
3100:   processor's portion of the vector data.  You MUST call `VecRestoreArray3d()`
3101:   when you no longer need access to the array.

3103:   Logically Collective

3105:   Input Parameters:
3106: + x      - the vector
3107: . m      - first dimension of three dimensional array
3108: . n      - second dimension of three dimensional array
3109: . p      - third dimension of three dimensional array
3110: . mstart - first index you will use in first coordinate direction (often 0)
3111: . nstart - first index in the second coordinate direction (often 0)
3112: - pstart - first index in the third coordinate direction (often 0)

3114:   Output Parameter:
3115: . a - location to put pointer to the array

3117:   Level: developer

3119:   Notes:
3120:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3121:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3122:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3123:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3d()`.

3125:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3127: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3128:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecRestoreArray3d()`,
3129:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3130: @*/
3131: PetscErrorCode VecGetArray3d(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3132: {
3133:   PetscInt     i, N, j;
3134:   PetscScalar *aa, **b;

3136:   PetscFunctionBegin;
3138:   PetscAssertPointer(a, 8);
3140:   PetscCall(VecGetLocalSize(x, &N));
3141:   PetscCheck(m * n * p == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 3d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p);
3142:   PetscCall(VecGetArray(x, &aa));

3144:   PetscCall(PetscMalloc(m * sizeof(PetscScalar **) + m * n * sizeof(PetscScalar *), a));
3145:   b = (PetscScalar **)((*a) + m);
3146:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3147:   for (i = 0; i < m; i++)
3148:     for (j = 0; j < n; j++) b[i * n + j] = PetscSafePointerPlusOffset(aa, i * n * p + j * p - pstart);
3149:   *a -= mstart;
3150:   PetscFunctionReturn(PETSC_SUCCESS);
3151: }

3153: /*@C
3154:   VecGetArray3dWrite - Returns a pointer to a 3d contiguous array that will contain this
3155:   processor's portion of the vector data.  You MUST call `VecRestoreArray3dWrite()`
3156:   when you no longer need access to the array.

3158:   Logically Collective

3160:   Input Parameters:
3161: + x      - the vector
3162: . m      - first dimension of three dimensional array
3163: . n      - second dimension of three dimensional array
3164: . p      - third dimension of three dimensional array
3165: . mstart - first index you will use in first coordinate direction (often 0)
3166: . nstart - first index in the second coordinate direction (often 0)
3167: - pstart - first index in the third coordinate direction (often 0)

3169:   Output Parameter:
3170: . a - location to put pointer to the array

3172:   Level: developer

3174:   Notes:
3175:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3176:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3177:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3178:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3d()`.

3180:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3182: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3183:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3184:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3185: @*/
3186: PetscErrorCode VecGetArray3dWrite(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3187: {
3188:   PetscInt     i, N, j;
3189:   PetscScalar *aa, **b;

3191:   PetscFunctionBegin;
3193:   PetscAssertPointer(a, 8);
3195:   PetscCall(VecGetLocalSize(x, &N));
3196:   PetscCheck(m * n * p == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 3d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p);
3197:   PetscCall(VecGetArrayWrite(x, &aa));

3199:   PetscCall(PetscMalloc(m * sizeof(PetscScalar **) + m * n * sizeof(PetscScalar *), a));
3200:   b = (PetscScalar **)((*a) + m);
3201:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3202:   for (i = 0; i < m; i++)
3203:     for (j = 0; j < n; j++) b[i * n + j] = aa + i * n * p + j * p - pstart;

3205:   *a -= mstart;
3206:   PetscFunctionReturn(PETSC_SUCCESS);
3207: }

3209: /*@C
3210:   VecRestoreArray3d - Restores a vector after `VecGetArray3d()` has been called.

3212:   Logically Collective

3214:   Input Parameters:
3215: + x      - the vector
3216: . m      - first dimension of three dimensional array
3217: . n      - second dimension of the three dimensional array
3218: . p      - third dimension of the three dimensional array
3219: . mstart - first index you will use in first coordinate direction (often 0)
3220: . nstart - first index in the second coordinate direction (often 0)
3221: . pstart - first index in the third coordinate direction (often 0)
3222: - a      - location of pointer to array obtained from VecGetArray3d()

3224:   Level: developer

3226:   Notes:
3227:   For regular PETSc vectors this routine does not involve any copies. For
3228:   any special vectors that do not store local vector data in a contiguous
3229:   array, this routine will copy the data back into the underlying
3230:   vector data structure from the array obtained with `VecGetArray()`.

3232:   This routine actually zeros out the `a` pointer.

3234: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3235:           `VecGetArray2d()`, `VecGetArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3236:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`, `VecGet`
3237: @*/
3238: PetscErrorCode VecRestoreArray3d(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3239: {
3240:   void *dummy;

3242:   PetscFunctionBegin;
3244:   PetscAssertPointer(a, 8);
3246:   dummy = (void *)(*a + mstart);
3247:   PetscCall(PetscFree(dummy));
3248:   PetscCall(VecRestoreArray(x, NULL));
3249:   *a = NULL;
3250:   PetscFunctionReturn(PETSC_SUCCESS);
3251: }

3253: /*@C
3254:   VecRestoreArray3dWrite - Restores a vector after `VecGetArray3dWrite()` has been called.

3256:   Logically Collective

3258:   Input Parameters:
3259: + x      - the vector
3260: . m      - first dimension of three dimensional array
3261: . n      - second dimension of the three dimensional array
3262: . p      - third dimension of the three dimensional array
3263: . mstart - first index you will use in first coordinate direction (often 0)
3264: . nstart - first index in the second coordinate direction (often 0)
3265: . pstart - first index in the third coordinate direction (often 0)
3266: - a      - location of pointer to array obtained from VecGetArray3d()

3268:   Level: developer

3270:   Notes:
3271:   For regular PETSc vectors this routine does not involve any copies. For
3272:   any special vectors that do not store local vector data in a contiguous
3273:   array, this routine will copy the data back into the underlying
3274:   vector data structure from the array obtained with `VecGetArray()`.

3276:   This routine actually zeros out the `a` pointer.

3278: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3279:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3280:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`, `VecGet`
3281: @*/
3282: PetscErrorCode VecRestoreArray3dWrite(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3283: {
3284:   void *dummy;

3286:   PetscFunctionBegin;
3288:   PetscAssertPointer(a, 8);
3290:   dummy = (void *)(*a + mstart);
3291:   PetscCall(PetscFree(dummy));
3292:   PetscCall(VecRestoreArrayWrite(x, NULL));
3293:   *a = NULL;
3294:   PetscFunctionReturn(PETSC_SUCCESS);
3295: }

3297: /*@C
3298:   VecGetArray4d - Returns a pointer to a 4d contiguous array that contains this
3299:   processor's portion of the vector data.  You MUST call `VecRestoreArray4d()`
3300:   when you no longer need access to the array.

3302:   Logically Collective

3304:   Input Parameters:
3305: + x      - the vector
3306: . m      - first dimension of four dimensional array
3307: . n      - second dimension of four dimensional array
3308: . p      - third dimension of four dimensional array
3309: . q      - fourth dimension of four dimensional array
3310: . mstart - first index you will use in first coordinate direction (often 0)
3311: . nstart - first index in the second coordinate direction (often 0)
3312: . pstart - first index in the third coordinate direction (often 0)
3313: - qstart - first index in the fourth coordinate direction (often 0)

3315:   Output Parameter:
3316: . a - location to put pointer to the array

3318:   Level: developer

3320:   Notes:
3321:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3322:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3323:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3324:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3d()`.

3326:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3328: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3329:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3330:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecRestoreArray4d()`
3331: @*/
3332: PetscErrorCode VecGetArray4d(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3333: {
3334:   PetscInt     i, N, j, k;
3335:   PetscScalar *aa, ***b, **c;

3337:   PetscFunctionBegin;
3339:   PetscAssertPointer(a, 10);
3341:   PetscCall(VecGetLocalSize(x, &N));
3342:   PetscCheck(m * n * p * q == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 4d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p, q);
3343:   PetscCall(VecGetArray(x, &aa));

3345:   PetscCall(PetscMalloc(m * sizeof(PetscScalar ***) + m * n * sizeof(PetscScalar **) + m * n * p * sizeof(PetscScalar *), a));
3346:   b = (PetscScalar ***)((*a) + m);
3347:   c = (PetscScalar **)(b + m * n);
3348:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3349:   for (i = 0; i < m; i++)
3350:     for (j = 0; j < n; j++) b[i * n + j] = c + i * n * p + j * p - pstart;
3351:   for (i = 0; i < m; i++)
3352:     for (j = 0; j < n; j++)
3353:       for (k = 0; k < p; k++) c[i * n * p + j * p + k] = aa + i * n * p * q + j * p * q + k * q - qstart;
3354:   *a -= mstart;
3355:   PetscFunctionReturn(PETSC_SUCCESS);
3356: }

3358: /*@C
3359:   VecGetArray4dWrite - Returns a pointer to a 4d contiguous array that will contain this
3360:   processor's portion of the vector data.  You MUST call `VecRestoreArray4dWrite()`
3361:   when you no longer need access to the array.

3363:   Logically Collective

3365:   Input Parameters:
3366: + x      - the vector
3367: . m      - first dimension of four dimensional array
3368: . n      - second dimension of four dimensional array
3369: . p      - third dimension of four dimensional array
3370: . q      - fourth dimension of four dimensional array
3371: . mstart - first index you will use in first coordinate direction (often 0)
3372: . nstart - first index in the second coordinate direction (often 0)
3373: . pstart - first index in the third coordinate direction (often 0)
3374: - qstart - first index in the fourth coordinate direction (often 0)

3376:   Output Parameter:
3377: . a - location to put pointer to the array

3379:   Level: developer

3381:   Notes:
3382:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3383:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3384:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3385:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3d()`.

3387:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3389: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3390:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3391:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3392: @*/
3393: PetscErrorCode VecGetArray4dWrite(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3394: {
3395:   PetscInt     i, N, j, k;
3396:   PetscScalar *aa, ***b, **c;

3398:   PetscFunctionBegin;
3400:   PetscAssertPointer(a, 10);
3402:   PetscCall(VecGetLocalSize(x, &N));
3403:   PetscCheck(m * n * p * q == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 4d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p, q);
3404:   PetscCall(VecGetArrayWrite(x, &aa));

3406:   PetscCall(PetscMalloc(m * sizeof(PetscScalar ***) + m * n * sizeof(PetscScalar **) + m * n * p * sizeof(PetscScalar *), a));
3407:   b = (PetscScalar ***)((*a) + m);
3408:   c = (PetscScalar **)(b + m * n);
3409:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3410:   for (i = 0; i < m; i++)
3411:     for (j = 0; j < n; j++) b[i * n + j] = c + i * n * p + j * p - pstart;
3412:   for (i = 0; i < m; i++)
3413:     for (j = 0; j < n; j++)
3414:       for (k = 0; k < p; k++) c[i * n * p + j * p + k] = aa + i * n * p * q + j * p * q + k * q - qstart;
3415:   *a -= mstart;
3416:   PetscFunctionReturn(PETSC_SUCCESS);
3417: }

3419: /*@C
3420:   VecRestoreArray4d - Restores a vector after `VecGetArray4d()` has been called.

3422:   Logically Collective

3424:   Input Parameters:
3425: + x      - the vector
3426: . m      - first dimension of four dimensional array
3427: . n      - second dimension of the four dimensional array
3428: . p      - third dimension of the four dimensional array
3429: . q      - fourth dimension of the four dimensional array
3430: . mstart - first index you will use in first coordinate direction (often 0)
3431: . nstart - first index in the second coordinate direction (often 0)
3432: . pstart - first index in the third coordinate direction (often 0)
3433: . qstart - first index in the fourth coordinate direction (often 0)
3434: - a      - location of pointer to array obtained from VecGetArray4d()

3436:   Level: developer

3438:   Notes:
3439:   For regular PETSc vectors this routine does not involve any copies. For
3440:   any special vectors that do not store local vector data in a contiguous
3441:   array, this routine will copy the data back into the underlying
3442:   vector data structure from the array obtained with `VecGetArray()`.

3444:   This routine actually zeros out the `a` pointer.

3446: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3447:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3448:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecGet`
3449: @*/
3450: PetscErrorCode VecRestoreArray4d(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3451: {
3452:   void *dummy;

3454:   PetscFunctionBegin;
3456:   PetscAssertPointer(a, 10);
3458:   dummy = (void *)(*a + mstart);
3459:   PetscCall(PetscFree(dummy));
3460:   PetscCall(VecRestoreArray(x, NULL));
3461:   *a = NULL;
3462:   PetscFunctionReturn(PETSC_SUCCESS);
3463: }

3465: /*@C
3466:   VecRestoreArray4dWrite - Restores a vector after `VecGetArray4dWrite()` has been called.

3468:   Logically Collective

3470:   Input Parameters:
3471: + x      - the vector
3472: . m      - first dimension of four dimensional array
3473: . n      - second dimension of the four dimensional array
3474: . p      - third dimension of the four dimensional array
3475: . q      - fourth dimension of the four dimensional array
3476: . mstart - first index you will use in first coordinate direction (often 0)
3477: . nstart - first index in the second coordinate direction (often 0)
3478: . pstart - first index in the third coordinate direction (often 0)
3479: . qstart - first index in the fourth coordinate direction (often 0)
3480: - a      - location of pointer to array obtained from `VecGetArray4d()`

3482:   Level: developer

3484:   Notes:
3485:   For regular PETSc vectors this routine does not involve any copies. For
3486:   any special vectors that do not store local vector data in a contiguous
3487:   array, this routine will copy the data back into the underlying
3488:   vector data structure from the array obtained with `VecGetArray()`.

3490:   This routine actually zeros out the `a` pointer.

3492: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3493:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3494:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`, `VecGet`
3495: @*/
3496: PetscErrorCode VecRestoreArray4dWrite(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3497: {
3498:   void *dummy;

3500:   PetscFunctionBegin;
3502:   PetscAssertPointer(a, 10);
3504:   dummy = (void *)(*a + mstart);
3505:   PetscCall(PetscFree(dummy));
3506:   PetscCall(VecRestoreArrayWrite(x, NULL));
3507:   *a = NULL;
3508:   PetscFunctionReturn(PETSC_SUCCESS);
3509: }

3511: /*@C
3512:   VecGetArray2dRead - Returns a pointer to a 2d contiguous array that contains this
3513:   processor's portion of the vector data.  You MUST call `VecRestoreArray2dRead()`
3514:   when you no longer need access to the array.

3516:   Logically Collective

3518:   Input Parameters:
3519: + x      - the vector
3520: . m      - first dimension of two dimensional array
3521: . n      - second dimension of two dimensional array
3522: . mstart - first index you will use in first coordinate direction (often 0)
3523: - nstart - first index in the second coordinate direction (often 0)

3525:   Output Parameter:
3526: . a - location to put pointer to the array

3528:   Level: developer

3530:   Notes:
3531:   For a vector obtained from `DMCreateLocalVector()` `mstart` and `nstart` are likely
3532:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3533:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3534:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray2d()`.

3536:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3538: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3539:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3540:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3541: @*/
3542: PetscErrorCode VecGetArray2dRead(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
3543: {
3544:   PetscInt           i, N;
3545:   const PetscScalar *aa;

3547:   PetscFunctionBegin;
3549:   PetscAssertPointer(a, 6);
3551:   PetscCall(VecGetLocalSize(x, &N));
3552:   PetscCheck(m * n == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 2d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n);
3553:   PetscCall(VecGetArrayRead(x, &aa));

3555:   PetscCall(PetscMalloc1(m, a));
3556:   for (i = 0; i < m; i++) (*a)[i] = (PetscScalar *)aa + i * n - nstart;
3557:   *a -= mstart;
3558:   PetscFunctionReturn(PETSC_SUCCESS);
3559: }

3561: /*@C
3562:   VecRestoreArray2dRead - Restores a vector after `VecGetArray2dRead()` has been called.

3564:   Logically Collective

3566:   Input Parameters:
3567: + x      - the vector
3568: . m      - first dimension of two dimensional array
3569: . n      - second dimension of the two dimensional array
3570: . mstart - first index you will use in first coordinate direction (often 0)
3571: . nstart - first index in the second coordinate direction (often 0)
3572: - a      - location of pointer to array obtained from VecGetArray2d()

3574:   Level: developer

3576:   Notes:
3577:   For regular PETSc vectors this routine does not involve any copies. For
3578:   any special vectors that do not store local vector data in a contiguous
3579:   array, this routine will copy the data back into the underlying
3580:   vector data structure from the array obtained with `VecGetArray()`.

3582:   This routine actually zeros out the `a` pointer.

3584: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3585:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3586:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3587: @*/
3588: PetscErrorCode VecRestoreArray2dRead(Vec x, PetscInt m, PetscInt n, PetscInt mstart, PetscInt nstart, PetscScalar **a[])
3589: {
3590:   void *dummy;

3592:   PetscFunctionBegin;
3594:   PetscAssertPointer(a, 6);
3596:   dummy = (void *)(*a + mstart);
3597:   PetscCall(PetscFree(dummy));
3598:   PetscCall(VecRestoreArrayRead(x, NULL));
3599:   *a = NULL;
3600:   PetscFunctionReturn(PETSC_SUCCESS);
3601: }

3603: /*@C
3604:   VecGetArray1dRead - Returns a pointer to a 1d contiguous array that contains this
3605:   processor's portion of the vector data.  You MUST call `VecRestoreArray1dRead()`
3606:   when you no longer need access to the array.

3608:   Logically Collective

3610:   Input Parameters:
3611: + x      - the vector
3612: . m      - first dimension of two dimensional array
3613: - mstart - first index you will use in first coordinate direction (often 0)

3615:   Output Parameter:
3616: . a - location to put pointer to the array

3618:   Level: developer

3620:   Notes:
3621:   For a vector obtained from `DMCreateLocalVector()` `mstart` is likely
3622:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3623:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`.

3625:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3627: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3628:           `VecRestoreArray2d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3629:           `VecGetArray2d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3630: @*/
3631: PetscErrorCode VecGetArray1dRead(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
3632: {
3633:   PetscInt N;

3635:   PetscFunctionBegin;
3637:   PetscAssertPointer(a, 4);
3639:   PetscCall(VecGetLocalSize(x, &N));
3640:   PetscCheck(m == N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Local array size %" PetscInt_FMT " does not match 1d array dimensions %" PetscInt_FMT, N, m);
3641:   PetscCall(VecGetArrayRead(x, (const PetscScalar **)a));
3642:   *a -= mstart;
3643:   PetscFunctionReturn(PETSC_SUCCESS);
3644: }

3646: /*@C
3647:   VecRestoreArray1dRead - Restores a vector after `VecGetArray1dRead()` has been called.

3649:   Logically Collective

3651:   Input Parameters:
3652: + x      - the vector
3653: . m      - first dimension of two dimensional array
3654: . mstart - first index you will use in first coordinate direction (often 0)
3655: - a      - location of pointer to array obtained from `VecGetArray1dRead()`

3657:   Level: developer

3659:   Notes:
3660:   For regular PETSc vectors this routine does not involve any copies. For
3661:   any special vectors that do not store local vector data in a contiguous
3662:   array, this routine will copy the data back into the underlying
3663:   vector data structure from the array obtained with `VecGetArray1dRead()`.

3665:   This routine actually zeros out the `a` pointer.

3667: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3668:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3669:           `VecGetArray1d()`, `VecRestoreArray2d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3670: @*/
3671: PetscErrorCode VecRestoreArray1dRead(Vec x, PetscInt m, PetscInt mstart, PetscScalar *a[])
3672: {
3673:   PetscFunctionBegin;
3676:   PetscCall(VecRestoreArrayRead(x, NULL));
3677:   *a = NULL;
3678:   PetscFunctionReturn(PETSC_SUCCESS);
3679: }

3681: /*@C
3682:   VecGetArray3dRead - Returns a pointer to a 3d contiguous array that contains this
3683:   processor's portion of the vector data.  You MUST call `VecRestoreArray3dRead()`
3684:   when you no longer need access to the array.

3686:   Logically Collective

3688:   Input Parameters:
3689: + x      - the vector
3690: . m      - first dimension of three dimensional array
3691: . n      - second dimension of three dimensional array
3692: . p      - third dimension of three dimensional array
3693: . mstart - first index you will use in first coordinate direction (often 0)
3694: . nstart - first index in the second coordinate direction (often 0)
3695: - pstart - first index in the third coordinate direction (often 0)

3697:   Output Parameter:
3698: . a - location to put pointer to the array

3700:   Level: developer

3702:   Notes:
3703:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3704:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3705:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3706:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3dRead()`.

3708:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3710: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3711:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3712:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3713: @*/
3714: PetscErrorCode VecGetArray3dRead(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3715: {
3716:   PetscInt           i, N, j;
3717:   const PetscScalar *aa;
3718:   PetscScalar      **b;

3720:   PetscFunctionBegin;
3722:   PetscAssertPointer(a, 8);
3724:   PetscCall(VecGetLocalSize(x, &N));
3725:   PetscCheck(m * n * p == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 3d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p);
3726:   PetscCall(VecGetArrayRead(x, &aa));

3728:   PetscCall(PetscMalloc(m * sizeof(PetscScalar **) + m * n * sizeof(PetscScalar *), a));
3729:   b = (PetscScalar **)((*a) + m);
3730:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3731:   for (i = 0; i < m; i++)
3732:     for (j = 0; j < n; j++) b[i * n + j] = PetscSafePointerPlusOffset((PetscScalar *)aa, i * n * p + j * p - pstart);
3733:   *a -= mstart;
3734:   PetscFunctionReturn(PETSC_SUCCESS);
3735: }

3737: /*@C
3738:   VecRestoreArray3dRead - Restores a vector after `VecGetArray3dRead()` has been called.

3740:   Logically Collective

3742:   Input Parameters:
3743: + x      - the vector
3744: . m      - first dimension of three dimensional array
3745: . n      - second dimension of the three dimensional array
3746: . p      - third dimension of the three dimensional array
3747: . mstart - first index you will use in first coordinate direction (often 0)
3748: . nstart - first index in the second coordinate direction (often 0)
3749: . pstart - first index in the third coordinate direction (often 0)
3750: - a      - location of pointer to array obtained from `VecGetArray3dRead()`

3752:   Level: developer

3754:   Notes:
3755:   For regular PETSc vectors this routine does not involve any copies. For
3756:   any special vectors that do not store local vector data in a contiguous
3757:   array, this routine will copy the data back into the underlying
3758:   vector data structure from the array obtained with `VecGetArray()`.

3760:   This routine actually zeros out the `a` pointer.

3762: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3763:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3764:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`, `VecGet`
3765: @*/
3766: PetscErrorCode VecRestoreArray3dRead(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscScalar ***a[])
3767: {
3768:   void *dummy;

3770:   PetscFunctionBegin;
3772:   PetscAssertPointer(a, 8);
3774:   dummy = (void *)(*a + mstart);
3775:   PetscCall(PetscFree(dummy));
3776:   PetscCall(VecRestoreArrayRead(x, NULL));
3777:   *a = NULL;
3778:   PetscFunctionReturn(PETSC_SUCCESS);
3779: }

3781: /*@C
3782:   VecGetArray4dRead - Returns a pointer to a 4d contiguous array that contains this
3783:   processor's portion of the vector data.  You MUST call `VecRestoreArray4dRead()`
3784:   when you no longer need access to the array.

3786:   Logically Collective

3788:   Input Parameters:
3789: + x      - the vector
3790: . m      - first dimension of four dimensional array
3791: . n      - second dimension of four dimensional array
3792: . p      - third dimension of four dimensional array
3793: . q      - fourth dimension of four dimensional array
3794: . mstart - first index you will use in first coordinate direction (often 0)
3795: . nstart - first index in the second coordinate direction (often 0)
3796: . pstart - first index in the third coordinate direction (often 0)
3797: - qstart - first index in the fourth coordinate direction (often 0)

3799:   Output Parameter:
3800: . a - location to put pointer to the array

3802:   Level: beginner

3804:   Notes:
3805:   For a vector obtained from `DMCreateLocalVector()` `mstart`, `nstart`, and `pstart` are likely
3806:   obtained from the corner indices obtained from `DMDAGetGhostCorners()` while for
3807:   `DMCreateGlobalVector()` they are the corner indices from `DMDAGetCorners()`. In both cases
3808:   the arguments from `DMDAGet[Ghost]Corners()` are reversed in the call to `VecGetArray3d()`.

3810:   For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.

3812: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecGetArrays()`, `VecGetArrayF90()`, `VecPlaceArray()`,
3813:           `VecRestoreArray2d()`, `DMDAVecGetarray()`, `DMDAVecRestoreArray()`, `VecGetArray3d()`, `VecRestoreArray3d()`,
3814:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`
3815: @*/
3816: PetscErrorCode VecGetArray4dRead(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3817: {
3818:   PetscInt           i, N, j, k;
3819:   const PetscScalar *aa;
3820:   PetscScalar     ***b, **c;

3822:   PetscFunctionBegin;
3824:   PetscAssertPointer(a, 10);
3826:   PetscCall(VecGetLocalSize(x, &N));
3827:   PetscCheck(m * n * p * q == N, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local array size %" PetscInt_FMT " does not match 4d array dimensions %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT " by %" PetscInt_FMT, N, m, n, p, q);
3828:   PetscCall(VecGetArrayRead(x, &aa));

3830:   PetscCall(PetscMalloc(m * sizeof(PetscScalar ***) + m * n * sizeof(PetscScalar **) + m * n * p * sizeof(PetscScalar *), a));
3831:   b = (PetscScalar ***)((*a) + m);
3832:   c = (PetscScalar **)(b + m * n);
3833:   for (i = 0; i < m; i++) (*a)[i] = b + i * n - nstart;
3834:   for (i = 0; i < m; i++)
3835:     for (j = 0; j < n; j++) b[i * n + j] = c + i * n * p + j * p - pstart;
3836:   for (i = 0; i < m; i++)
3837:     for (j = 0; j < n; j++)
3838:       for (k = 0; k < p; k++) c[i * n * p + j * p + k] = (PetscScalar *)aa + i * n * p * q + j * p * q + k * q - qstart;
3839:   *a -= mstart;
3840:   PetscFunctionReturn(PETSC_SUCCESS);
3841: }

3843: /*@C
3844:   VecRestoreArray4dRead - Restores a vector after `VecGetArray4d()` has been called.

3846:   Logically Collective

3848:   Input Parameters:
3849: + x      - the vector
3850: . m      - first dimension of four dimensional array
3851: . n      - second dimension of the four dimensional array
3852: . p      - third dimension of the four dimensional array
3853: . q      - fourth dimension of the four dimensional array
3854: . mstart - first index you will use in first coordinate direction (often 0)
3855: . nstart - first index in the second coordinate direction (often 0)
3856: . pstart - first index in the third coordinate direction (often 0)
3857: . qstart - first index in the fourth coordinate direction (often 0)
3858: - a      - location of pointer to array obtained from `VecGetArray4dRead()`

3860:   Level: beginner

3862:   Notes:
3863:   For regular PETSc vectors this routine does not involve any copies. For
3864:   any special vectors that do not store local vector data in a contiguous
3865:   array, this routine will copy the data back into the underlying
3866:   vector data structure from the array obtained with `VecGetArray()`.

3868:   This routine actually zeros out the `a` pointer.

3870: .seealso: [](ch_vectors), `Vec`, `VecGetArray()`, `VecRestoreArray()`, `VecRestoreArrays()`, `VecRestoreArrayF90()`, `VecPlaceArray()`,
3871:           `VecGetArray2d()`, `VecGetArray3d()`, `VecRestoreArray3d()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`
3872:           `VecGetArray1d()`, `VecRestoreArray1d()`, `VecGetArray4d()`, `VecRestoreArray4d()`, `VecGet`
3873: @*/
3874: PetscErrorCode VecRestoreArray4dRead(Vec x, PetscInt m, PetscInt n, PetscInt p, PetscInt q, PetscInt mstart, PetscInt nstart, PetscInt pstart, PetscInt qstart, PetscScalar ****a[])
3875: {
3876:   void *dummy;

3878:   PetscFunctionBegin;
3880:   PetscAssertPointer(a, 10);
3882:   dummy = (void *)(*a + mstart);
3883:   PetscCall(PetscFree(dummy));
3884:   PetscCall(VecRestoreArrayRead(x, NULL));
3885:   *a = NULL;
3886:   PetscFunctionReturn(PETSC_SUCCESS);
3887: }

3889: #if defined(PETSC_USE_DEBUG)

3891: /*@
3892:   VecLockGet  - Gets the current lock status of a vector

3894:   Logically Collective

3896:   Input Parameter:
3897: . x - the vector

3899:   Output Parameter:
3900: . state - greater than zero indicates the vector is locked for read; less than zero indicates the vector is
3901:            locked for write; equal to zero means the vector is unlocked, that is, it is free to read or write.

3903:   Level: advanced

3905: .seealso: [](ch_vectors), `Vec`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecLockReadPush()`, `VecLockReadPop()`
3906: @*/
3907: PetscErrorCode VecLockGet(Vec x, PetscInt *state)
3908: {
3909:   PetscFunctionBegin;
3911:   *state = x->lock;
3912:   PetscFunctionReturn(PETSC_SUCCESS);
3913: }

3915: PetscErrorCode VecLockGetLocation(Vec x, const char *file[], const char *func[], int *line)
3916: {
3917:   PetscFunctionBegin;
3919:   PetscAssertPointer(file, 2);
3920:   PetscAssertPointer(func, 3);
3921:   PetscAssertPointer(line, 4);
3922:   #if !PetscDefined(HAVE_THREADSAFETY)
3923:   {
3924:     const int index = x->lockstack.currentsize - 1;

3926:     PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Corrupted vec lock stack, have negative index %d", index);
3927:     *file = x->lockstack.file[index];
3928:     *func = x->lockstack.function[index];
3929:     *line = x->lockstack.line[index];
3930:   }
3931:   #else
3932:   *file = NULL;
3933:   *func = NULL;
3934:   *line = 0;
3935:   #endif
3936:   PetscFunctionReturn(PETSC_SUCCESS);
3937: }

3939: /*@
3940:   VecLockReadPush  - Pushes a read-only lock on a vector to prevent it from being written to

3942:   Logically Collective

3944:   Input Parameter:
3945: . x - the vector

3947:   Level: intermediate

3949:   Notes:
3950:   If this is set then calls to `VecGetArray()` or `VecSetValues()` or any other routines that change the vectors values will generate an error.

3952:   The call can be nested, i.e., called multiple times on the same vector, but each `VecLockReadPush()` has to have one matching
3953:   `VecLockReadPop()`, which removes the latest read-only lock.

3955: .seealso: [](ch_vectors), `Vec`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecLockReadPop()`, `VecLockGet()`
3956: @*/
3957: PetscErrorCode VecLockReadPush(Vec x)
3958: {
3959:   PetscFunctionBegin;
3961:   PetscCheck(x->lock++ >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector is already locked for exclusive write access but you want to read it");
3962:   #if !PetscDefined(HAVE_THREADSAFETY)
3963:   {
3964:     const char *file, *func;
3965:     int         index, line;

3967:     if ((index = petscstack.currentsize - 2) == -1) {
3968:       // vec was locked "outside" of petsc, either in user-land or main. the error message will
3969:       // now show this function as the culprit, but it will include the stacktrace
3970:       file = "unknown user-file";
3971:       func = "unknown_user_function";
3972:       line = 0;
3973:     } else {
3974:       PetscCheck(index >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected petscstack, have negative index %d", index);
3975:       file = petscstack.file[index];
3976:       func = petscstack.function[index];
3977:       line = petscstack.line[index];
3978:     }
3979:     PetscStackPush_Private(x->lockstack, file, func, line, petscstack.petscroutine[index], PETSC_FALSE);
3980:   }
3981:   #endif
3982:   PetscFunctionReturn(PETSC_SUCCESS);
3983: }

3985: /*@
3986:   VecLockReadPop  - Pops a read-only lock from a vector

3988:   Logically Collective

3990:   Input Parameter:
3991: . x - the vector

3993:   Level: intermediate

3995: .seealso: [](ch_vectors), `Vec`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecLockReadPush()`, `VecLockGet()`
3996: @*/
3997: PetscErrorCode VecLockReadPop(Vec x)
3998: {
3999:   PetscFunctionBegin;
4001:   PetscCheck(--x->lock >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector has been unlocked from read-only access too many times");
4002:   #if !PetscDefined(HAVE_THREADSAFETY)
4003:   {
4004:     const char *previous = x->lockstack.function[x->lockstack.currentsize - 1];

4006:     PetscStackPop_Private(x->lockstack, previous);
4007:   }
4008:   #endif
4009:   PetscFunctionReturn(PETSC_SUCCESS);
4010: }

4012: /*@C
4013:   VecLockWriteSet  - Lock or unlock a vector for exclusive read/write access

4015:   Logically Collective

4017:   Input Parameters:
4018: + x   - the vector
4019: - flg - `PETSC_TRUE` to lock the vector for exclusive read/write access; `PETSC_FALSE` to unlock it.

4021:   Level: intermediate

4023:   Notes:
4024:   The function is useful in split-phase computations, which usually have a begin phase and an end phase.
4025:   One can call `VecLockWriteSet`(x,`PETSC_TRUE`) in the begin phase to lock a vector for exclusive
4026:   access, and call `VecLockWriteSet`(x,`PETSC_FALSE`) in the end phase to unlock the vector from exclusive
4027:   access. In this way, one is ensured no other operations can access the vector in between. The code may like

4029: .vb
4030:        VecGetArray(x,&xdata); // begin phase
4031:        VecLockWriteSet(v,PETSC_TRUE);

4033:        Other operations, which can not access x anymore (they can access xdata, of course)

4035:        VecRestoreArray(x,&vdata); // end phase
4036:        VecLockWriteSet(v,PETSC_FALSE);
4037: .ve

4039:   The call can not be nested on the same vector, in other words, one can not call `VecLockWriteSet`(x,`PETSC_TRUE`)
4040:   again before calling `VecLockWriteSet`(v,`PETSC_FALSE`).

4042: .seealso: [](ch_vectors), `Vec`, `VecRestoreArray()`, `VecGetArrayRead()`, `VecLockReadPush()`, `VecLockReadPop()`, `VecLockGet()`
4043: @*/
4044: PetscErrorCode VecLockWriteSet(Vec x, PetscBool flg)
4045: {
4046:   PetscFunctionBegin;
4048:   if (flg) {
4049:     PetscCheck(x->lock <= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector is already locked for read-only access but you want to write it");
4050:     PetscCheck(x->lock >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector is already locked for exclusive write access but you want to write it");
4051:     x->lock = -1;
4052:   } else {
4053:     PetscCheck(x->lock == -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vector is not locked for exclusive write access but you want to unlock it from that");
4054:     x->lock = 0;
4055:   }
4056:   PetscFunctionReturn(PETSC_SUCCESS);
4057: }

4059: // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
4060: /*@
4061:   VecLockPush  - Pushes a read-only lock on a vector to prevent it from being written to

4063:   Level: deprecated

4065: .seealso: [](ch_vectors), `Vec`, `VecLockReadPush()`
4066: @*/
4067: PetscErrorCode VecLockPush(Vec x)
4068: {
4069:   PetscFunctionBegin;
4070:   PetscCall(VecLockReadPush(x));
4071:   PetscFunctionReturn(PETSC_SUCCESS);
4072: }

4074: // PetscClangLinter pragma disable: -fdoc-param-list-func-parameter-documentation
4075: /*@
4076:   VecLockPop  - Pops a read-only lock from a vector

4078:   Level: deprecated

4080: .seealso: [](ch_vectors), `Vec`, `VecLockReadPop()`
4081: @*/
4082: PetscErrorCode VecLockPop(Vec x)
4083: {
4084:   PetscFunctionBegin;
4085:   PetscCall(VecLockReadPop(x));
4086:   PetscFunctionReturn(PETSC_SUCCESS);
4087: }

4089: #endif