Actual source code: rvector.c

petsc-master 2018-06-24
Report Typos and Errors

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

  9: PETSC_EXTERN PetscErrorCode VecValidValues(Vec vec,PetscInt argnum,PetscBool begin)
 10: {
 11: #if defined(PETSC_USE_DEBUG)
 12:   PetscErrorCode    ierr;
 13:   PetscInt          n,i;
 14:   const PetscScalar *x;

 17: #if defined(PETSC_HAVE_VECCUDA) || defined(PETSC_HAVE_VIENNACL)
 18:   if ((vec->petscnative || vec->ops->getarray) && (vec->valid_GPU_array == PETSC_OFFLOAD_CPU || vec->valid_GPU_array == PETSC_OFFLOAD_BOTH)) {
 19: #else
 20:   if (vec->petscnative || vec->ops->getarray) {
 21: #endif
 22:     VecGetLocalSize(vec,&n);
 23:     VecGetArrayRead(vec,&x);
 24:     for (i=0; i<n; i++) {
 25:       if (begin) {
 26:         if (PetscIsInfOrNanScalar(x[i])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FP,"Vec entry at local location %D is not-a-number or infinite at beginning of function: Parameter number %D",i,argnum);
 27:       } else {
 28:         if (PetscIsInfOrNanScalar(x[i])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FP,"Vec entry at local location %D is not-a-number or infinite at end of function: Parameter number %D",i,argnum);
 29:       }
 30:     }
 31:     VecRestoreArrayRead(vec,&x);
 32:   }
 33:   return(0);
 34: #else
 35:   return 0;
 36: #endif
 37: }

 39: /*@
 40:    VecMaxPointwiseDivide - Computes the maximum of the componentwise division max = max_i abs(x_i/y_i).

 42:    Logically Collective on Vec

 44:    Input Parameters:
 45: .  x, y  - the vectors

 47:    Output Parameter:
 48: .  max - the result

 50:    Level: advanced

 52:    Notes:
 53:     x and y may be the same vector
 54:           if a particular y_i is zero, it is treated as 1 in the above formula

 56: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMax(), VecPointwiseMin(), VecPointwiseMaxAbs()
 57: @*/
 58: PetscErrorCode  VecMaxPointwiseDivide(Vec x,Vec y,PetscReal *max)
 59: {

 69:   VecCheckSameSize(x,1,y,2);
 70:   (*x->ops->maxpointwisedivide)(x,y,max);
 71:   return(0);
 72: }

 74: /*@
 75:    VecDot - Computes the vector dot product.

 77:    Collective on Vec

 79:    Input Parameters:
 80: .  x, y - the vectors

 82:    Output Parameter:
 83: .  val - the dot product

 85:    Performance Issues:
 86: $    per-processor memory bandwidth
 87: $    interprocessor latency
 88: $    work load inbalance that causes certain processes to arrive much earlier than others

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

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

101:    Level: intermediate

103:    Concepts: inner product
104:    Concepts: vector^inner product

106: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDotRealPart()
107: @*/
108: PetscErrorCode  VecDot(Vec x,Vec y,PetscScalar *val)
109: {

119:   VecCheckSameSize(x,1,y,2);

121:   PetscLogEventBegin(VEC_Dot,x,y,0,0);
122:   (*x->ops->dot)(x,y,val);
123:   PetscLogEventEnd(VEC_Dot,x,y,0,0);
124:   return(0);
125: }

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

130:    Collective on Vec

132:    Input Parameters:
133: .  x, y - the vectors

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

138:    Performance Issues:
139: $    per-processor memory bandwidth
140: $    interprocessor latency
141: $    work load inbalance that causes certain processes to arrive much earlier than others

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

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

148:      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
149:      the space R^{2n} (that is a vector of 2n components with the real or imaginary part of the complex numbers for components)

151:    Developer Note: This is not currently optimized to compute only the real part of the dot product.

153:    Level: intermediate

155:    Concepts: inner product
156:    Concepts: vector^inner product

158: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDot(), VecDotNorm2()
159: @*/
160: PetscErrorCode  VecDotRealPart(Vec x,Vec y,PetscReal *val)
161: {
163:   PetscScalar    fdot;

166:   VecDot(x,y,&fdot);
167:   *val = PetscRealPart(fdot);
168:   return(0);
169: }

171: /*@
172:    VecNorm  - Computes the vector norm.

174:    Collective on Vec

176:    Input Parameters:
177: +  x - the vector
178: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
179:           NORM_1_AND_2, which computes both norms and stores them
180:           in a two element array.

182:    Output Parameter:
183: .  val - the norm

185:    Notes:
186: $     NORM_1 denotes sum_i |x_i|
187: $     NORM_2 denotes sqrt(sum_i |x_i|^2)
188: $     NORM_INFINITY denotes max_i |x_i|

190:       For complex numbers NORM_1 will return the traditional 1 norm of the 2 norm of the complex numbers; that is the 1
191:       norm of the absolutely values of the complex entries. In PETSc 3.6 and earlier releases it returned the 1 norm of
192:       the 1 norm of the complex entries (what is returned by the BLAS routine asum()). Both are valid norms but most
193:       people expect the former.

195:    Level: intermediate

197:    Performance Issues:
198: $    per-processor memory bandwidth
199: $    interprocessor latency
200: $    work load inbalance that causes certain processes to arrive much earlier than others

202:    Concepts: norm
203:    Concepts: vector^norm

205: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNormAvailable(),
206:           VecNormBegin(), VecNormEnd()

208: @*/
209: PetscErrorCode  VecNorm(Vec x,NormType type,PetscReal *val)
210: {
211:   PetscBool      flg;


219:   /*
220:    * Cached data?
221:    */
222:   if (type!=NORM_1_AND_2) {
223:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,flg);
224:     if (flg) return(0);
225:   }
226:   PetscLogEventBegin(VEC_Norm,x,0,0,0);
227:   (*x->ops->norm)(x,type,val);
228:   PetscLogEventEnd(VEC_Norm,x,0,0,0);

230:   if (type!=NORM_1_AND_2) {
231:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[type],*val);
232:   }
233:   return(0);
234: }

236: /*@
237:    VecNormAvailable  - Returns the vector norm if it is already known.

239:    Not Collective

241:    Input Parameters:
242: +  x - the vector
243: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
244:           NORM_1_AND_2, which computes both norms and stores them
245:           in a two element array.

247:    Output Parameter:
248: +  available - PETSC_TRUE if the val returned is valid
249: -  val - the norm

251:    Notes:
252: $     NORM_1 denotes sum_i |x_i|
253: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
254: $     NORM_INFINITY denotes max_i |x_i|

256:    Level: intermediate

258:    Performance Issues:
259: $    per-processor memory bandwidth
260: $    interprocessor latency
261: $    work load inbalance that causes certain processes to arrive much earlier than others

263:    Compile Option:
264:    PETSC_HAVE_SLOW_BLAS_NORM2 will cause a C (loop unrolled) version of the norm to be used, rather
265:  than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines
266:  (as opposed to vendor provided) because the FORTRAN BLAS NRM2() routine is very slow.

268:    Concepts: norm
269:    Concepts: vector^norm

271: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNorm()
272:           VecNormBegin(), VecNormEnd()

274: @*/
275: PetscErrorCode  VecNormAvailable(Vec x,NormType type,PetscBool  *available,PetscReal *val)
276: {


284:   *available = PETSC_FALSE;
285:   if (type!=NORM_1_AND_2) {
286:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,*available);
287:   }
288:   return(0);
289: }

291: /*@
292:    VecNormalize - Normalizes a vector by 2-norm.

294:    Collective on Vec

296:    Input Parameters:
297: +  x - the vector

299:    Output Parameter:
300: .  x - the normalized vector
301: -  val - the vector norm before normalization

303:    Level: intermediate

305:    Concepts: vector^normalizing
306:    Concepts: normalizing^vector

308: @*/
309: PetscErrorCode  VecNormalize(Vec x,PetscReal *val)
310: {
312:   PetscReal      norm;

317:   PetscLogEventBegin(VEC_Normalize,x,0,0,0);
318:   VecNorm(x,NORM_2,&norm);
319:   if (norm == 0.0) {
320:     PetscInfo(x,"Vector of zero norm can not be normalized; Returning only the zero norm\n");
321:   } else if (norm != 1.0) {
322:     PetscScalar tmp = 1.0/norm;
323:     VecScale(x,tmp);
324:   }
325:   if (val) *val = norm;
326:   PetscLogEventEnd(VEC_Normalize,x,0,0,0);
327:   return(0);
328: }

330: /*@C
331:    VecMax - Determines the maximum vector component and its location.

333:    Collective on Vec

335:    Input Parameter:
336: .  x - the vector

338:    Output Parameters:
339: +  val - the maximum component
340: -  p - the location of val (pass NULL if you don't want this)

342:    Notes:
343:    Returns the value PETSC_MIN_REAL and p = -1 if the vector is of length 0.

345:    Returns the smallest index with the maximum value
346:    Level: intermediate

348:    Concepts: maximum^of vector
349:    Concepts: vector^maximum value

351: .seealso: VecNorm(), VecMin()
352: @*/
353: PetscErrorCode  VecMax(Vec x,PetscInt *p,PetscReal *val)
354: {

361:   PetscLogEventBegin(VEC_Max,x,0,0,0);
362:   (*x->ops->max)(x,p,val);
363:   PetscLogEventEnd(VEC_Max,x,0,0,0);
364:   return(0);
365: }

367: /*@C
368:    VecMin - Determines the minimum vector component and its location.

370:    Collective on Vec

372:    Input Parameters:
373: .  x - the vector

375:    Output Parameter:
376: +  val - the minimum component
377: -  p - the location of val (pass NULL if you don't want this location)

379:    Level: intermediate

381:    Notes:
382:    Returns the value PETSC_MAX_REAL and p = -1 if the vector is of length 0.

384:    This returns the smallest index with the minumum value

386:    Concepts: minimum^of vector
387:    Concepts: vector^minimum entry

389: .seealso: VecMax()
390: @*/
391: PetscErrorCode  VecMin(Vec x,PetscInt *p,PetscReal *val)
392: {

399:   PetscLogEventBegin(VEC_Min,x,0,0,0);
400:   (*x->ops->min)(x,p,val);
401:   PetscLogEventEnd(VEC_Min,x,0,0,0);
402:   return(0);
403: }

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

409:    Collective on Vec

411:    Input Parameters:
412: .  x, y - the vectors

414:    Output Parameter:
415: .  val - the dot product

417:    Notes for Users of Complex Numbers:
418:    For complex vectors, VecTDot() computes the indefinite form
419: $     val = (x,y) = y^T x,
420:    where y^T denotes the transpose of y.

422:    Use VecDot() for the inner product
423: $     val = (x,y) = y^H x,
424:    where y^H denotes the conjugate transpose of y.

426:    Level: intermediate

428:    Concepts: inner product^non-Hermitian
429:    Concepts: vector^inner product
430:    Concepts: non-Hermitian inner product

432: .seealso: VecDot(), VecMTDot()
433: @*/
434: PetscErrorCode  VecTDot(Vec x,Vec y,PetscScalar *val)
435: {

445:   VecCheckSameSize(x,1,y,2);

447:   PetscLogEventBegin(VEC_TDot,x,y,0,0);
448:   (*x->ops->tdot)(x,y,val);
449:   PetscLogEventEnd(VEC_TDot,x,y,0,0);
450:   return(0);
451: }

453: /*@
454:    VecScale - Scales a vector.

456:    Not collective on Vec

458:    Input Parameters:
459: +  x - the vector
460: -  alpha - the scalar

462:    Output Parameter:
463: .  x - the scaled vector

465:    Note:
466:    For a vector with n components, VecScale() computes
467: $      x[i] = alpha * x[i], for i=1,...,n.

469:    Level: intermediate

471:    Concepts: vector^scaling
472:    Concepts: scaling^vector

474: @*/
475: PetscErrorCode  VecScale(Vec x, PetscScalar alpha)
476: {
477:   PetscReal      norms[4] = {0.0,0.0,0.0, 0.0};
478:   PetscBool      flgs[4];
480:   PetscInt       i;

485:   if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
486:   PetscLogEventBegin(VEC_Scale,x,0,0,0);
487:   if (alpha != (PetscScalar)1.0) {
488:     /* get current stashed norms */
489:     for (i=0; i<4; i++) {
490:       PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],norms[i],flgs[i]);
491:     }
492:     (*x->ops->scale)(x,alpha);
493:     PetscObjectStateIncrease((PetscObject)x);
494:     /* put the scaled stashed norms back into the Vec */
495:     for (i=0; i<4; i++) {
496:       if (flgs[i]) {
497:         PetscObjectComposedDataSetReal((PetscObject)x,NormIds[i],PetscAbsScalar(alpha)*norms[i]);
498:       }
499:     }
500:   }
501:   PetscLogEventEnd(VEC_Scale,x,0,0,0);
502:   return(0);
503: }

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

508:    Logically Collective on Vec

510:    Input Parameters:
511: +  x  - the vector
512: -  alpha - the scalar

514:    Output Parameter:
515: .  x  - the vector

517:    Note:
518:    For a vector of dimension n, VecSet() computes
519: $     x[i] = alpha, for i=1,...,n,
520:    so that all vector entries then equal the identical
521:    scalar value, alpha.  Use the more general routine
522:    VecSetValues() to set different vector entries.

524:    You CANNOT call this after you have called VecSetValues() but before you call
525:    VecAssemblyBegin/End().

527:    Level: beginner

529: .seealso VecSetValues(), VecSetValuesBlocked(), VecSetRandom()

531:    Concepts: vector^setting to constant

533: @*/
534: PetscErrorCode  VecSet(Vec x,PetscScalar alpha)
535: {
536:   PetscReal      val;

542:   if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"You cannot call this after you have called VecSetValues() but\n before you have called VecAssemblyBegin/End()");
544:   VecLocked(x,1);

546:   PetscLogEventBegin(VEC_Set,x,0,0,0);
547:   (*x->ops->set)(x,alpha);
548:   PetscLogEventEnd(VEC_Set,x,0,0,0);
549:   PetscObjectStateIncrease((PetscObject)x);

551:   /*  norms can be simply set (if |alpha|*N not too large) */
552:   val  = PetscAbsScalar(alpha);
553:   if (x->map->N == 0) {
554:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],0.0l);
555:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],0.0);
556:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],0.0);
557:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],0.0);
558:   } else if (val > PETSC_MAX_REAL/x->map->N) {
559:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
560:   } else {
561:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],x->map->N * val);
562:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
563:     val  = PetscSqrtReal((PetscReal)x->map->N) * val;
564:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],val);
565:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],val);
566:   }
567:   return(0);
568: }


571: /*@
572:    VecAXPY - Computes y = alpha x + y.

574:    Logically Collective on Vec

576:    Input Parameters:
577: +  alpha - the scalar
578: -  x, y  - the vectors

580:    Output Parameter:
581: .  y - output vector

583:    Level: intermediate

585:    Notes:
586:     x and y MUST be different vectors

588:    Concepts: vector^BLAS
589:    Concepts: BLAS

591: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
592: @*/
593: PetscErrorCode  VecAXPY(Vec y,PetscScalar alpha,Vec x)
594: {

603:   VecCheckSameSize(x,1,y,3);
604:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");
606:   VecLocked(y,1);

608:   VecLockPush(x);
609:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
610:   (*y->ops->axpy)(y,alpha,x);
611:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
612:   VecLockPop(x);
613:   PetscObjectStateIncrease((PetscObject)y);
614:   return(0);
615: }

617: /*@
618:    VecAXPBY - Computes y = alpha x + beta y.

620:    Logically Collective on Vec

622:    Input Parameters:
623: +  alpha,beta - the scalars
624: -  x, y  - the vectors

626:    Output Parameter:
627: .  y - output vector

629:    Level: intermediate

631:    Notes:
632:     x and y MUST be different vectors

634:    Concepts: BLAS
635:    Concepts: vector^BLAS

637: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
638: @*/
639: PetscErrorCode  VecAXPBY(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
640: {

649:   VecCheckSameSize(x,1,y,4);
650:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");

654:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
655:   (*y->ops->axpby)(y,alpha,beta,x);
656:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
657:   PetscObjectStateIncrease((PetscObject)y);
658:   return(0);
659: }

661: /*@
662:    VecAXPBYPCZ - Computes z = alpha x + beta y + gamma z

664:    Logically Collective on Vec

666:    Input Parameters:
667: +  alpha,beta, gamma - the scalars
668: -  x, y, z  - the vectors

670:    Output Parameter:
671: .  z - output vector

673:    Level: intermediate

675:    Notes:
676:     x, y and z must be different vectors

678:    Developer Note:   alpha = 1 or gamma = 1 or gamma = 0.0 are handled as special cases

680:    Concepts: BLAS
681:    Concepts: vector^BLAS

683: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
684: @*/
685: PetscErrorCode  VecAXPBYPCZ(Vec z,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec x,Vec y)
686: {

698:   VecCheckSameSize(x,1,y,5);
699:   VecCheckSameSize(x,1,z,6);
700:   if (x == y || x == z) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");
701:   if (y == z) SETERRQ(PetscObjectComm((PetscObject)y),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");

706:   PetscLogEventBegin(VEC_AXPBYPCZ,x,y,z,0);
707:   (*y->ops->axpbypcz)(z,alpha,beta,gamma,x,y);
708:   PetscLogEventEnd(VEC_AXPBYPCZ,x,y,z,0);
709:   PetscObjectStateIncrease((PetscObject)z);
710:   return(0);
711: }

713: /*@
714:    VecAYPX - Computes y = x + alpha y.

716:    Logically Collective on Vec

718:    Input Parameters:
719: +  alpha - the scalar
720: -  x, y  - the vectors

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

725:    Level: intermediate

727:    Notes:
728:     x and y MUST be different vectors

730:    Concepts: vector^BLAS
731:    Concepts: BLAS

733: .seealso: VecAXPY(), VecWAXPY()
734: @*/
735: PetscErrorCode  VecAYPX(Vec y,PetscScalar alpha,Vec x)
736: {

745:   VecCheckSameSize(x,1,y,3);
746:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y must be different vectors");

749:   PetscLogEventBegin(VEC_AYPX,x,y,0,0);
750:    (*y->ops->aypx)(y,alpha,x);
751:   PetscLogEventEnd(VEC_AYPX,x,y,0,0);
752:   PetscObjectStateIncrease((PetscObject)y);
753:   return(0);
754: }


757: /*@
758:    VecWAXPY - Computes w = alpha x + y.

760:    Logically Collective on Vec

762:    Input Parameters:
763: +  alpha - the scalar
764: -  x, y  - the vectors

766:    Output Parameter:
767: .  w - the result

769:    Level: intermediate

771:    Notes:
772:     w cannot be either x or y, but x and y can be the same

774:    Concepts: vector^BLAS
775:    Concepts: BLAS

777: .seealso: VecAXPY(), VecAYPX(), VecAXPBY()
778: @*/
779: PetscErrorCode  VecWAXPY(Vec w,PetscScalar alpha,Vec x,Vec y)
780: {

792:   VecCheckSameSize(x,3,y,4);
793:   VecCheckSameSize(x,3,w,1);
794:   if (w == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector y, suggest VecAXPY()");
795:   if (w == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector x, suggest VecAYPX()");

798:   PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
799:    (*w->ops->waxpy)(w,alpha,x,y);
800:   PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
801:   PetscObjectStateIncrease((PetscObject)w);
802:   return(0);
803: }


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

809:    Not Collective

811:    Input Parameters:
812: +  x - vector to insert in
813: .  ni - number of elements to add
814: .  ix - indices where to add
815: .  y - array of values
816: -  iora - either INSERT_VALUES or ADD_VALUES, where
817:    ADD_VALUES adds values to any existing entries, and
818:    INSERT_VALUES replaces existing entries with new values

820:    Notes:
821:    VecSetValues() sets x[ix[i]] = y[i], for i=0,...,ni-1.

823:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
824:    options cannot be mixed without intervening calls to the assembly
825:    routines.

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

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

832:    If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
833:    negative indices may be passed in ix. These rows are
834:    simply ignored. This allows easily inserting element load matrices
835:    with homogeneous Dirchlet boundary conditions that you don't want represented
836:    in the vector.

838:    Level: beginner

840:    Concepts: vector^setting values

842: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
843:            VecSetValue(), VecSetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecGetValues()
844: @*/
845: PetscErrorCode  VecSetValues(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
846: {

851:   if (!ni) return(0);
855:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
856:   (*x->ops->setvalues)(x,ni,ix,y,iora);
857:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
858:   PetscObjectStateIncrease((PetscObject)x);
859:   return(0);
860: }

862: /*@
863:    VecGetValues - Gets values from certain locations of a vector. Currently
864:           can only get values on the same processor

866:     Not Collective

868:    Input Parameters:
869: +  x - vector to get values from
870: .  ni - number of elements to get
871: -  ix - indices where to get them from (in global 1d numbering)

873:    Output Parameter:
874: .   y - array of values

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

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

881:    VecAssemblyBegin() and VecAssemblyEnd()  MUST be called before calling this

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

885:    If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
886:    negative indices may be passed in ix. These rows are
887:    simply ignored.

889:    Level: beginner

891:    Concepts: vector^getting values

893: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues()
894: @*/
895: PetscErrorCode  VecGetValues(Vec x,PetscInt ni,const PetscInt ix[],PetscScalar y[])
896: {

901:   if (!ni) return(0);
905:   (*x->ops->getvalues)(x,ni,ix,y);
906:   return(0);
907: }

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

912:    Not Collective

914:    Input Parameters:
915: +  x - vector to insert in
916: .  ni - number of blocks to add
917: .  ix - indices where to add in block count, rather than element count
918: .  y - array of values
919: -  iora - either INSERT_VALUES or ADD_VALUES, where
920:    ADD_VALUES adds values to any existing entries, and
921:    INSERT_VALUES replaces existing entries with new values

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

927:    Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
928:    options cannot be mixed without intervening calls to the assembly
929:    routines.

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

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

936:    Negative indices may be passed in ix, these rows are
937:    simply ignored. This allows easily inserting element load matrices
938:    with homogeneous Dirchlet boundary conditions that you don't want represented
939:    in the vector.

941:    Level: intermediate

943:    Concepts: vector^setting values blocked

945: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
946:            VecSetValues()
947: @*/
948: PetscErrorCode  VecSetValuesBlocked(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
949: {

957:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
958:   (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
959:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
960:   PetscObjectStateIncrease((PetscObject)x);
961:   return(0);
962: }


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

969:    Not Collective

971:    Input Parameters:
972: +  x - vector to insert in
973: .  ni - number of elements to add
974: .  ix - indices where to add
975: .  y - array of values
976: -  iora - either INSERT_VALUES or ADD_VALUES, where
977:    ADD_VALUES adds values to any existing entries, and
978:    INSERT_VALUES replaces existing entries with new values

980:    Level: intermediate

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

985:    Calls to VecSetValues() 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 VecSetValuesLocal() have been completed.

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

994:    Concepts: vector^setting values with local numbering

996: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
997:            VecSetValuesBlockedLocal()
998: @*/
999: PetscErrorCode  VecSetValuesLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1000: {
1002:   PetscInt       lixp[128],*lix = lixp;

1006:   if (!ni) return(0);

1011:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1012:   if (!x->ops->setvalueslocal) {
1013:     if (!x->map->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
1014:     if (ni > 128) {
1015:       PetscMalloc1(ni,&lix);
1016:     }
1017:     ISLocalToGlobalMappingApply(x->map->mapping,ni,(PetscInt*)ix,lix);
1018:     (*x->ops->setvalues)(x,ni,lix,y,iora);
1019:     if (ni > 128) {
1020:       PetscFree(lix);
1021:     }
1022:   } else {
1023:     (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
1024:   }
1025:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1026:   PetscObjectStateIncrease((PetscObject)x);
1027:   return(0);
1028: }

1030: /*@
1031:    VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1032:    using a local ordering of the nodes.

1034:    Not Collective

1036:    Input Parameters:
1037: +  x - vector to insert in
1038: .  ni - number of blocks to add
1039: .  ix - indices where to add in block count, not element count
1040: .  y - array of values
1041: -  iora - either INSERT_VALUES or ADD_VALUES, where
1042:    ADD_VALUES adds values to any existing entries, and
1043:    INSERT_VALUES replaces existing entries with new values

1045:    Level: intermediate

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

1051:    Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1052:    options cannot be mixed without intervening calls to the assembly
1053:    routines.

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

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


1061:    Concepts: vector^setting values blocked with local numbering

1063: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1064:            VecSetLocalToGlobalMapping()
1065: @*/
1066: PetscErrorCode  VecSetValuesBlockedLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1067: {
1069:   PetscInt       lixp[128],*lix = lixp;

1076:   if (ni > 128) {
1077:     PetscMalloc1(ni,&lix);
1078:   }

1080:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1081:   ISLocalToGlobalMappingApplyBlock(x->map->mapping,ni,(PetscInt*)ix,lix);
1082:   (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1083:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1084:   if (ni > 128) {
1085:     PetscFree(lix);
1086:   }
1087:   PetscObjectStateIncrease((PetscObject)x);
1088:   return(0);
1089: }

1091: /*@
1092:    VecMTDot - Computes indefinite vector multiple dot products.
1093:    That is, it does NOT use the complex conjugate.

1095:    Collective on Vec

1097:    Input Parameters:
1098: +  x - one vector
1099: .  nv - number of vectors
1100: -  y - array of vectors.  Note that vectors are pointers

1102:    Output Parameter:
1103: .  val - array of the dot products

1105:    Notes for Users of Complex Numbers:
1106:    For complex vectors, VecMTDot() computes the indefinite form
1107: $      val = (x,y) = y^T x,
1108:    where y^T denotes the transpose of y.

1110:    Use VecMDot() for the inner product
1111: $      val = (x,y) = y^H x,
1112:    where y^H denotes the conjugate transpose of y.

1114:    Level: intermediate

1116:    Concepts: inner product^multiple
1117:    Concepts: vector^multiple inner products

1119: .seealso: VecMDot(), VecTDot()
1120: @*/
1121: PetscErrorCode  VecMTDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1122: {

1133:   VecCheckSameSize(x,1,*y,3);

1135:   PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1136:   (*x->ops->mtdot)(x,nv,y,val);
1137:   PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1138:   return(0);
1139: }

1141: /*@
1142:    VecMDot - Computes vector multiple dot products.

1144:    Collective on Vec

1146:    Input Parameters:
1147: +  x - one vector
1148: .  nv - number of vectors
1149: -  y - array of vectors.

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

1154:    Notes for Users of Complex Numbers:
1155:    For complex vectors, VecMDot() computes
1156: $     val = (x,y) = y^H x,
1157:    where y^H denotes the conjugate transpose of y.

1159:    Use VecMTDot() for the indefinite form
1160: $     val = (x,y) = y^T x,
1161:    where y^T denotes the transpose of y.

1163:    Level: intermediate

1165:    Concepts: inner product^multiple
1166:    Concepts: vector^multiple inner products

1168: .seealso: VecMTDot(), VecDot()
1169: @*/
1170: PetscErrorCode  VecMDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1171: {

1176:   if (!nv) return(0);
1177:   if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1184:   VecCheckSameSize(x,1,*y,3);

1186:   PetscLogEventBegin(VEC_MDot,x,*y,0,0);
1187:   (*x->ops->mdot)(x,nv,y,val);
1188:   PetscLogEventEnd(VEC_MDot,x,*y,0,0);
1189:   return(0);
1190: }

1192: /*@
1193:    VecMAXPY - Computes y = y + sum alpha[j] x[j]

1195:    Logically Collective on Vec

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

1203:    Level: intermediate

1205:    Notes:
1206:     y cannot be any of the x vectors

1208:    Concepts: BLAS

1210: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1211: @*/
1212: PetscErrorCode  VecMAXPY(Vec y,PetscInt nv,const PetscScalar alpha[],Vec x[])
1213: {
1215:   PetscInt       i;

1219:   if (!nv) return(0);
1220:   if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1227:   VecCheckSameSize(y,1,*x,4);

1230:   PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1231:   (*y->ops->maxpy)(y,nv,alpha,x);
1232:   PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1233:   PetscObjectStateIncrease((PetscObject)y);
1234:   return(0);
1235: }

1237: /*@
1238:    VecGetSubVector - Gets a vector representing part of another vector

1240:    Collective on IS (and Vec if nonlocal entries are needed)

1242:    Input Arguments:
1243: + X - vector from which to extract a subvector
1244: - is - index set representing portion of X to extract

1246:    Output Arguments:
1247: . Y - subvector corresponding to is

1249:    Level: advanced

1251:    Notes:
1252:    The subvector Y should be returned with VecRestoreSubVector().

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

1257: .seealso: MatCreateSubMatrix()
1258: @*/
1259: PetscErrorCode  VecGetSubVector(Vec X,IS is,Vec *Y)
1260: {
1261:   PetscErrorCode   ierr;
1262:   Vec              Z;

1268:   if (X->ops->getsubvector) {
1269:     (*X->ops->getsubvector)(X,is,&Z);
1270:   } else {                      /* Default implementation currently does no caching */
1271:     PetscInt  gstart,gend,start;
1272:     PetscBool contiguous,gcontiguous;
1273:     VecGetOwnershipRange(X,&gstart,&gend);
1274:     ISContiguousLocal(is,gstart,gend,&start,&contiguous);
1275:     MPIU_Allreduce(&contiguous,&gcontiguous,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)is));
1276:     if (gcontiguous) {          /* We can do a no-copy implementation */
1277:       PetscInt    n,N,bs;
1278:       PetscMPIInt size;
1279:       PetscInt    state;

1281:       ISGetLocalSize(is,&n);
1282:       VecGetBlockSize(X,&bs);
1283:       if (n%bs || bs == 1 || !n) bs = -1; /* Do not decide block size if we do not have to */
1284:       MPI_Comm_size(PetscObjectComm((PetscObject)X),&size);
1285:       VecLockGet(X,&state);
1286:       if (state) {
1287:         const PetscScalar *x;
1288:         VecGetArrayRead(X,&x);
1289:         if (size == 1) {
1290:           VecCreateSeqWithArray(PetscObjectComm((PetscObject)X),bs,n,(PetscScalar*)x+start,&Z);
1291:         } else {
1292:           ISGetSize(is,&N);
1293:           VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),bs,n,N,(PetscScalar*)x+start,&Z);
1294:         }
1295:         VecRestoreArrayRead(X,&x);
1296:         VecLockPush(Z);
1297:       } else {
1298:         PetscScalar *x;
1299:         VecGetArray(X,&x);
1300:         if (size == 1) {
1301:           VecCreateSeqWithArray(PetscObjectComm((PetscObject)X),bs,n,x+start,&Z);
1302:         } else {
1303:           ISGetSize(is,&N);
1304:           VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),bs,n,N,x+start,&Z);
1305:         }
1306:         VecRestoreArray(X,&x);
1307:       }
1308:     } else {                    /* Have to create a scatter and do a copy */
1309:       VecScatter scatter;
1310:       PetscInt   n,N;
1311:       ISGetLocalSize(is,&n);
1312:       ISGetSize(is,&N);
1313:       VecCreate(PetscObjectComm((PetscObject)is),&Z);
1314:       VecSetSizes(Z,n,N);
1315:       VecSetType(Z,((PetscObject)X)->type_name);
1316:       VecScatterCreate(X,is,Z,NULL,&scatter);
1317:       VecScatterBegin(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1318:       VecScatterEnd(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1319:       PetscObjectCompose((PetscObject)Z,"VecGetSubVector_Scatter",(PetscObject)scatter);
1320:       VecScatterDestroy(&scatter);
1321:     }
1322:   }
1323:   /* Record the state when the subvector was gotten so we know whether its values need to be put back */
1324:   if (VecGetSubVectorSavedStateId < 0) {PetscObjectComposedDataRegister(&VecGetSubVectorSavedStateId);}
1325:   PetscObjectComposedDataSetInt((PetscObject)Z,VecGetSubVectorSavedStateId,1);
1326:   *Y   = Z;
1327:   return(0);
1328: }

1330: /*@
1331:    VecRestoreSubVector - Restores a subvector extracted using VecGetSubVector()

1333:    Collective on IS (and Vec if nonlocal entries need to be written)

1335:    Input Arguments:
1336: + X - vector from which subvector was obtained
1337: . is - index set representing the subset of X
1338: - Y - subvector being restored

1340:    Level: advanced

1342: .seealso: VecGetSubVector()
1343: @*/
1344: PetscErrorCode  VecRestoreSubVector(Vec X,IS is,Vec *Y)
1345: {

1353:   if (X->ops->restoresubvector) {
1354:     (*X->ops->restoresubvector)(X,is,Y);
1355:   } else {
1356:     PETSC_UNUSED PetscObjectState dummystate = 0;
1357:     PetscBool valid;
1358:     PetscObjectComposedDataGetInt((PetscObject)*Y,VecGetSubVectorSavedStateId,dummystate,valid);
1359:     if (!valid) {
1360:       VecScatter scatter;

1362:       PetscObjectQuery((PetscObject)*Y,"VecGetSubVector_Scatter",(PetscObject*)&scatter);
1363:       if (scatter) {
1364:         VecScatterBegin(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1365:         VecScatterEnd(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1366:       }
1367:     }
1368:     VecDestroy(Y);
1369:   }
1370:   return(0);
1371: }

1373: /*@
1374:    VecGetLocalVectorRead - Maps the local portion of a vector into a
1375:    vector.  You must call VecRestoreLocalVectorRead() when the local
1376:    vector is no longer needed.

1378:    Not collective.

1380:    Input parameter:
1381: .  v - The vector for which the local vector is desired.

1383:    Output parameter:
1384: .  w - Upon exit this contains the local vector.

1386:    Level: beginner

1388:    Notes:
1389:    This function is similar to VecGetArrayRead() which maps the local
1390:    portion into a raw pointer.  VecGetLocalVectorRead() is usually
1391:    almost as efficient as VecGetArrayRead() but in certain circumstances
1392:    VecGetLocalVectorRead() can be much more efficient than
1393:    VecGetArrayRead().  This is because the construction of a contiguous
1394:    array representing the vector data required by VecGetArrayRead() can
1395:    be an expensive operation for certain vector types.  For example, for
1396:    GPU vectors VecGetArrayRead() requires that the data between device
1397:    and host is synchronized.

1399:    Unlike VecGetLocalVector(), this routine is not collective and
1400:    preserves cached information.

1402: .seealso: VecRestoreLocalVectorRead(), VecGetLocalVector(), VecGetArrayRead(), VecGetArray()
1403: @*/
1404: PetscErrorCode VecGetLocalVectorRead(Vec v,Vec w)
1405: {
1407:   PetscScalar    *a;

1412:   VecCheckSameLocalSize(v,1,w,2);
1413:   if (v->ops->getlocalvectorread) {
1414:     (*v->ops->getlocalvectorread)(v,w);
1415:   } else {
1416:     VecGetArrayRead(v,(const PetscScalar**)&a);
1417:     VecPlaceArray(w,a);
1418:   }
1419:   return(0);
1420: }

1422: /*@
1423:    VecRestoreLocalVectorRead - Unmaps the local portion of a vector
1424:    previously mapped into a vector using VecGetLocalVectorRead().

1426:    Not collective.

1428:    Input parameter:
1429: .  v - The local portion of this vector was previously mapped into w using VecGetLocalVectorRead().
1430: .  w - The vector into which the local portion of v was mapped.

1432:    Level: beginner

1434: .seealso: VecGetLocalVectorRead(), VecGetLocalVector(), VecGetArrayRead(), VecGetArray()
1435: @*/
1436: PetscErrorCode VecRestoreLocalVectorRead(Vec v,Vec w)
1437: {
1439:   PetscScalar    *a;

1444:   if (v->ops->restorelocalvectorread) {
1445:     (*v->ops->restorelocalvectorread)(v,w);
1446:   } else {
1447:     VecGetArrayRead(w,(const PetscScalar**)&a);
1448:     VecRestoreArrayRead(v,(const PetscScalar**)&a);
1449:     VecResetArray(w);
1450:   }
1451:   return(0);
1452: }

1454: /*@
1455:    VecGetLocalVector - Maps the local portion of a vector into a
1456:    vector.

1458:    Collective on v, not collective on w.

1460:    Input parameter:
1461: .  v - The vector for which the local vector is desired.

1463:    Output parameter:
1464: .  w - Upon exit this contains the local vector.

1466:    Level: beginner

1468:    Notes:
1469:    This function is similar to VecGetArray() which maps the local
1470:    portion into a raw pointer.  VecGetLocalVector() is usually about as
1471:    efficient as VecGetArray() but in certain circumstances
1472:    VecGetLocalVector() can be much more efficient than VecGetArray().
1473:    This is because the construction of a contiguous array representing
1474:    the vector data required by VecGetArray() can be an expensive
1475:    operation for certain vector types.  For example, for GPU vectors
1476:    VecGetArray() requires that the data between device and host is
1477:    synchronized.

1479: .seealso: VecRestoreLocalVector(), VecGetLocalVectorRead(), VecGetArrayRead(), VecGetArray()
1480: @*/
1481: PetscErrorCode VecGetLocalVector(Vec v,Vec w)
1482: {
1484:   PetscScalar    *a;

1489:   VecCheckSameLocalSize(v,1,w,2);
1490:   if (v->ops->getlocalvector) {
1491:     (*v->ops->getlocalvector)(v,w);
1492:   } else {
1493:     VecGetArray(v,&a);
1494:     VecPlaceArray(w,a);
1495:   }
1496:   return(0);
1497: }

1499: /*@
1500:    VecRestoreLocalVector - Unmaps the local portion of a vector
1501:    previously mapped into a vector using VecGetLocalVector().

1503:    Logically collective.

1505:    Input parameter:
1506: .  v - The local portion of this vector was previously mapped into w using VecGetLocalVector().
1507: .  w - The vector into which the local portion of v was mapped.

1509:    Level: beginner

1511: .seealso: VecGetLocalVector(), VecGetLocalVectorRead(), VecRestoreLocalVectorRead(), LocalVectorRead(), VecGetArrayRead(), VecGetArray()
1512: @*/
1513: PetscErrorCode VecRestoreLocalVector(Vec v,Vec w)
1514: {
1516:   PetscScalar    *a;

1521:   if (v->ops->restorelocalvector) {
1522:     (*v->ops->restorelocalvector)(v,w);
1523:   } else {
1524:     VecGetArray(w,&a);
1525:     VecRestoreArray(v,&a);
1526:     VecResetArray(w);
1527:   }
1528:   return(0);
1529: }

1531: /*@C
1532:    VecGetArray - Returns a pointer to a contiguous array that contains this
1533:    processor's portion of the vector data. For the standard PETSc
1534:    vectors, VecGetArray() returns a pointer to the local data array and
1535:    does not use any copies. If the underlying vector data is not stored
1536:    in a contiguous array this routine will copy the data to a contiguous
1537:    array and return a pointer to that. You MUST call VecRestoreArray()
1538:    when you no longer need access to the array.

1540:    Logically Collective on Vec

1542:    Input Parameter:
1543: .  x - the vector

1545:    Output Parameter:
1546: .  a - location to put pointer to the array

1548:    Fortran Note:
1549:    This routine is used differently from Fortran 77
1550: $    Vec         x
1551: $    PetscScalar x_array(1)
1552: $    PetscOffset i_x
1553: $    PetscErrorCode ierr
1554: $       call VecGetArray(x,x_array,i_x,ierr)
1555: $
1556: $   Access first local entry in vector with
1557: $      value = x_array(i_x + 1)
1558: $
1559: $      ...... other code
1560: $       call VecRestoreArray(x,x_array,i_x,ierr)
1561:    For Fortran 90 see VecGetArrayF90()

1563:    See the Fortran chapter of the users manual and
1564:    petsc/src/snes/examples/tutorials/ex5f.F for details.

1566:    Level: beginner

1568:    Concepts: vector^accessing local values

1570: .seealso: VecRestoreArray(), VecGetArrayRead(), VecGetArrays(), VecGetArrayF90(), VecGetArrayReadF90(), VecPlaceArray(), VecGetArray2d(),
1571:           VecGetArrayPair(), VecRestoreArrayPair()
1572: @*/
1573: PetscErrorCode VecGetArray(Vec x,PetscScalar **a)
1574: {
1576: #if defined(PETSC_HAVE_VIENNACL)
1577:   PetscBool      is_viennacltype = PETSC_FALSE;
1578: #endif

1582:   VecLocked(x,1);
1583:   if (x->petscnative) {
1584: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1585:     if (x->valid_GPU_array == PETSC_OFFLOAD_GPU) {
1586: #if defined(PETSC_HAVE_VIENNACL)
1587:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1588:       if (is_viennacltype) {
1589:         VecViennaCLCopyFromGPU(x);
1590:       } else
1591: #endif
1592:       {
1593: #if defined(PETSC_HAVE_VECCUDA)
1594:         VecCUDACopyFromGPU(x);
1595: #endif
1596:       }
1597:     } else if (x->valid_GPU_array == PETSC_OFFLOAD_UNALLOCATED) {
1598: #if defined(PETSC_HAVE_VIENNACL)
1599:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1600:       if (is_viennacltype) {
1601:         VecViennaCLAllocateCheckHost(x);
1602:       } else
1603: #endif
1604:       {
1605: #if defined(PETSC_HAVE_VECCUDA)
1606:         VecCUDAAllocateCheckHost(x);
1607: #endif
1608:       }
1609:     }
1610: #endif
1611:     *a = *((PetscScalar**)x->data);
1612:   } else {
1613:     (*x->ops->getarray)(x,a);
1614:   }
1615:   return(0);
1616: }

1618: /*@C
1619:    VecGetArrayRead - Get read-only pointer to contiguous array containing this processor's portion of the vector data.

1621:    Not Collective

1623:    Input Parameters:
1624: .  x - the vector

1626:    Output Parameter:
1627: .  a - the array

1629:    Level: beginner

1631:    Notes:
1632:    The array must be returned using a matching call to VecRestoreArrayRead().

1634:    Unlike VecGetArray(), this routine is not collective and preserves cached information like vector norms.

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

1640: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrayPair(), VecRestoreArrayPair()
1641: @*/
1642: PetscErrorCode VecGetArrayRead(Vec x,const PetscScalar **a)
1643: {
1645: #if defined(PETSC_HAVE_VIENNACL)
1646:   PetscBool      is_viennacltype = PETSC_FALSE;
1647: #endif

1651:   if (x->petscnative) {
1652: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1653:     if (x->valid_GPU_array == PETSC_OFFLOAD_GPU) {
1654: #if defined(PETSC_HAVE_VIENNACL)
1655:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1656:       if (is_viennacltype) {
1657:         VecViennaCLCopyFromGPU(x);
1658:       } else
1659: #endif
1660:       {
1661: #if defined(PETSC_HAVE_VECCUDA)
1662:         VecCUDACopyFromGPU(x);
1663: #endif
1664:       }
1665:     }
1666: #endif
1667:     *a = *((PetscScalar **)x->data);
1668:   } else if (x->ops->getarrayread) {
1669:     (*x->ops->getarrayread)(x,a);
1670:   } else {
1671:     (*x->ops->getarray)(x,(PetscScalar**)a);
1672:   }
1673:   return(0);
1674: }

1676: /*@C
1677:    VecGetArrays - Returns a pointer to the arrays in a set of vectors
1678:    that were created by a call to VecDuplicateVecs().  You MUST call
1679:    VecRestoreArrays() when you no longer need access to the array.

1681:    Logically Collective on Vec

1683:    Input Parameter:
1684: +  x - the vectors
1685: -  n - the number of vectors

1687:    Output Parameter:
1688: .  a - location to put pointer to the array

1690:    Fortran Note:
1691:    This routine is not supported in Fortran.

1693:    Level: intermediate

1695: .seealso: VecGetArray(), VecRestoreArrays()
1696: @*/
1697: PetscErrorCode  VecGetArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1698: {
1700:   PetscInt       i;
1701:   PetscScalar    **q;

1707:   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %D",n);
1708:   PetscMalloc1(n,&q);
1709:   for (i=0; i<n; ++i) {
1710:     VecGetArray(x[i],&q[i]);
1711:   }
1712:   *a = q;
1713:   return(0);
1714: }

1716: /*@C
1717:    VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1718:    has been called.

1720:    Logically Collective on Vec

1722:    Input Parameters:
1723: +  x - the vector
1724: .  n - the number of vectors
1725: -  a - location of pointer to arrays obtained from VecGetArrays()

1727:    Notes:
1728:    For regular PETSc vectors this routine does not involve any copies. For
1729:    any special vectors that do not store local vector data in a contiguous
1730:    array, this routine will copy the data back into the underlying
1731:    vector data structure from the arrays obtained with VecGetArrays().

1733:    Fortran Note:
1734:    This routine is not supported in Fortran.

1736:    Level: intermediate

1738: .seealso: VecGetArrays(), VecRestoreArray()
1739: @*/
1740: PetscErrorCode  VecRestoreArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1741: {
1743:   PetscInt       i;
1744:   PetscScalar    **q = *a;


1751:   for (i=0; i<n; ++i) {
1752:     VecRestoreArray(x[i],&q[i]);
1753:   }
1754:   PetscFree(q);
1755:   return(0);
1756: }

1758: /*@C
1759:    VecRestoreArray - Restores a vector after VecGetArray() has been called.

1761:    Logically Collective on Vec

1763:    Input Parameters:
1764: +  x - the vector
1765: -  a - location of pointer to array obtained from VecGetArray()

1767:    Level: beginner

1769:    Notes:
1770:    For regular PETSc vectors this routine does not involve any copies. For
1771:    any special vectors that do not store local vector data in a contiguous
1772:    array, this routine will copy the data back into the underlying
1773:    vector data structure from the array obtained with VecGetArray().

1775:    This routine actually zeros out the a pointer. This is to prevent accidental
1776:    us of the array after it has been restored. If you pass null for a it will
1777:    not zero the array pointer a.

1779:    Fortran Note:
1780:    This routine is used differently from Fortran 77
1781: $    Vec         x
1782: $    PetscScalar x_array(1)
1783: $    PetscOffset i_x
1784: $    PetscErrorCode ierr
1785: $       call VecGetArray(x,x_array,i_x,ierr)
1786: $
1787: $   Access first local entry in vector with
1788: $      value = x_array(i_x + 1)
1789: $
1790: $      ...... other code
1791: $       call VecRestoreArray(x,x_array,i_x,ierr)

1793:    See the Fortran chapter of the users manual and
1794:    petsc/src/snes/examples/tutorials/ex5f.F for details.
1795:    For Fortran 90 see VecRestoreArrayF90()

1797: .seealso: VecGetArray(), VecRestoreArrayRead(), VecRestoreArrays(), VecRestoreArrayF90(), VecRestoreArrayReadF90(), VecPlaceArray(), VecRestoreArray2d(),
1798:           VecGetArrayPair(), VecRestoreArrayPair()
1799: @*/
1800: PetscErrorCode VecRestoreArray(Vec x,PetscScalar **a)
1801: {

1806:   if (x->petscnative) {
1807: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
1808:     x->valid_GPU_array = PETSC_OFFLOAD_CPU;
1809: #endif
1810:   } else {
1811:     (*x->ops->restorearray)(x,a);
1812:   }
1813:   if (a) *a = NULL;
1814:   PetscObjectStateIncrease((PetscObject)x);
1815:   return(0);
1816: }

1818: /*@C
1819:    VecRestoreArrayRead - Restore array obtained with VecGetArrayRead()

1821:    Not Collective

1823:    Input Parameters:
1824: +  vec - the vector
1825: -  array - the array

1827:    Level: beginner

1829: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrayPair(), VecRestoreArrayPair()
1830: @*/
1831: PetscErrorCode VecRestoreArrayRead(Vec x,const PetscScalar **a)
1832: {

1837:   if (x->petscnative) {
1838:     /* nothing */
1839:   } else if (x->ops->restorearrayread) {
1840:     (*x->ops->restorearrayread)(x,a);
1841:   } else {
1842:     (*x->ops->restorearray)(x,(PetscScalar**)a);
1843:   }
1844:   if (a) *a = NULL;
1845:   return(0);
1846: }

1848: /*@
1849:    VecPlaceArray - Allows one to replace the array in a vector with an
1850:    array provided by the user. This is useful to avoid copying an array
1851:    into a vector.

1853:    Not Collective

1855:    Input Parameters:
1856: +  vec - the vector
1857: -  array - the array

1859:    Notes:
1860:    You can return to the original array with a call to VecResetArray()

1862:    Level: developer

1864: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()

1866: @*/
1867: PetscErrorCode  VecPlaceArray(Vec vec,const PetscScalar array[])
1868: {

1875:   if (vec->ops->placearray) {
1876:     (*vec->ops->placearray)(vec,array);
1877:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot place array in this type of vector");
1878:   PetscObjectStateIncrease((PetscObject)vec);
1879:   return(0);
1880: }

1882: /*@C
1883:    VecReplaceArray - Allows one to replace the array in a vector with an
1884:    array provided by the user. This is useful to avoid copying an array
1885:    into a vector.

1887:    Not Collective

1889:    Input Parameters:
1890: +  vec - the vector
1891: -  array - the array

1893:    Notes:
1894:    This permanently replaces the array and frees the memory associated
1895:    with the old array.

1897:    The memory passed in MUST be obtained with PetscMalloc() and CANNOT be
1898:    freed by the user. It will be freed when the vector is destroy.

1900:    Not supported from Fortran

1902:    Level: developer

1904: .seealso: VecGetArray(), VecRestoreArray(), VecPlaceArray(), VecResetArray()

1906: @*/
1907: PetscErrorCode  VecReplaceArray(Vec vec,const PetscScalar array[])
1908: {

1914:   if (vec->ops->replacearray) {
1915:     (*vec->ops->replacearray)(vec,array);
1916:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot replace array in this type of vector");
1917:   PetscObjectStateIncrease((PetscObject)vec);
1918:   return(0);
1919: }

1921: /*MC
1922:     VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
1923:     and makes them accessible via a Fortran90 pointer.

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

1928:     Collective on Vec

1930:     Input Parameters:
1931: +   x - a vector to mimic
1932: -   n - the number of vectors to obtain

1934:     Output Parameters:
1935: +   y - Fortran90 pointer to the array of vectors
1936: -   ierr - error code

1938:     Example of Usage:
1939: .vb
1940: #include <petsc/finclude/petscvec.h>
1941:     use petscvec

1943:     Vec x
1944:     Vec, pointer :: y(:)
1945:     ....
1946:     call VecDuplicateVecsF90(x,2,y,ierr)
1947:     call VecSet(y(2),alpha,ierr)
1948:     call VecSet(y(2),alpha,ierr)
1949:     ....
1950:     call VecDestroyVecsF90(2,y,ierr)
1951: .ve

1953:     Notes:
1954:     Not yet supported for all F90 compilers

1956:     Use VecDestroyVecsF90() to free the space.

1958:     Level: beginner

1960: .seealso:  VecDestroyVecsF90(), VecDuplicateVecs()

1962: M*/

1964: /*MC
1965:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
1966:     VecGetArrayF90().

1968:     Synopsis:
1969:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

1971:     Logically Collective on Vec

1973:     Input Parameters:
1974: +   x - vector
1975: -   xx_v - the Fortran90 pointer to the array

1977:     Output Parameter:
1978: .   ierr - error code

1980:     Example of Usage:
1981: .vb
1982: #include <petsc/finclude/petscvec.h>
1983:     use petscvec

1985:     PetscScalar, pointer :: xx_v(:)
1986:     ....
1987:     call VecGetArrayF90(x,xx_v,ierr)
1988:     xx_v(3) = a
1989:     call VecRestoreArrayF90(x,xx_v,ierr)
1990: .ve

1992:     Level: beginner

1994: .seealso:  VecGetArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran, VecRestoreArrayReadF90()

1996: M*/

1998: /*MC
1999:     VecDestroyVecsF90 - Frees a block of vectors obtained with VecDuplicateVecsF90().

2001:     Synopsis:
2002:     VecDestroyVecsF90(PetscInt n,{Vec, pointer :: x(:)},PetscErrorCode ierr)

2004:     Collective on Vec

2006:     Input Parameters:
2007: +   n - the number of vectors previously obtained
2008: -   x - pointer to array of vector pointers

2010:     Output Parameter:
2011: .   ierr - error code

2013:     Notes:
2014:     Not yet supported for all F90 compilers

2016:     Level: beginner

2018: .seealso:  VecDestroyVecs(), VecDuplicateVecsF90()

2020: M*/

2022: /*MC
2023:     VecGetArrayF90 - Accesses a vector array from Fortran90. For default PETSc
2024:     vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
2025:     this routine is implementation dependent. You MUST call VecRestoreArrayF90()
2026:     when you no longer need access to the array.

2028:     Synopsis:
2029:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2031:     Logically Collective on Vec

2033:     Input Parameter:
2034: .   x - vector

2036:     Output Parameters:
2037: +   xx_v - the Fortran90 pointer to the array
2038: -   ierr - error code

2040:     Example of Usage:
2041: .vb
2042: #include <petsc/finclude/petscvec.h>
2043:     use petscvec

2045:     PetscScalar, pointer :: xx_v(:)
2046:     ....
2047:     call VecGetArrayF90(x,xx_v,ierr)
2048:     xx_v(3) = a
2049:     call VecRestoreArrayF90(x,xx_v,ierr)
2050: .ve

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

2054:     Level: beginner

2056: .seealso:  VecRestoreArrayF90(), VecGetArray(), VecRestoreArray(), VecGetArrayReadF90(), UsingFortran

2058: M*/

2060:  /*MC
2061:     VecGetArrayReadF90 - Accesses a read only array from Fortran90. For default PETSc
2062:     vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
2063:     this routine is implementation dependent. You MUST call VecRestoreArrayReadF90()
2064:     when you no longer need access to the array.

2066:     Synopsis:
2067:     VecGetArrayReadF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2069:     Logically Collective on Vec

2071:     Input Parameter:
2072: .   x - vector

2074:     Output Parameters:
2075: +   xx_v - the Fortran90 pointer to the array
2076: -   ierr - error code

2078:     Example of Usage:
2079: .vb
2080: #include <petsc/finclude/petscvec.h>
2081:     use petscvec

2083:     PetscScalar, pointer :: xx_v(:)
2084:     ....
2085:     call VecGetArrayReadF90(x,xx_v,ierr)
2086:     a = xx_v(3)
2087:     call VecRestoreArrayReadF90(x,xx_v,ierr)
2088: .ve

2090:     If you intend to write entries into the array you must use VecGetArrayF90().

2092:     Level: beginner

2094: .seealso:  VecRestoreArrayReadF90(), VecGetArray(), VecRestoreArray(), VecGetArrayRead(), VecRestoreArrayRead(), VecGetArrayF90(), UsingFortran

2096: M*/

2098: /*MC
2099:     VecRestoreArrayReadF90 - Restores a readonly vector to a usable state after a call to
2100:     VecGetArrayReadF90().

2102:     Synopsis:
2103:     VecRestoreArrayReadF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2105:     Logically Collective on Vec

2107:     Input Parameters:
2108: +   x - vector
2109: -   xx_v - the Fortran90 pointer to the array

2111:     Output Parameter:
2112: .   ierr - error code

2114:     Example of Usage:
2115: .vb
2116: #include <petsc/finclude/petscvec.h>
2117:     use petscvec

2119:     PetscScalar, pointer :: xx_v(:)
2120:     ....
2121:     call VecGetArrayReadF90(x,xx_v,ierr)
2122:     a = xx_v(3)
2123:     call VecRestoreArrayReadF90(x,xx_v,ierr)
2124: .ve

2126:     Level: beginner

2128: .seealso:  VecGetArrayReadF90(), VecGetArray(), VecRestoreArray(), VecGetArrayRead(), VecRestoreArrayRead(),UsingFortran, VecRestoreArrayF90()

2130: M*/

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

2137:    Logically Collective

2139:    Input Parameter:
2140: +  x - the vector
2141: .  m - first dimension of two dimensional array
2142: .  n - second dimension of two dimensional array
2143: .  mstart - first index you will use in first coordinate direction (often 0)
2144: -  nstart - first index in the second coordinate direction (often 0)

2146:    Output Parameter:
2147: .  a - location to put pointer to the array

2149:    Level: developer

2151:   Notes:
2152:    For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
2153:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2154:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2155:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().

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

2159:    Concepts: vector^accessing local values as 2d array

2161: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2162:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2163:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2164: @*/
2165: PetscErrorCode  VecGetArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2166: {
2168:   PetscInt       i,N;
2169:   PetscScalar    *aa;

2175:   VecGetLocalSize(x,&N);
2176:   if (m*n != N) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 2d array dimensions %D by %D",N,m,n);
2177:   VecGetArray(x,&aa);

2179:   PetscMalloc1(m,a);
2180:   for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
2181:   *a -= mstart;
2182:   return(0);
2183: }

2185: /*@C
2186:    VecRestoreArray2d - Restores a vector after VecGetArray2d() has been called.

2188:    Logically Collective

2190:    Input Parameters:
2191: +  x - the vector
2192: .  m - first dimension of two dimensional array
2193: .  n - second dimension of the two dimensional array
2194: .  mstart - first index you will use in first coordinate direction (often 0)
2195: .  nstart - first index in the second coordinate direction (often 0)
2196: -  a - location of pointer to array obtained from VecGetArray2d()

2198:    Level: developer

2200:    Notes:
2201:    For regular PETSc vectors this routine does not involve any copies. For
2202:    any special vectors that do not store local vector data in a contiguous
2203:    array, this routine will copy the data back into the underlying
2204:    vector data structure from the array obtained with VecGetArray().

2206:    This routine actually zeros out the a pointer.

2208: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2209:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2210:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2211: @*/
2212: PetscErrorCode  VecRestoreArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2213: {
2215:   void           *dummy;

2221:   dummy = (void*)(*a + mstart);
2222:   PetscFree(dummy);
2223:   VecRestoreArray(x,NULL);
2224:   return(0);
2225: }

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

2232:    Logically Collective

2234:    Input Parameter:
2235: +  x - the vector
2236: .  m - first dimension of two dimensional array
2237: -  mstart - first index you will use in first coordinate direction (often 0)

2239:    Output Parameter:
2240: .  a - location to put pointer to the array

2242:    Level: developer

2244:   Notes:
2245:    For a vector obtained from DMCreateLocalVector() mstart are likely
2246:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2247:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().

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

2251: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2252:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2253:           VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2254: @*/
2255: PetscErrorCode  VecGetArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2256: {
2258:   PetscInt       N;

2264:   VecGetLocalSize(x,&N);
2265:   if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2266:   VecGetArray(x,a);
2267:   *a  -= mstart;
2268:   return(0);
2269: }

2271: /*@C
2272:    VecRestoreArray1d - Restores a vector after VecGetArray1d() has been called.

2274:    Logically Collective

2276:    Input Parameters:
2277: +  x - the vector
2278: .  m - first dimension of two dimensional array
2279: .  mstart - first index you will use in first coordinate direction (often 0)
2280: -  a - location of pointer to array obtained from VecGetArray21()

2282:    Level: developer

2284:    Notes:
2285:    For regular PETSc vectors this routine does not involve any copies. For
2286:    any special vectors that do not store local vector data in a contiguous
2287:    array, this routine will copy the data back into the underlying
2288:    vector data structure from the array obtained with VecGetArray1d().

2290:    This routine actually zeros out the a pointer.

2292:    Concepts: vector^accessing local values as 1d array

2294: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2295:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2296:           VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2297: @*/
2298: PetscErrorCode  VecRestoreArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2299: {

2305:   VecRestoreArray(x,NULL);
2306:   return(0);
2307: }


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

2315:    Logically Collective

2317:    Input Parameter:
2318: +  x - the vector
2319: .  m - first dimension of three dimensional array
2320: .  n - second dimension of three dimensional array
2321: .  p - third dimension of three dimensional array
2322: .  mstart - first index you will use in first coordinate direction (often 0)
2323: .  nstart - first index in the second coordinate direction (often 0)
2324: -  pstart - first index in the third coordinate direction (often 0)

2326:    Output Parameter:
2327: .  a - location to put pointer to the array

2329:    Level: developer

2331:   Notes:
2332:    For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2333:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2334:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2335:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().

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

2339:    Concepts: vector^accessing local values as 3d array

2341: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2342:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2343:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2344: @*/
2345: PetscErrorCode  VecGetArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2346: {
2348:   PetscInt       i,N,j;
2349:   PetscScalar    *aa,**b;

2355:   VecGetLocalSize(x,&N);
2356:   if (m*n*p != N) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 3d array dimensions %D by %D by %D",N,m,n,p);
2357:   VecGetArray(x,&aa);

2359:   PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2360:   b    = (PetscScalar**)((*a) + m);
2361:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2362:   for (i=0; i<m; i++)
2363:     for (j=0; j<n; j++)
2364:       b[i*n+j] = aa + i*n*p + j*p - pstart;

2366:   *a -= mstart;
2367:   return(0);
2368: }

2370: /*@C
2371:    VecRestoreArray3d - Restores a vector after VecGetArray3d() has been called.

2373:    Logically Collective

2375:    Input Parameters:
2376: +  x - the vector
2377: .  m - first dimension of three dimensional array
2378: .  n - second dimension of the three dimensional array
2379: .  p - third dimension of the three dimensional array
2380: .  mstart - first index you will use in first coordinate direction (often 0)
2381: .  nstart - first index in the second coordinate direction (often 0)
2382: .  pstart - first index in the third coordinate direction (often 0)
2383: -  a - location of pointer to array obtained from VecGetArray3d()

2385:    Level: developer

2387:    Notes:
2388:    For regular PETSc vectors this routine does not involve any copies. For
2389:    any special vectors that do not store local vector data in a contiguous
2390:    array, this routine will copy the data back into the underlying
2391:    vector data structure from the array obtained with VecGetArray().

2393:    This routine actually zeros out the a pointer.

2395: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2396:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2397:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2398: @*/
2399: PetscErrorCode  VecRestoreArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2400: {
2402:   void           *dummy;

2408:   dummy = (void*)(*a + mstart);
2409:   PetscFree(dummy);
2410:   VecRestoreArray(x,NULL);
2411:   return(0);
2412: }

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

2419:    Logically Collective

2421:    Input Parameter:
2422: +  x - the vector
2423: .  m - first dimension of four dimensional array
2424: .  n - second dimension of four dimensional array
2425: .  p - third dimension of four dimensional array
2426: .  q - fourth dimension of four dimensional array
2427: .  mstart - first index you will use in first coordinate direction (often 0)
2428: .  nstart - first index in the second coordinate direction (often 0)
2429: .  pstart - first index in the third coordinate direction (often 0)
2430: -  qstart - first index in the fourth coordinate direction (often 0)

2432:    Output Parameter:
2433: .  a - location to put pointer to the array

2435:    Level: beginner

2437:   Notes:
2438:    For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2439:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2440:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2441:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().

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

2445:    Concepts: vector^accessing local values as 3d array

2447: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2448:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2449:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2450: @*/
2451: PetscErrorCode  VecGetArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2452: {
2454:   PetscInt       i,N,j,k;
2455:   PetscScalar    *aa,***b,**c;

2461:   VecGetLocalSize(x,&N);
2462:   if (m*n*p*q != N) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 4d array dimensions %D by %D by %D by %D",N,m,n,p,q);
2463:   VecGetArray(x,&aa);

2465:   PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2466:   b    = (PetscScalar***)((*a) + m);
2467:   c    = (PetscScalar**)(b + m*n);
2468:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2469:   for (i=0; i<m; i++)
2470:     for (j=0; j<n; j++)
2471:       b[i*n+j] = c + i*n*p + j*p - pstart;
2472:   for (i=0; i<m; i++)
2473:     for (j=0; j<n; j++)
2474:       for (k=0; k<p; k++)
2475:         c[i*n*p+j*p+k] = aa + i*n*p*q + j*p*q + k*q - qstart;
2476:   *a -= mstart;
2477:   return(0);
2478: }

2480: /*@C
2481:    VecRestoreArray4d - Restores a vector after VecGetArray3d() has been called.

2483:    Logically Collective

2485:    Input Parameters:
2486: +  x - the vector
2487: .  m - first dimension of four dimensional array
2488: .  n - second dimension of the four dimensional array
2489: .  p - third dimension of the four dimensional array
2490: .  q - fourth dimension of the four dimensional array
2491: .  mstart - first index you will use in first coordinate direction (often 0)
2492: .  nstart - first index in the second coordinate direction (often 0)
2493: .  pstart - first index in the third coordinate direction (often 0)
2494: .  qstart - first index in the fourth coordinate direction (often 0)
2495: -  a - location of pointer to array obtained from VecGetArray4d()

2497:    Level: beginner

2499:    Notes:
2500:    For regular PETSc vectors this routine does not involve any copies. For
2501:    any special vectors that do not store local vector data in a contiguous
2502:    array, this routine will copy the data back into the underlying
2503:    vector data structure from the array obtained with VecGetArray().

2505:    This routine actually zeros out the a pointer.

2507: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2508:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2509:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2510: @*/
2511: PetscErrorCode  VecRestoreArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2512: {
2514:   void           *dummy;

2520:   dummy = (void*)(*a + mstart);
2521:   PetscFree(dummy);
2522:   VecRestoreArray(x,NULL);
2523:   return(0);
2524: }

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

2531:    Logically Collective

2533:    Input Parameter:
2534: +  x - the vector
2535: .  m - first dimension of two dimensional array
2536: .  n - second dimension of two dimensional array
2537: .  mstart - first index you will use in first coordinate direction (often 0)
2538: -  nstart - first index in the second coordinate direction (often 0)

2540:    Output Parameter:
2541: .  a - location to put pointer to the array

2543:    Level: developer

2545:   Notes:
2546:    For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
2547:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2548:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2549:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().

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

2553:    Concepts: vector^accessing local values as 2d array

2555: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2556:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2557:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2558: @*/
2559: PetscErrorCode  VecGetArray2dRead(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2560: {
2561:   PetscErrorCode    ierr;
2562:   PetscInt          i,N;
2563:   const PetscScalar *aa;

2569:   VecGetLocalSize(x,&N);
2570:   if (m*n != N) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 2d array dimensions %D by %D",N,m,n);
2571:   VecGetArrayRead(x,&aa);

2573:   PetscMalloc1(m,a);
2574:   for (i=0; i<m; i++) (*a)[i] = (PetscScalar*) aa + i*n - nstart;
2575:   *a -= mstart;
2576:   return(0);
2577: }

2579: /*@C
2580:    VecRestoreArray2dRead - Restores a vector after VecGetArray2dRead() has been called.

2582:    Logically Collective

2584:    Input Parameters:
2585: +  x - the vector
2586: .  m - first dimension of two dimensional array
2587: .  n - second dimension of the two dimensional array
2588: .  mstart - first index you will use in first coordinate direction (often 0)
2589: .  nstart - first index in the second coordinate direction (often 0)
2590: -  a - location of pointer to array obtained from VecGetArray2d()

2592:    Level: developer

2594:    Notes:
2595:    For regular PETSc vectors this routine does not involve any copies. For
2596:    any special vectors that do not store local vector data in a contiguous
2597:    array, this routine will copy the data back into the underlying
2598:    vector data structure from the array obtained with VecGetArray().

2600:    This routine actually zeros out the a pointer.

2602: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2603:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2604:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2605: @*/
2606: PetscErrorCode  VecRestoreArray2dRead(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2607: {
2609:   void           *dummy;

2615:   dummy = (void*)(*a + mstart);
2616:   PetscFree(dummy);
2617:   VecRestoreArrayRead(x,NULL);
2618:   return(0);
2619: }

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

2626:    Logically Collective

2628:    Input Parameter:
2629: +  x - the vector
2630: .  m - first dimension of two dimensional array
2631: -  mstart - first index you will use in first coordinate direction (often 0)

2633:    Output Parameter:
2634: .  a - location to put pointer to the array

2636:    Level: developer

2638:   Notes:
2639:    For a vector obtained from DMCreateLocalVector() mstart are likely
2640:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2641:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().

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

2645: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2646:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2647:           VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2648: @*/
2649: PetscErrorCode  VecGetArray1dRead(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2650: {
2652:   PetscInt       N;

2658:   VecGetLocalSize(x,&N);
2659:   if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2660:   VecGetArrayRead(x,(const PetscScalar**)a);
2661:   *a  -= mstart;
2662:   return(0);
2663: }

2665: /*@C
2666:    VecRestoreArray1dRead - Restores a vector after VecGetArray1dRead() has been called.

2668:    Logically Collective

2670:    Input Parameters:
2671: +  x - the vector
2672: .  m - first dimension of two dimensional array
2673: .  mstart - first index you will use in first coordinate direction (often 0)
2674: -  a - location of pointer to array obtained from VecGetArray21()

2676:    Level: developer

2678:    Notes:
2679:    For regular PETSc vectors this routine does not involve any copies. For
2680:    any special vectors that do not store local vector data in a contiguous
2681:    array, this routine will copy the data back into the underlying
2682:    vector data structure from the array obtained with VecGetArray1dRead().

2684:    This routine actually zeros out the a pointer.

2686:    Concepts: vector^accessing local values as 1d array

2688: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2689:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2690:           VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2691: @*/
2692: PetscErrorCode  VecRestoreArray1dRead(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2693: {

2699:   VecRestoreArrayRead(x,NULL);
2700:   return(0);
2701: }


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

2709:    Logically Collective

2711:    Input Parameter:
2712: +  x - the vector
2713: .  m - first dimension of three dimensional array
2714: .  n - second dimension of three dimensional array
2715: .  p - third dimension of three dimensional array
2716: .  mstart - first index you will use in first coordinate direction (often 0)
2717: .  nstart - first index in the second coordinate direction (often 0)
2718: -  pstart - first index in the third coordinate direction (often 0)

2720:    Output Parameter:
2721: .  a - location to put pointer to the array

2723:    Level: developer

2725:   Notes:
2726:    For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2727:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2728:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2729:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3dRead().

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

2733:    Concepts: vector^accessing local values as 3d array

2735: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2736:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2737:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2738: @*/
2739: PetscErrorCode  VecGetArray3dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2740: {
2741:   PetscErrorCode    ierr;
2742:   PetscInt          i,N,j;
2743:   const PetscScalar *aa;
2744:   PetscScalar       **b;

2750:   VecGetLocalSize(x,&N);
2751:   if (m*n*p != N) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 3d array dimensions %D by %D by %D",N,m,n,p);
2752:   VecGetArrayRead(x,&aa);

2754:   PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2755:   b    = (PetscScalar**)((*a) + m);
2756:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2757:   for (i=0; i<m; i++)
2758:     for (j=0; j<n; j++)
2759:       b[i*n+j] = (PetscScalar *)aa + i*n*p + j*p - pstart;

2761:   *a -= mstart;
2762:   return(0);
2763: }

2765: /*@C
2766:    VecRestoreArray3dRead - Restores a vector after VecGetArray3dRead() has been called.

2768:    Logically Collective

2770:    Input Parameters:
2771: +  x - the vector
2772: .  m - first dimension of three dimensional array
2773: .  n - second dimension of the three dimensional array
2774: .  p - third dimension of the three dimensional array
2775: .  mstart - first index you will use in first coordinate direction (often 0)
2776: .  nstart - first index in the second coordinate direction (often 0)
2777: .  pstart - first index in the third coordinate direction (often 0)
2778: -  a - location of pointer to array obtained from VecGetArray3dRead()

2780:    Level: developer

2782:    Notes:
2783:    For regular PETSc vectors this routine does not involve any copies. For
2784:    any special vectors that do not store local vector data in a contiguous
2785:    array, this routine will copy the data back into the underlying
2786:    vector data structure from the array obtained with VecGetArray().

2788:    This routine actually zeros out the a pointer.

2790: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2791:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2792:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2793: @*/
2794: PetscErrorCode  VecRestoreArray3dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2795: {
2797:   void           *dummy;

2803:   dummy = (void*)(*a + mstart);
2804:   PetscFree(dummy);
2805:   VecRestoreArrayRead(x,NULL);
2806:   return(0);
2807: }

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

2814:    Logically Collective

2816:    Input Parameter:
2817: +  x - the vector
2818: .  m - first dimension of four dimensional array
2819: .  n - second dimension of four dimensional array
2820: .  p - third dimension of four dimensional array
2821: .  q - fourth dimension of four dimensional array
2822: .  mstart - first index you will use in first coordinate direction (often 0)
2823: .  nstart - first index in the second coordinate direction (often 0)
2824: .  pstart - first index in the third coordinate direction (often 0)
2825: -  qstart - first index in the fourth coordinate direction (often 0)

2827:    Output Parameter:
2828: .  a - location to put pointer to the array

2830:    Level: beginner

2832:   Notes:
2833:    For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2834:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2835:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2836:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().

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

2840:    Concepts: vector^accessing local values as 3d array

2842: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2843:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2844:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2845: @*/
2846: PetscErrorCode  VecGetArray4dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2847: {
2848:   PetscErrorCode    ierr;
2849:   PetscInt          i,N,j,k;
2850:   const PetscScalar *aa;
2851:   PetscScalar       ***b,**c;

2857:   VecGetLocalSize(x,&N);
2858:   if (m*n*p*q != N) SETERRQ5(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 4d array dimensions %D by %D by %D by %D",N,m,n,p,q);
2859:   VecGetArrayRead(x,&aa);

2861:   PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2862:   b    = (PetscScalar***)((*a) + m);
2863:   c    = (PetscScalar**)(b + m*n);
2864:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2865:   for (i=0; i<m; i++)
2866:     for (j=0; j<n; j++)
2867:       b[i*n+j] = c + i*n*p + j*p - pstart;
2868:   for (i=0; i<m; i++)
2869:     for (j=0; j<n; j++)
2870:       for (k=0; k<p; k++)
2871:         c[i*n*p+j*p+k] = (PetscScalar*) aa + i*n*p*q + j*p*q + k*q - qstart;
2872:   *a -= mstart;
2873:   return(0);
2874: }

2876: /*@C
2877:    VecRestoreArray4dRead - Restores a vector after VecGetArray3d() has been called.

2879:    Logically Collective

2881:    Input Parameters:
2882: +  x - the vector
2883: .  m - first dimension of four dimensional array
2884: .  n - second dimension of the four dimensional array
2885: .  p - third dimension of the four dimensional array
2886: .  q - fourth dimension of the four dimensional array
2887: .  mstart - first index you will use in first coordinate direction (often 0)
2888: .  nstart - first index in the second coordinate direction (often 0)
2889: .  pstart - first index in the third coordinate direction (often 0)
2890: .  qstart - first index in the fourth coordinate direction (often 0)
2891: -  a - location of pointer to array obtained from VecGetArray4dRead()

2893:    Level: beginner

2895:    Notes:
2896:    For regular PETSc vectors this routine does not involve any copies. For
2897:    any special vectors that do not store local vector data in a contiguous
2898:    array, this routine will copy the data back into the underlying
2899:    vector data structure from the array obtained with VecGetArray().

2901:    This routine actually zeros out the a pointer.

2903: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2904:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2905:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2906: @*/
2907: PetscErrorCode  VecRestoreArray4dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2908: {
2910:   void           *dummy;

2916:   dummy = (void*)(*a + mstart);
2917:   PetscFree(dummy);
2918:   VecRestoreArrayRead(x,NULL);
2919:   return(0);
2920: }

2922: #if defined(PETSC_USE_DEBUG)

2924: /*@
2925:    VecLockGet  - Gets the current lock status of a vector

2927:    Logically Collective on Vec

2929:    Input Parameter:
2930: .  x - the vector

2932:    Output Parameter:
2933: .  state - greater than zero indicates the vector is still locked

2935:    Level: beginner

2937:    Concepts: vector^accessing local values

2939: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockPush(), VecLockGet()
2940: @*/
2941: PetscErrorCode VecLockGet(Vec x,PetscInt *state)
2942: {
2945:   *state = x->lock;
2946:   return(0);
2947: }

2949: /*@
2950:    VecLockPush  - Lock a vector from writing

2952:    Logically Collective on Vec

2954:    Input Parameter:
2955: .  x - the vector

2957:    Notes:
2958:     If this is set then calls to VecGetArray() or VecSetValues() or any other routines that change the vectors values will fail.

2960:     Call VecLockPop() to remove the latest lock

2962:    Level: beginner

2964:    Concepts: vector^accessing local values

2966: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockPop(), VecLockGet()
2967: @*/
2968: PetscErrorCode VecLockPush(Vec x)
2969: {
2972:   x->lock++;
2973:   return(0);
2974: }

2976: /*@
2977:    VecLockPop  - Unlock a vector from writing

2979:    Logically Collective on Vec

2981:    Input Parameter:
2982: .  x - the vector

2984:    Level: beginner

2986:    Concepts: vector^accessing local values

2988: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockPush(), VecLockGet()
2989: @*/
2990: PetscErrorCode VecLockPop(Vec x)
2991: {
2994:   x->lock--;
2995:   if (x->lock < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector has been unlocked too many times");
2996:   return(0);
2997: }

2999: #endif