Actual source code: rvector.c

petsc-master 2019-06-15
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_CUDA) || 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


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

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

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

125: /*@
126:    VecDotRealPart - Computes the real part of the vector dot product.

128:    Collective on Vec

130:    Input Parameters:
131: .  x, y - the vectors

133:    Output Parameter:
134: .  val - the real part of the dot product;

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

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

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

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

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

151:    Level: intermediate


154: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDot(), VecDotNorm2()
155: @*/
156: PetscErrorCode  VecDotRealPart(Vec x,Vec y,PetscReal *val)
157: {
159:   PetscScalar    fdot;

162:   VecDot(x,y,&fdot);
163:   *val = PetscRealPart(fdot);
164:   return(0);
165: }

167: /*@
168:    VecNorm  - Computes the vector norm.

170:    Collective on Vec

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

178:    Output Parameter:
179: .  val - the norm

181:    Notes:
182: $     NORM_1 denotes sum_i |x_i|
183: $     NORM_2 denotes sqrt(sum_i |x_i|^2)
184: $     NORM_INFINITY denotes max_i |x_i|

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

191:    Level: intermediate

193:    Performance Issues:
194: $    per-processor memory bandwidth
195: $    interprocessor latency
196: $    work load inbalance that causes certain processes to arrive much earlier than others


199: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNormAvailable(),
200:           VecNormBegin(), VecNormEnd()

202: @*/
203: PetscErrorCode  VecNorm(Vec x,NormType type,PetscReal *val)
204: {
205:   PetscBool      flg;


213:   /*
214:    * Cached data?
215:    */
216:   if (type!=NORM_1_AND_2) {
217:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,flg);
218:     if (flg) return(0);
219:   }
220:   PetscLogEventBegin(VEC_Norm,x,0,0,0);
221:   (*x->ops->norm)(x,type,val);
222:   PetscLogEventEnd(VEC_Norm,x,0,0,0);

224:   if (type!=NORM_1_AND_2) {
225:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[type],*val);
226:   }
227:   return(0);
228: }

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

233:    Not Collective

235:    Input Parameters:
236: +  x - the vector
237: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
238:           NORM_1_AND_2, which computes both norms and stores them
239:           in a two element array.

241:    Output Parameter:
242: +  available - PETSC_TRUE if the val returned is valid
243: -  val - the norm

245:    Notes:
246: $     NORM_1 denotes sum_i |x_i|
247: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
248: $     NORM_INFINITY denotes max_i |x_i|

250:    Level: intermediate

252:    Performance Issues:
253: $    per-processor memory bandwidth
254: $    interprocessor latency
255: $    work load inbalance that causes certain processes to arrive much earlier than others

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


263: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNorm()
264:           VecNormBegin(), VecNormEnd()

266: @*/
267: PetscErrorCode  VecNormAvailable(Vec x,NormType type,PetscBool  *available,PetscReal *val)
268: {


276:   *available = PETSC_FALSE;
277:   if (type!=NORM_1_AND_2) {
278:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,*available);
279:   }
280:   return(0);
281: }

283: /*@
284:    VecNormalize - Normalizes a vector by 2-norm.

286:    Collective on Vec

288:    Input Parameters:
289: +  x - the vector

291:    Output Parameter:
292: .  x - the normalized vector
293: -  val - the vector norm before normalization

295:    Level: intermediate


298: @*/
299: PetscErrorCode  VecNormalize(Vec x,PetscReal *val)
300: {
302:   PetscReal      norm;

307:   PetscLogEventBegin(VEC_Normalize,x,0,0,0);
308:   VecNorm(x,NORM_2,&norm);
309:   if (norm == 0.0) {
310:     PetscInfo(x,"Vector of zero norm can not be normalized; Returning only the zero norm\n");
311:   } else if (norm != 1.0) {
312:     PetscScalar tmp = 1.0/norm;
313:     VecScale(x,tmp);
314:   }
315:   if (val) *val = norm;
316:   PetscLogEventEnd(VEC_Normalize,x,0,0,0);
317:   return(0);
318: }

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

323:    Collective on Vec

325:    Input Parameter:
326: .  x - the vector

328:    Output Parameters:
329: +  p - the location of val (pass NULL if you don't want this)
330: -  val - the maximum component

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

335:    Returns the smallest index with the maximum value
336:    Level: intermediate


339: .seealso: VecNorm(), VecMin()
340: @*/
341: PetscErrorCode  VecMax(Vec x,PetscInt *p,PetscReal *val)
342: {

349:   PetscLogEventBegin(VEC_Max,x,0,0,0);
350:   (*x->ops->max)(x,p,val);
351:   PetscLogEventEnd(VEC_Max,x,0,0,0);
352:   return(0);
353: }

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

358:    Collective on Vec

360:    Input Parameters:
361: .  x - the vector

363:    Output Parameter:
364: +  p - the location of val (pass NULL if you don't want this location)
365: -  val - the minimum component

367:    Level: intermediate

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

372:    This returns the smallest index with the minumum value


375: .seealso: VecMax()
376: @*/
377: PetscErrorCode  VecMin(Vec x,PetscInt *p,PetscReal *val)
378: {

385:   PetscLogEventBegin(VEC_Min,x,0,0,0);
386:   (*x->ops->min)(x,p,val);
387:   PetscLogEventEnd(VEC_Min,x,0,0,0);
388:   return(0);
389: }

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

395:    Collective on Vec

397:    Input Parameters:
398: .  x, y - the vectors

400:    Output Parameter:
401: .  val - the dot product

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

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

412:    Level: intermediate

414: .seealso: VecDot(), VecMTDot()
415: @*/
416: PetscErrorCode  VecTDot(Vec x,Vec y,PetscScalar *val)
417: {

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

429:   PetscLogEventBegin(VEC_TDot,x,y,0,0);
430:   (*x->ops->tdot)(x,y,val);
431:   PetscLogEventEnd(VEC_TDot,x,y,0,0);
432:   return(0);
433: }

435: /*@
436:    VecScale - Scales a vector.

438:    Not collective on Vec

440:    Input Parameters:
441: +  x - the vector
442: -  alpha - the scalar

444:    Output Parameter:
445: .  x - the scaled vector

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

451:    Level: intermediate


454: @*/
455: PetscErrorCode  VecScale(Vec x, PetscScalar alpha)
456: {
457:   PetscReal      norms[4] = {0.0,0.0,0.0, 0.0};
458:   PetscBool      flgs[4];
460:   PetscInt       i;

465:   if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
466:   PetscLogEventBegin(VEC_Scale,x,0,0,0);
467:   if (alpha != (PetscScalar)1.0) {
468:     /* get current stashed norms */
469:     for (i=0; i<4; i++) {
470:       PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],norms[i],flgs[i]);
471:     }
472:     (*x->ops->scale)(x,alpha);
473:     PetscObjectStateIncrease((PetscObject)x);
474:     /* put the scaled stashed norms back into the Vec */
475:     for (i=0; i<4; i++) {
476:       if (flgs[i]) {
477:         PetscObjectComposedDataSetReal((PetscObject)x,NormIds[i],PetscAbsScalar(alpha)*norms[i]);
478:       }
479:     }
480:   }
481:   PetscLogEventEnd(VEC_Scale,x,0,0,0);
482:   return(0);
483: }

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

488:    Logically Collective on Vec

490:    Input Parameters:
491: +  x  - the vector
492: -  alpha - the scalar

494:    Output Parameter:
495: .  x  - the vector

497:    Note:
498:    For a vector of dimension n, VecSet() computes
499: $     x[i] = alpha, for i=1,...,n,
500:    so that all vector entries then equal the identical
501:    scalar value, alpha.  Use the more general routine
502:    VecSetValues() to set different vector entries.

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

507:    Level: beginner

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

511: @*/
512: PetscErrorCode  VecSet(Vec x,PetscScalar alpha)
513: {
514:   PetscReal      val;

520:   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()");
522:   VecSetErrorIfLocked(x,1);

524:   PetscLogEventBegin(VEC_Set,x,0,0,0);
525:   (*x->ops->set)(x,alpha);
526:   PetscLogEventEnd(VEC_Set,x,0,0,0);
527:   PetscObjectStateIncrease((PetscObject)x);

529:   /*  norms can be simply set (if |alpha|*N not too large) */
530:   val  = PetscAbsScalar(alpha);
531:   if (x->map->N == 0) {
532:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],0.0l);
533:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],0.0);
534:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],0.0);
535:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],0.0);
536:   } else if (val > PETSC_MAX_REAL/x->map->N) {
537:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
538:   } else {
539:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],x->map->N * val);
540:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
541:     val  = PetscSqrtReal((PetscReal)x->map->N) * val;
542:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],val);
543:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],val);
544:   }
545:   return(0);
546: }


549: /*@
550:    VecAXPY - Computes y = alpha x + y.

552:    Logically Collective on Vec

554:    Input Parameters:
555: +  alpha - the scalar
556: -  x, y  - the vectors

558:    Output Parameter:
559: .  y - output vector

561:    Level: intermediate

563:    Notes:
564:     x and y MUST be different vectors
565:     This routine is optimized for alpha of 0.0, otherwise it calls the BLAS routine

567: $    VecAXPY(y,alpha,x)                   y = alpha x           +      y
568: $    VecAYPX(y,beta,x)                    y =       x           + beta y
569: $    VecAXPBY(y,alpha,beta,x)             y = alpha x           + beta y
570: $    VecWAXPY(w,alpha,x,y)                w = alpha x           +      y
571: $    VecAXPBYPCZ(w,alpha,beta,gamma,x,y)  z = alpha x           + beta y + gamma z
572: $    VecMAXPY(y,nv,alpha[],x[])           y = sum alpha[i] x[i] +      y


575: .seealso:  VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPBYPCZ(), VecAXPBY()
576: @*/
577: PetscErrorCode  VecAXPY(Vec y,PetscScalar alpha,Vec x)
578: {

587:   VecCheckSameSize(x,1,y,3);
588:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");
590:   VecSetErrorIfLocked(y,1);

592:   VecLockReadPush(x);
593:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
594:   (*y->ops->axpy)(y,alpha,x);
595:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
596:   VecLockReadPop(x);
597:   PetscObjectStateIncrease((PetscObject)y);
598:   return(0);
599: }

601: /*@
602:    VecAXPBY - Computes y = alpha x + beta y.

604:    Logically Collective on Vec

606:    Input Parameters:
607: +  alpha,beta - the scalars
608: -  x, y  - the vectors

610:    Output Parameter:
611: .  y - output vector

613:    Level: intermediate

615:    Notes:
616:     x and y MUST be different vectors
617:     The implementation is optimized for alpha and/or beta values of 0.0 and 1.0


620: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY(), VecAXPBYPCZ()
621: @*/
622: PetscErrorCode  VecAXPBY(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
623: {

632:   VecCheckSameSize(y,1,x,4);
633:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");
636:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
637:   (*y->ops->axpby)(y,alpha,beta,x);
638:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
639:   PetscObjectStateIncrease((PetscObject)y);
640:   return(0);
641: }

643: /*@
644:    VecAXPBYPCZ - Computes z = alpha x + beta y + gamma z

646:    Logically Collective on Vec

648:    Input Parameters:
649: +  alpha,beta, gamma - the scalars
650: -  x, y, z  - the vectors

652:    Output Parameter:
653: .  z - output vector

655:    Level: intermediate

657:    Notes:
658:     x, y and z must be different vectors
659:     The implementation is optimized for alpha of 1.0 and gamma of 1.0 or 0.0


662: .seealso:  VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY(), VecAXPBY()
663: @*/
664: PetscErrorCode  VecAXPBYPCZ(Vec z,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec x,Vec y)
665: {

677:   VecCheckSameSize(x,1,y,5);
678:   VecCheckSameSize(x,1,z,6);
679:   if (x == y || x == z) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");
680:   if (y == z) SETERRQ(PetscObjectComm((PetscObject)y),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");

685:   PetscLogEventBegin(VEC_AXPBYPCZ,x,y,z,0);
686:   (*y->ops->axpbypcz)(z,alpha,beta,gamma,x,y);
687:   PetscLogEventEnd(VEC_AXPBYPCZ,x,y,z,0);
688:   PetscObjectStateIncrease((PetscObject)z);
689:   return(0);
690: }

692: /*@
693:    VecAYPX - Computes y = x + beta y.

695:    Logically Collective on Vec

697:    Input Parameters:
698: +  beta - the scalar
699: -  x, y  - the vectors

701:    Output Parameter:
702: .  y - output vector

704:    Level: intermediate

706:    Notes:
707:     x and y MUST be different vectors
708:     The implementation is optimized for beta of -1.0, 0.0, and 1.0


711: .seealso:  VecMAXPY(), VecWAXPY(), VecAXPY(), VecAXPBYPCZ(), VecAXPBY()
712: @*/
713: PetscErrorCode  VecAYPX(Vec y,PetscScalar beta,Vec x)
714: {

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

727:   PetscLogEventBegin(VEC_AYPX,x,y,0,0);
728:    (*y->ops->aypx)(y,beta,x);
729:   PetscLogEventEnd(VEC_AYPX,x,y,0,0);
730:   PetscObjectStateIncrease((PetscObject)y);
731:   return(0);
732: }


735: /*@
736:    VecWAXPY - Computes w = alpha x + y.

738:    Logically Collective on Vec

740:    Input Parameters:
741: +  alpha - the scalar
742: -  x, y  - the vectors

744:    Output Parameter:
745: .  w - the result

747:    Level: intermediate

749:    Notes:
750:     w cannot be either x or y, but x and y can be the same
751:     The implementation is optimzed for alpha of -1.0, 0.0, and 1.0


754: .seealso: VecAXPY(), VecAYPX(), VecAXPBY(), VecMAXPY(), VecAXPBYPCZ()
755: @*/
756: PetscErrorCode  VecWAXPY(Vec w,PetscScalar alpha,Vec x,Vec y)
757: {

769:   VecCheckSameSize(x,3,y,4);
770:   VecCheckSameSize(x,3,w,1);
771:   if (w == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector y, suggest VecAXPY()");
772:   if (w == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector x, suggest VecAYPX()");

775:   PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
776:    (*w->ops->waxpy)(w,alpha,x,y);
777:   PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
778:   PetscObjectStateIncrease((PetscObject)w);
779:   return(0);
780: }


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

786:    Not Collective

788:    Input Parameters:
789: +  x - vector to insert in
790: .  ni - number of elements to add
791: .  ix - indices where to add
792: .  y - array of values
793: -  iora - either INSERT_VALUES or ADD_VALUES, where
794:    ADD_VALUES adds values to any existing entries, and
795:    INSERT_VALUES replaces existing entries with new values

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

800:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
801:    options cannot be mixed without intervening calls to the assembly
802:    routines.

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

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

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

815:    Level: beginner

817: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
818:            VecSetValue(), VecSetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecGetValues()
819: @*/
820: PetscErrorCode  VecSetValues(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
821: {

826:   if (!ni) return(0);
830:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
831:   (*x->ops->setvalues)(x,ni,ix,y,iora);
832:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
833:   PetscObjectStateIncrease((PetscObject)x);
834:   return(0);
835: }

837: /*@
838:    VecGetValues - Gets values from certain locations of a vector. Currently
839:           can only get values on the same processor

841:     Not Collective

843:    Input Parameters:
844: +  x - vector to get values from
845: .  ni - number of elements to get
846: -  ix - indices where to get them from (in global 1d numbering)

848:    Output Parameter:
849: .   y - array of values

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

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

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

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

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

864:    Level: beginner

866: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues()
867: @*/
868: PetscErrorCode  VecGetValues(Vec x,PetscInt ni,const PetscInt ix[],PetscScalar y[])
869: {

874:   if (!ni) return(0);
878:   (*x->ops->getvalues)(x,ni,ix,y);
879:   return(0);
880: }

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

885:    Not Collective

887:    Input Parameters:
888: +  x - vector to insert in
889: .  ni - number of blocks to add
890: .  ix - indices where to add in block count, rather than element count
891: .  y - array of values
892: -  iora - either INSERT_VALUES or ADD_VALUES, where
893:    ADD_VALUES adds values to any existing entries, and
894:    INSERT_VALUES replaces existing entries with new values

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

900:    Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
901:    options cannot be mixed without intervening calls to the assembly
902:    routines.

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

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

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

914:    Level: intermediate

916: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
917:            VecSetValues()
918: @*/
919: PetscErrorCode  VecSetValuesBlocked(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
920: {

925:   if (!ni) return(0);
929:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
930:   (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
931:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
932:   PetscObjectStateIncrease((PetscObject)x);
933:   return(0);
934: }


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

941:    Not Collective

943:    Input Parameters:
944: +  x - vector to insert in
945: .  ni - number of elements to add
946: .  ix - indices where to add
947: .  y - array of values
948: -  iora - either INSERT_VALUES or ADD_VALUES, where
949:    ADD_VALUES adds values to any existing entries, and
950:    INSERT_VALUES replaces existing entries with new values

952:    Level: intermediate

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

957:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
958:    options cannot be mixed without intervening calls to the assembly
959:    routines.

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

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

966: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
967:            VecSetValuesBlockedLocal()
968: @*/
969: PetscErrorCode  VecSetValuesLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
970: {
972:   PetscInt       lixp[128],*lix = lixp;

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

981:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
982:   if (!x->ops->setvalueslocal) {
983:     if (!x->map->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
984:     if (ni > 128) {
985:       PetscMalloc1(ni,&lix);
986:     }
987:     ISLocalToGlobalMappingApply(x->map->mapping,ni,(PetscInt*)ix,lix);
988:     (*x->ops->setvalues)(x,ni,lix,y,iora);
989:     if (ni > 128) {
990:       PetscFree(lix);
991:     }
992:   } else {
993:     (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
994:   }
995:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
996:   PetscObjectStateIncrease((PetscObject)x);
997:   return(0);
998: }

1000: /*@
1001:    VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1002:    using a local ordering of the nodes.

1004:    Not Collective

1006:    Input Parameters:
1007: +  x - vector to insert in
1008: .  ni - number of blocks to add
1009: .  ix - indices where to add in block count, not element count
1010: .  y - array of values
1011: -  iora - either INSERT_VALUES or ADD_VALUES, where
1012:    ADD_VALUES adds values to any existing entries, and
1013:    INSERT_VALUES replaces existing entries with new values

1015:    Level: intermediate

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

1021:    Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1022:    options cannot be mixed without intervening calls to the assembly
1023:    routines.

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

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


1031: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1032:            VecSetLocalToGlobalMapping()
1033: @*/
1034: PetscErrorCode  VecSetValuesBlockedLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1035: {
1037:   PetscInt       lixp[128],*lix = lixp;

1041:   if (!ni) return(0);
1045:   if (ni > 128) {
1046:     PetscMalloc1(ni,&lix);
1047:   }

1049:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1050:   ISLocalToGlobalMappingApplyBlock(x->map->mapping,ni,(PetscInt*)ix,lix);
1051:   (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1052:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1053:   if (ni > 128) {
1054:     PetscFree(lix);
1055:   }
1056:   PetscObjectStateIncrease((PetscObject)x);
1057:   return(0);
1058: }

1060: /*@
1061:    VecMTDot - Computes indefinite vector multiple dot products.
1062:    That is, it does NOT use the complex conjugate.

1064:    Collective on Vec

1066:    Input Parameters:
1067: +  x - one vector
1068: .  nv - number of vectors
1069: -  y - array of vectors.  Note that vectors are pointers

1071:    Output Parameter:
1072: .  val - array of the dot products

1074:    Notes for Users of Complex Numbers:
1075:    For complex vectors, VecMTDot() computes the indefinite form
1076: $      val = (x,y) = y^T x,
1077:    where y^T denotes the transpose of y.

1079:    Use VecMDot() for the inner product
1080: $      val = (x,y) = y^H x,
1081:    where y^H denotes the conjugate transpose of y.

1083:    Level: intermediate


1086: .seealso: VecMDot(), VecTDot()
1087: @*/
1088: PetscErrorCode  VecMTDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1089: {

1095:   if (!nv) return(0);
1102:   VecCheckSameSize(x,1,*y,3);

1104:   PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1105:   (*x->ops->mtdot)(x,nv,y,val);
1106:   PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1107:   return(0);
1108: }

1110: /*@
1111:    VecMDot - Computes vector multiple dot products.

1113:    Collective on Vec

1115:    Input Parameters:
1116: +  x - one vector
1117: .  nv - number of vectors
1118: -  y - array of vectors.

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

1123:    Notes for Users of Complex Numbers:
1124:    For complex vectors, VecMDot() computes
1125: $     val = (x,y) = y^H x,
1126:    where y^H denotes the conjugate transpose of y.

1128:    Use VecMTDot() for the indefinite form
1129: $     val = (x,y) = y^T x,
1130:    where y^T denotes the transpose of y.

1132:    Level: intermediate


1135: .seealso: VecMTDot(), VecDot()
1136: @*/
1137: PetscErrorCode  VecMDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1138: {

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

1154:   PetscLogEventBegin(VEC_MDot,x,*y,0,0);
1155:   (*x->ops->mdot)(x,nv,y,val);
1156:   PetscLogEventEnd(VEC_MDot,x,*y,0,0);
1157:   return(0);
1158: }

1160: /*@
1161:    VecMAXPY - Computes y = y + sum alpha[i] x[i]

1163:    Logically Collective on Vec

1165:    Input Parameters:
1166: +  nv - number of scalars and x-vectors
1167: .  alpha - array of scalars
1168: .  y - one vector
1169: -  x - array of vectors

1171:    Level: intermediate

1173:    Notes:
1174:     y cannot be any of the x vectors

1176: .seealso:  VecAYPX(), VecWAXPY(), VecAXPY(), VecAXPBYPCZ(), VecAXPBY()
1177: @*/
1178: PetscErrorCode  VecMAXPY(Vec y,PetscInt nv,const PetscScalar alpha[],Vec x[])
1179: {
1181:   PetscInt       i;

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

1197:   PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1198:   (*y->ops->maxpy)(y,nv,alpha,x);
1199:   PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1200:   PetscObjectStateIncrease((PetscObject)y);
1201:   return(0);
1202: }

1204: /*@
1205:    VecGetSubVector - Gets a vector representing part of another vector

1207:    Collective on IS

1209:    Input Arguments:
1210: + X - vector from which to extract a subvector
1211: - is - index set representing portion of X to extract

1213:    Output Arguments:
1214: . Y - subvector corresponding to is

1216:    Level: advanced

1218:    Notes:
1219:    The subvector Y should be returned with VecRestoreSubVector().

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

1224: .seealso: MatCreateSubMatrix()
1225: @*/
1226: PetscErrorCode  VecGetSubVector(Vec X,IS is,Vec *Y)
1227: {
1228:   PetscErrorCode   ierr;
1229:   Vec              Z;

1235:   if (X->ops->getsubvector) {
1236:     (*X->ops->getsubvector)(X,is,&Z);
1237:   } else { /* Default implementation currently does no caching */
1238:     PetscInt  gstart,gend,start;
1239:     PetscBool contiguous,gcontiguous;
1240:     VecGetOwnershipRange(X,&gstart,&gend);
1241:     ISContiguousLocal(is,gstart,gend,&start,&contiguous);
1242:     MPIU_Allreduce(&contiguous,&gcontiguous,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)is));
1243:     if (gcontiguous) { /* We can do a no-copy implementation */
1244:       PetscInt n,N,bs;
1245:       PetscInt state;

1247:       ISGetSize(is,&N);
1248:       ISGetLocalSize(is,&n);
1249:       VecGetBlockSize(X,&bs);
1250:       if (n%bs || bs == 1 || !n) bs = -1; /* Do not decide block size if we do not have to */
1251:       VecLockGet(X,&state);
1252:       if (state) {
1253:         const PetscScalar *x;
1254:         VecGetArrayRead(X,&x);
1255:         VecCreate(PetscObjectComm((PetscObject)X),&Z);
1256:         VecSetType(Z,((PetscObject)X)->type_name);
1257:         VecSetSizes(Z,n,N);
1258:         VecSetBlockSize(Z,bs);
1259:         VecPlaceArray(Z,(PetscScalar*)x+start);
1260:         VecLockReadPush(Z);
1261:         VecRestoreArrayRead(X,&x);
1262:       } else {
1263:         PetscScalar *x;
1264:         VecGetArray(X,&x);
1265:         VecCreate(PetscObjectComm((PetscObject)X),&Z);
1266:         VecSetType(Z,((PetscObject)X)->type_name);
1267:         VecSetSizes(Z,n,N);
1268:         VecSetBlockSize(Z,bs);
1269:         VecPlaceArray(Z,(PetscScalar*)x+start);
1270:         VecRestoreArray(X,&x);
1271:       }
1272:     } else { /* Have to create a scatter and do a copy */
1273:       VecScatter scatter;
1274:       PetscInt   n,N;
1275:       ISGetLocalSize(is,&n);
1276:       ISGetSize(is,&N);
1277:       VecCreate(PetscObjectComm((PetscObject)is),&Z);
1278:       VecSetSizes(Z,n,N);
1279:       VecSetType(Z,((PetscObject)X)->type_name);
1280:       VecScatterCreate(X,is,Z,NULL,&scatter);
1281:       VecScatterBegin(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1282:       VecScatterEnd(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1283:       PetscObjectCompose((PetscObject)Z,"VecGetSubVector_Scatter",(PetscObject)scatter);
1284:       VecScatterDestroy(&scatter);
1285:     }
1286:   }
1287:   /* Record the state when the subvector was gotten so we know whether its values need to be put back */
1288:   if (VecGetSubVectorSavedStateId < 0) {PetscObjectComposedDataRegister(&VecGetSubVectorSavedStateId);}
1289:   PetscObjectComposedDataSetInt((PetscObject)Z,VecGetSubVectorSavedStateId,1);
1290:   *Y   = Z;
1291:   return(0);
1292: }

1294: /*@
1295:    VecRestoreSubVector - Restores a subvector extracted using VecGetSubVector()

1297:    Collective on IS

1299:    Input Arguments:
1300: + X - vector from which subvector was obtained
1301: . is - index set representing the subset of X
1302: - Y - subvector being restored

1304:    Level: advanced

1306: .seealso: VecGetSubVector()
1307: @*/
1308: PetscErrorCode  VecRestoreSubVector(Vec X,IS is,Vec *Y)
1309: {

1317:   if (X->ops->restoresubvector) {
1318:     (*X->ops->restoresubvector)(X,is,Y);
1319:   } else {
1320:     PETSC_UNUSED PetscObjectState dummystate = 0;
1321:     PetscBool valid;
1322:     PetscObjectComposedDataGetInt((PetscObject)*Y,VecGetSubVectorSavedStateId,dummystate,valid);
1323:     if (!valid) {
1324:       VecScatter scatter;

1326:       PetscObjectQuery((PetscObject)*Y,"VecGetSubVector_Scatter",(PetscObject*)&scatter);
1327:       if (scatter) {
1328:         VecScatterBegin(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1329:         VecScatterEnd(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1330:       }
1331:     }
1332:     VecDestroy(Y);
1333:   }
1334:   return(0);
1335: }

1337: /*@
1338:    VecGetLocalVectorRead - Maps the local portion of a vector into a
1339:    vector.  You must call VecRestoreLocalVectorRead() when the local
1340:    vector is no longer needed.

1342:    Not collective.

1344:    Input parameter:
1345: .  v - The vector for which the local vector is desired.

1347:    Output parameter:
1348: .  w - Upon exit this contains the local vector.

1350:    Level: beginner

1352:    Notes:
1353:    This function is similar to VecGetArrayRead() which maps the local
1354:    portion into a raw pointer.  VecGetLocalVectorRead() is usually
1355:    almost as efficient as VecGetArrayRead() but in certain circumstances
1356:    VecGetLocalVectorRead() can be much more efficient than
1357:    VecGetArrayRead().  This is because the construction of a contiguous
1358:    array representing the vector data required by VecGetArrayRead() can
1359:    be an expensive operation for certain vector types.  For example, for
1360:    GPU vectors VecGetArrayRead() requires that the data between device
1361:    and host is synchronized.

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

1366: .seealso: VecRestoreLocalVectorRead(), VecGetLocalVector(), VecGetArrayRead(), VecGetArray()
1367: @*/
1368: PetscErrorCode VecGetLocalVectorRead(Vec v,Vec w)
1369: {
1371:   PetscScalar    *a;

1376:   VecCheckSameLocalSize(v,1,w,2);
1377:   if (v->ops->getlocalvectorread) {
1378:     (*v->ops->getlocalvectorread)(v,w);
1379:   } else {
1380:     VecGetArrayRead(v,(const PetscScalar**)&a);
1381:     VecPlaceArray(w,a);
1382:   }
1383:   return(0);
1384: }

1386: /*@
1387:    VecRestoreLocalVectorRead - Unmaps the local portion of a vector
1388:    previously mapped into a vector using VecGetLocalVectorRead().

1390:    Not collective.

1392:    Input parameter:
1393: .  v - The local portion of this vector was previously mapped into w using VecGetLocalVectorRead().
1394: .  w - The vector into which the local portion of v was mapped.

1396:    Level: beginner

1398: .seealso: VecGetLocalVectorRead(), VecGetLocalVector(), VecGetArrayRead(), VecGetArray()
1399: @*/
1400: PetscErrorCode VecRestoreLocalVectorRead(Vec v,Vec w)
1401: {
1403:   PetscScalar    *a;

1408:   if (v->ops->restorelocalvectorread) {
1409:     (*v->ops->restorelocalvectorread)(v,w);
1410:   } else {
1411:     VecGetArrayRead(w,(const PetscScalar**)&a);
1412:     VecRestoreArrayRead(v,(const PetscScalar**)&a);
1413:     VecResetArray(w);
1414:   }
1415:   return(0);
1416: }

1418: /*@
1419:    VecGetLocalVector - Maps the local portion of a vector into a
1420:    vector.

1422:    Collective on v, not collective on w.

1424:    Input parameter:
1425: .  v - The vector for which the local vector is desired.

1427:    Output parameter:
1428: .  w - Upon exit this contains the local vector.

1430:    Level: beginner

1432:    Notes:
1433:    This function is similar to VecGetArray() which maps the local
1434:    portion into a raw pointer.  VecGetLocalVector() is usually about as
1435:    efficient as VecGetArray() but in certain circumstances
1436:    VecGetLocalVector() can be much more efficient than VecGetArray().
1437:    This is because the construction of a contiguous array representing
1438:    the vector data required by VecGetArray() can be an expensive
1439:    operation for certain vector types.  For example, for GPU vectors
1440:    VecGetArray() requires that the data between device and host is
1441:    synchronized.

1443: .seealso: VecRestoreLocalVector(), VecGetLocalVectorRead(), VecGetArrayRead(), VecGetArray()
1444: @*/
1445: PetscErrorCode VecGetLocalVector(Vec v,Vec w)
1446: {
1448:   PetscScalar    *a;

1453:   VecCheckSameLocalSize(v,1,w,2);
1454:   if (v->ops->getlocalvector) {
1455:     (*v->ops->getlocalvector)(v,w);
1456:   } else {
1457:     VecGetArray(v,&a);
1458:     VecPlaceArray(w,a);
1459:   }
1460:   return(0);
1461: }

1463: /*@
1464:    VecRestoreLocalVector - Unmaps the local portion of a vector
1465:    previously mapped into a vector using VecGetLocalVector().

1467:    Logically collective.

1469:    Input parameter:
1470: .  v - The local portion of this vector was previously mapped into w using VecGetLocalVector().
1471: .  w - The vector into which the local portion of v was mapped.

1473:    Level: beginner

1475: .seealso: VecGetLocalVector(), VecGetLocalVectorRead(), VecRestoreLocalVectorRead(), LocalVectorRead(), VecGetArrayRead(), VecGetArray()
1476: @*/
1477: PetscErrorCode VecRestoreLocalVector(Vec v,Vec w)
1478: {
1480:   PetscScalar    *a;

1485:   if (v->ops->restorelocalvector) {
1486:     (*v->ops->restorelocalvector)(v,w);
1487:   } else {
1488:     VecGetArray(w,&a);
1489:     VecRestoreArray(v,&a);
1490:     VecResetArray(w);
1491:   }
1492:   return(0);
1493: }

1495: /*@C
1496:    VecGetArray - Returns a pointer to a contiguous array that contains this
1497:    processor's portion of the vector data. For the standard PETSc
1498:    vectors, VecGetArray() returns a pointer to the local data array and
1499:    does not use any copies. If the underlying vector data is not stored
1500:    in a contiguous array this routine will copy the data to a contiguous
1501:    array and return a pointer to that. You MUST call VecRestoreArray()
1502:    when you no longer need access to the array.

1504:    Logically Collective on Vec

1506:    Input Parameter:
1507: .  x - the vector

1509:    Output Parameter:
1510: .  a - location to put pointer to the array

1512:    Fortran Note:
1513:    This routine is used differently from Fortran 77
1514: $    Vec         x
1515: $    PetscScalar x_array(1)
1516: $    PetscOffset i_x
1517: $    PetscErrorCode ierr
1518: $       call VecGetArray(x,x_array,i_x,ierr)
1519: $
1520: $   Access first local entry in vector with
1521: $      value = x_array(i_x + 1)
1522: $
1523: $      ...... other code
1524: $       call VecRestoreArray(x,x_array,i_x,ierr)
1525:    For Fortran 90 see VecGetArrayF90()

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

1530:    Level: beginner

1532: .seealso: VecRestoreArray(), VecGetArrayRead(), VecGetArrays(), VecGetArrayF90(), VecGetArrayReadF90(), VecPlaceArray(), VecGetArray2d(),
1533:           VecGetArrayPair(), VecRestoreArrayPair()
1534: @*/
1535: PetscErrorCode VecGetArray(Vec x,PetscScalar **a)
1536: {
1538: #if defined(PETSC_HAVE_VIENNACL)
1539:   PetscBool      is_viennacltype = PETSC_FALSE;
1540: #endif

1544:   VecSetErrorIfLocked(x,1);
1545:   if (x->petscnative) {
1546: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1547:     if (x->valid_GPU_array == PETSC_OFFLOAD_GPU) {
1548: #if defined(PETSC_HAVE_VIENNACL)
1549:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1550:       if (is_viennacltype) {
1551:         VecViennaCLCopyFromGPU(x);
1552:       } else
1553: #endif
1554:       {
1555: #if defined(PETSC_HAVE_CUDA)
1556:         VecCUDACopyFromGPU(x);
1557: #endif
1558:       }
1559:     } else if (x->valid_GPU_array == PETSC_OFFLOAD_UNALLOCATED) {
1560: #if defined(PETSC_HAVE_VIENNACL)
1561:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1562:       if (is_viennacltype) {
1563:         VecViennaCLAllocateCheckHost(x);
1564:       } else
1565: #endif
1566:       {
1567: #if defined(PETSC_HAVE_CUDA)
1568:         VecCUDAAllocateCheckHost(x);
1569: #endif
1570:       }
1571:     }
1572: #endif
1573:     *a = *((PetscScalar**)x->data);
1574:   } else {
1575:     if (x->ops->getarray) {
1576:       (*x->ops->getarray)(x,a);
1577:     } else SETERRQ1(PetscObjectComm((PetscObject)x),PETSC_ERR_SUP,"Cannot get array for vector type \"%s\"",((PetscObject)x)->type_name);
1578:   }
1579:   return(0);
1580: }

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

1585:    Not Collective

1587:    Input Parameters:
1588: .  x - the vector

1590:    Output Parameter:
1591: .  a - the array

1593:    Level: beginner

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

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

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

1604: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrayPair(), VecRestoreArrayPair()
1605: @*/
1606: PetscErrorCode VecGetArrayRead(Vec x,const PetscScalar **a)
1607: {
1609: #if defined(PETSC_HAVE_VIENNACL)
1610:   PetscBool      is_viennacltype = PETSC_FALSE;
1611: #endif

1615:   if (x->petscnative) {
1616: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1617:     if (x->valid_GPU_array == PETSC_OFFLOAD_GPU) {
1618: #if defined(PETSC_HAVE_VIENNACL)
1619:       PetscObjectTypeCompareAny((PetscObject)x,&is_viennacltype,VECSEQVIENNACL,VECMPIVIENNACL,VECVIENNACL,"");
1620:       if (is_viennacltype) {
1621:         VecViennaCLCopyFromGPU(x);
1622:       } else
1623: #endif
1624:       {
1625: #if defined(PETSC_HAVE_CUDA)
1626:         VecCUDACopyFromGPU(x);
1627: #endif
1628:       }
1629:     }
1630: #endif
1631:     *a = *((PetscScalar **)x->data);
1632:   } else if (x->ops->getarrayread) {
1633:     (*x->ops->getarrayread)(x,a);
1634:   } else {
1635:     (*x->ops->getarray)(x,(PetscScalar**)a);
1636:   }
1637:   return(0);
1638: }

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

1645:    Logically Collective on Vec

1647:    Input Parameter:
1648: +  x - the vectors
1649: -  n - the number of vectors

1651:    Output Parameter:
1652: .  a - location to put pointer to the array

1654:    Fortran Note:
1655:    This routine is not supported in Fortran.

1657:    Level: intermediate

1659: .seealso: VecGetArray(), VecRestoreArrays()
1660: @*/
1661: PetscErrorCode  VecGetArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1662: {
1664:   PetscInt       i;
1665:   PetscScalar    **q;

1671:   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %D",n);
1672:   PetscMalloc1(n,&q);
1673:   for (i=0; i<n; ++i) {
1674:     VecGetArray(x[i],&q[i]);
1675:   }
1676:   *a = q;
1677:   return(0);
1678: }

1680: /*@C
1681:    VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1682:    has been called.

1684:    Logically Collective on Vec

1686:    Input Parameters:
1687: +  x - the vector
1688: .  n - the number of vectors
1689: -  a - location of pointer to arrays obtained from VecGetArrays()

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

1697:    Fortran Note:
1698:    This routine is not supported in Fortran.

1700:    Level: intermediate

1702: .seealso: VecGetArrays(), VecRestoreArray()
1703: @*/
1704: PetscErrorCode  VecRestoreArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1705: {
1707:   PetscInt       i;
1708:   PetscScalar    **q = *a;


1715:   for (i=0; i<n; ++i) {
1716:     VecRestoreArray(x[i],&q[i]);
1717:   }
1718:   PetscFree(q);
1719:   return(0);
1720: }

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

1725:    Logically Collective on Vec

1727:    Input Parameters:
1728: +  x - the vector
1729: -  a - location of pointer to array obtained from VecGetArray()

1731:    Level: beginner

1733:    Notes:
1734:    For regular PETSc vectors this routine does not involve any copies. For
1735:    any special vectors that do not store local vector data in a contiguous
1736:    array, this routine will copy the data back into the underlying
1737:    vector data structure from the array obtained with VecGetArray().

1739:    This routine actually zeros out the a pointer. This is to prevent accidental
1740:    use of the array after it has been restored. If you pass null for a it will
1741:    not zero the array pointer a.

1743:    Fortran Note:
1744:    This routine is used differently from Fortran 77
1745: $    Vec         x
1746: $    PetscScalar x_array(1)
1747: $    PetscOffset i_x
1748: $    PetscErrorCode ierr
1749: $       call VecGetArray(x,x_array,i_x,ierr)
1750: $
1751: $   Access first local entry in vector with
1752: $      value = x_array(i_x + 1)
1753: $
1754: $      ...... other code
1755: $       call VecRestoreArray(x,x_array,i_x,ierr)

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

1761: .seealso: VecGetArray(), VecRestoreArrayRead(), VecRestoreArrays(), VecRestoreArrayF90(), VecRestoreArrayReadF90(), VecPlaceArray(), VecRestoreArray2d(),
1762:           VecGetArrayPair(), VecRestoreArrayPair()
1763: @*/
1764: PetscErrorCode VecRestoreArray(Vec x,PetscScalar **a)
1765: {

1770:   if (x->petscnative) {
1771: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1772:     x->valid_GPU_array = PETSC_OFFLOAD_CPU;
1773: #endif
1774:   } else {
1775:     (*x->ops->restorearray)(x,a);
1776:   }
1777:   if (a) *a = NULL;
1778:   PetscObjectStateIncrease((PetscObject)x);
1779:   return(0);
1780: }

1782: /*@C
1783:    VecRestoreArrayRead - Restore array obtained with VecGetArrayRead()

1785:    Not Collective

1787:    Input Parameters:
1788: +  vec - the vector
1789: -  array - the array

1791:    Level: beginner

1793: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrayPair(), VecRestoreArrayPair()
1794: @*/
1795: PetscErrorCode VecRestoreArrayRead(Vec x,const PetscScalar **a)
1796: {

1801:   if (x->petscnative) {
1802:     /* nothing */
1803:   } else if (x->ops->restorearrayread) {
1804:     (*x->ops->restorearrayread)(x,a);
1805:   } else {
1806:     (*x->ops->restorearray)(x,(PetscScalar**)a);
1807:   }
1808:   if (a) *a = NULL;
1809:   return(0);
1810: }

1812: /*@
1813:    VecPlaceArray - Allows one to replace the array in a vector with an
1814:    array provided by the user. This is useful to avoid copying an array
1815:    into a vector.

1817:    Not Collective

1819:    Input Parameters:
1820: +  vec - the vector
1821: -  array - the array

1823:    Notes:
1824:    You can return to the original array with a call to VecResetArray()

1826:    Level: developer

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

1830: @*/
1831: PetscErrorCode  VecPlaceArray(Vec vec,const PetscScalar array[])
1832: {

1839:   if (vec->ops->placearray) {
1840:     (*vec->ops->placearray)(vec,array);
1841:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot place array in this type of vector");
1842:   PetscObjectStateIncrease((PetscObject)vec);
1843:   return(0);
1844: }

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

1851:    Not Collective

1853:    Input Parameters:
1854: +  vec - the vector
1855: -  array - the array

1857:    Notes:
1858:    This permanently replaces the array and frees the memory associated
1859:    with the old array.

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

1864:    Not supported from Fortran

1866:    Level: developer

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

1870: @*/
1871: PetscErrorCode  VecReplaceArray(Vec vec,const PetscScalar array[])
1872: {

1878:   if (vec->ops->replacearray) {
1879:     (*vec->ops->replacearray)(vec,array);
1880:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot replace array in this type of vector");
1881:   PetscObjectStateIncrease((PetscObject)vec);
1882:   return(0);
1883: }

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

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

1892:     Collective on Vec

1894:     Input Parameters:
1895: +   x - a vector to mimic
1896: -   n - the number of vectors to obtain

1898:     Output Parameters:
1899: +   y - Fortran90 pointer to the array of vectors
1900: -   ierr - error code

1902:     Example of Usage:
1903: .vb
1904: #include <petsc/finclude/petscvec.h>
1905:     use petscvec

1907:     Vec x
1908:     Vec, pointer :: y(:)
1909:     ....
1910:     call VecDuplicateVecsF90(x,2,y,ierr)
1911:     call VecSet(y(2),alpha,ierr)
1912:     call VecSet(y(2),alpha,ierr)
1913:     ....
1914:     call VecDestroyVecsF90(2,y,ierr)
1915: .ve

1917:     Notes:
1918:     Not yet supported for all F90 compilers

1920:     Use VecDestroyVecsF90() to free the space.

1922:     Level: beginner

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

1926: M*/

1928: /*MC
1929:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
1930:     VecGetArrayF90().

1932:     Synopsis:
1933:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

1935:     Logically Collective on Vec

1937:     Input Parameters:
1938: +   x - vector
1939: -   xx_v - the Fortran90 pointer to the array

1941:     Output Parameter:
1942: .   ierr - error code

1944:     Example of Usage:
1945: .vb
1946: #include <petsc/finclude/petscvec.h>
1947:     use petscvec

1949:     PetscScalar, pointer :: xx_v(:)
1950:     ....
1951:     call VecGetArrayF90(x,xx_v,ierr)
1952:     xx_v(3) = a
1953:     call VecRestoreArrayF90(x,xx_v,ierr)
1954: .ve

1956:     Level: beginner

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

1960: M*/

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

1965:     Synopsis:
1966:     VecDestroyVecsF90(PetscInt n,{Vec, pointer :: x(:)},PetscErrorCode ierr)

1968:     Collective on Vec

1970:     Input Parameters:
1971: +   n - the number of vectors previously obtained
1972: -   x - pointer to array of vector pointers

1974:     Output Parameter:
1975: .   ierr - error code

1977:     Notes:
1978:     Not yet supported for all F90 compilers

1980:     Level: beginner

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

1984: M*/

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

1992:     Synopsis:
1993:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

1995:     Logically Collective on Vec

1997:     Input Parameter:
1998: .   x - vector

2000:     Output Parameters:
2001: +   xx_v - the Fortran90 pointer to the array
2002: -   ierr - error code

2004:     Example of Usage:
2005: .vb
2006: #include <petsc/finclude/petscvec.h>
2007:     use petscvec

2009:     PetscScalar, pointer :: xx_v(:)
2010:     ....
2011:     call VecGetArrayF90(x,xx_v,ierr)
2012:     xx_v(3) = a
2013:     call VecRestoreArrayF90(x,xx_v,ierr)
2014: .ve

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

2018:     Level: beginner

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

2022: M*/

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

2030:     Synopsis:
2031:     VecGetArrayReadF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

2033:     Logically Collective on Vec

2035:     Input Parameter:
2036: .   x - vector

2038:     Output Parameters:
2039: +   xx_v - the Fortran90 pointer to the array
2040: -   ierr - error code

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

2047:     PetscScalar, pointer :: xx_v(:)
2048:     ....
2049:     call VecGetArrayReadF90(x,xx_v,ierr)
2050:     a = xx_v(3)
2051:     call VecRestoreArrayReadF90(x,xx_v,ierr)
2052: .ve

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

2056:     Level: beginner

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

2060: M*/

2062: /*MC
2063:     VecRestoreArrayReadF90 - Restores a readonly vector to a usable state after a call to
2064:     VecGetArrayReadF90().

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

2069:     Logically Collective on Vec

2071:     Input Parameters:
2072: +   x - vector
2073: -   xx_v - the Fortran90 pointer to the array

2075:     Output Parameter:
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:     Level: beginner

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

2094: M*/

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

2101:    Logically Collective

2103:    Input Parameter:
2104: +  x - the vector
2105: .  m - first dimension of two dimensional array
2106: .  n - second dimension of two dimensional array
2107: .  mstart - first index you will use in first coordinate direction (often 0)
2108: -  nstart - first index in the second coordinate direction (often 0)

2110:    Output Parameter:
2111: .  a - location to put pointer to the array

2113:    Level: developer

2115:   Notes:
2116:    For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
2117:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2118:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2119:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().

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

2123: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2124:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2125:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2126: @*/
2127: PetscErrorCode  VecGetArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2128: {
2130:   PetscInt       i,N;
2131:   PetscScalar    *aa;

2137:   VecGetLocalSize(x,&N);
2138:   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);
2139:   VecGetArray(x,&aa);

2141:   PetscMalloc1(m,a);
2142:   for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
2143:   *a -= mstart;
2144:   return(0);
2145: }

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

2150:    Logically Collective

2152:    Input Parameters:
2153: +  x - the vector
2154: .  m - first dimension of two dimensional array
2155: .  n - second dimension of the two dimensional array
2156: .  mstart - first index you will use in first coordinate direction (often 0)
2157: .  nstart - first index in the second coordinate direction (often 0)
2158: -  a - location of pointer to array obtained from VecGetArray2d()

2160:    Level: developer

2162:    Notes:
2163:    For regular PETSc vectors this routine does not involve any copies. For
2164:    any special vectors that do not store local vector data in a contiguous
2165:    array, this routine will copy the data back into the underlying
2166:    vector data structure from the array obtained with VecGetArray().

2168:    This routine actually zeros out the a pointer.

2170: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2171:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2172:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2173: @*/
2174: PetscErrorCode  VecRestoreArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2175: {
2177:   void           *dummy;

2183:   dummy = (void*)(*a + mstart);
2184:   PetscFree(dummy);
2185:   VecRestoreArray(x,NULL);
2186:   return(0);
2187: }

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

2194:    Logically Collective

2196:    Input Parameter:
2197: +  x - the vector
2198: .  m - first dimension of two dimensional array
2199: -  mstart - first index you will use in first coordinate direction (often 0)

2201:    Output Parameter:
2202: .  a - location to put pointer to the array

2204:    Level: developer

2206:   Notes:
2207:    For a vector obtained from DMCreateLocalVector() mstart are likely
2208:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2209:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().

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

2213: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2214:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2215:           VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2216: @*/
2217: PetscErrorCode  VecGetArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2218: {
2220:   PetscInt       N;

2226:   VecGetLocalSize(x,&N);
2227:   if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2228:   VecGetArray(x,a);
2229:   *a  -= mstart;
2230:   return(0);
2231: }

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

2236:    Logically Collective

2238:    Input Parameters:
2239: +  x - the vector
2240: .  m - first dimension of two dimensional array
2241: .  mstart - first index you will use in first coordinate direction (often 0)
2242: -  a - location of pointer to array obtained from VecGetArray21()

2244:    Level: developer

2246:    Notes:
2247:    For regular PETSc vectors this routine does not involve any copies. For
2248:    any special vectors that do not store local vector data in a contiguous
2249:    array, this routine will copy the data back into the underlying
2250:    vector data structure from the array obtained with VecGetArray1d().

2252:    This routine actually zeros out the a pointer.

2254: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2255:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2256:           VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2257: @*/
2258: PetscErrorCode  VecRestoreArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2259: {

2265:   VecRestoreArray(x,NULL);
2266:   return(0);
2267: }


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

2275:    Logically Collective

2277:    Input Parameter:
2278: +  x - the vector
2279: .  m - first dimension of three dimensional array
2280: .  n - second dimension of three dimensional array
2281: .  p - third dimension of three dimensional array
2282: .  mstart - first index you will use in first coordinate direction (often 0)
2283: .  nstart - first index in the second coordinate direction (often 0)
2284: -  pstart - first index in the third coordinate direction (often 0)

2286:    Output Parameter:
2287: .  a - location to put pointer to the array

2289:    Level: developer

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

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

2299: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2300:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2301:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2302: @*/
2303: PetscErrorCode  VecGetArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2304: {
2306:   PetscInt       i,N,j;
2307:   PetscScalar    *aa,**b;

2313:   VecGetLocalSize(x,&N);
2314:   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);
2315:   VecGetArray(x,&aa);

2317:   PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2318:   b    = (PetscScalar**)((*a) + m);
2319:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2320:   for (i=0; i<m; i++)
2321:     for (j=0; j<n; j++)
2322:       b[i*n+j] = aa + i*n*p + j*p - pstart;

2324:   *a -= mstart;
2325:   return(0);
2326: }

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

2331:    Logically Collective

2333:    Input Parameters:
2334: +  x - the vector
2335: .  m - first dimension of three dimensional array
2336: .  n - second dimension of the three dimensional array
2337: .  p - third dimension of the three dimensional array
2338: .  mstart - first index you will use in first coordinate direction (often 0)
2339: .  nstart - first index in the second coordinate direction (often 0)
2340: .  pstart - first index in the third coordinate direction (often 0)
2341: -  a - location of pointer to array obtained from VecGetArray3d()

2343:    Level: developer

2345:    Notes:
2346:    For regular PETSc vectors this routine does not involve any copies. For
2347:    any special vectors that do not store local vector data in a contiguous
2348:    array, this routine will copy the data back into the underlying
2349:    vector data structure from the array obtained with VecGetArray().

2351:    This routine actually zeros out the a pointer.

2353: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2354:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2355:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2356: @*/
2357: PetscErrorCode  VecRestoreArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2358: {
2360:   void           *dummy;

2366:   dummy = (void*)(*a + mstart);
2367:   PetscFree(dummy);
2368:   VecRestoreArray(x,NULL);
2369:   return(0);
2370: }

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

2377:    Logically Collective

2379:    Input Parameter:
2380: +  x - the vector
2381: .  m - first dimension of four dimensional array
2382: .  n - second dimension of four dimensional array
2383: .  p - third dimension of four dimensional array
2384: .  q - fourth dimension of four dimensional array
2385: .  mstart - first index you will use in first coordinate direction (often 0)
2386: .  nstart - first index in the second coordinate direction (often 0)
2387: .  pstart - first index in the third coordinate direction (often 0)
2388: -  qstart - first index in the fourth coordinate direction (often 0)

2390:    Output Parameter:
2391: .  a - location to put pointer to the array

2393:    Level: beginner

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

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

2403: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2404:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2405:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2406: @*/
2407: PetscErrorCode  VecGetArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2408: {
2410:   PetscInt       i,N,j,k;
2411:   PetscScalar    *aa,***b,**c;

2417:   VecGetLocalSize(x,&N);
2418:   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);
2419:   VecGetArray(x,&aa);

2421:   PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2422:   b    = (PetscScalar***)((*a) + m);
2423:   c    = (PetscScalar**)(b + m*n);
2424:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2425:   for (i=0; i<m; i++)
2426:     for (j=0; j<n; j++)
2427:       b[i*n+j] = c + i*n*p + j*p - pstart;
2428:   for (i=0; i<m; i++)
2429:     for (j=0; j<n; j++)
2430:       for (k=0; k<p; k++)
2431:         c[i*n*p+j*p+k] = aa + i*n*p*q + j*p*q + k*q - qstart;
2432:   *a -= mstart;
2433:   return(0);
2434: }

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

2439:    Logically Collective

2441:    Input Parameters:
2442: +  x - the vector
2443: .  m - first dimension of four dimensional array
2444: .  n - second dimension of the four dimensional array
2445: .  p - third dimension of the four dimensional array
2446: .  q - fourth dimension of the four dimensional array
2447: .  mstart - first index you will use in first coordinate direction (often 0)
2448: .  nstart - first index in the second coordinate direction (often 0)
2449: .  pstart - first index in the third coordinate direction (often 0)
2450: .  qstart - first index in the fourth coordinate direction (often 0)
2451: -  a - location of pointer to array obtained from VecGetArray4d()

2453:    Level: beginner

2455:    Notes:
2456:    For regular PETSc vectors this routine does not involve any copies. For
2457:    any special vectors that do not store local vector data in a contiguous
2458:    array, this routine will copy the data back into the underlying
2459:    vector data structure from the array obtained with VecGetArray().

2461:    This routine actually zeros out the a pointer.

2463: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2464:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2465:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2466: @*/
2467: PetscErrorCode  VecRestoreArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2468: {
2470:   void           *dummy;

2476:   dummy = (void*)(*a + mstart);
2477:   PetscFree(dummy);
2478:   VecRestoreArray(x,NULL);
2479:   return(0);
2480: }

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

2487:    Logically Collective

2489:    Input Parameter:
2490: +  x - the vector
2491: .  m - first dimension of two dimensional array
2492: .  n - second dimension of two dimensional array
2493: .  mstart - first index you will use in first coordinate direction (often 0)
2494: -  nstart - first index in the second coordinate direction (often 0)

2496:    Output Parameter:
2497: .  a - location to put pointer to the array

2499:    Level: developer

2501:   Notes:
2502:    For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
2503:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2504:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2505:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().

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

2509: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2510:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2511:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2512: @*/
2513: PetscErrorCode  VecGetArray2dRead(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2514: {
2515:   PetscErrorCode    ierr;
2516:   PetscInt          i,N;
2517:   const PetscScalar *aa;

2523:   VecGetLocalSize(x,&N);
2524:   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);
2525:   VecGetArrayRead(x,&aa);

2527:   PetscMalloc1(m,a);
2528:   for (i=0; i<m; i++) (*a)[i] = (PetscScalar*) aa + i*n - nstart;
2529:   *a -= mstart;
2530:   return(0);
2531: }

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

2536:    Logically Collective

2538:    Input Parameters:
2539: +  x - the vector
2540: .  m - first dimension of two dimensional array
2541: .  n - second dimension of the two dimensional array
2542: .  mstart - first index you will use in first coordinate direction (often 0)
2543: .  nstart - first index in the second coordinate direction (often 0)
2544: -  a - location of pointer to array obtained from VecGetArray2d()

2546:    Level: developer

2548:    Notes:
2549:    For regular PETSc vectors this routine does not involve any copies. For
2550:    any special vectors that do not store local vector data in a contiguous
2551:    array, this routine will copy the data back into the underlying
2552:    vector data structure from the array obtained with VecGetArray().

2554:    This routine actually zeros out the a pointer.

2556: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2557:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2558:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2559: @*/
2560: PetscErrorCode  VecRestoreArray2dRead(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
2561: {
2563:   void           *dummy;

2569:   dummy = (void*)(*a + mstart);
2570:   PetscFree(dummy);
2571:   VecRestoreArrayRead(x,NULL);
2572:   return(0);
2573: }

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

2580:    Logically Collective

2582:    Input Parameter:
2583: +  x - the vector
2584: .  m - first dimension of two dimensional array
2585: -  mstart - first index you will use in first coordinate direction (often 0)

2587:    Output Parameter:
2588: .  a - location to put pointer to the array

2590:    Level: developer

2592:   Notes:
2593:    For a vector obtained from DMCreateLocalVector() mstart are likely
2594:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2595:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().

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

2599: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2600:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2601:           VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2602: @*/
2603: PetscErrorCode  VecGetArray1dRead(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2604: {
2606:   PetscInt       N;

2612:   VecGetLocalSize(x,&N);
2613:   if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2614:   VecGetArrayRead(x,(const PetscScalar**)a);
2615:   *a  -= mstart;
2616:   return(0);
2617: }

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

2622:    Logically Collective

2624:    Input Parameters:
2625: +  x - the vector
2626: .  m - first dimension of two dimensional array
2627: .  mstart - first index you will use in first coordinate direction (often 0)
2628: -  a - location of pointer to array obtained from VecGetArray21()

2630:    Level: developer

2632:    Notes:
2633:    For regular PETSc vectors this routine does not involve any copies. For
2634:    any special vectors that do not store local vector data in a contiguous
2635:    array, this routine will copy the data back into the underlying
2636:    vector data structure from the array obtained with VecGetArray1dRead().

2638:    This routine actually zeros out the a pointer.

2640: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2641:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2642:           VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2643: @*/
2644: PetscErrorCode  VecRestoreArray1dRead(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2645: {

2651:   VecRestoreArrayRead(x,NULL);
2652:   return(0);
2653: }


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

2661:    Logically Collective

2663:    Input Parameter:
2664: +  x - the vector
2665: .  m - first dimension of three dimensional array
2666: .  n - second dimension of three dimensional array
2667: .  p - third dimension of three dimensional array
2668: .  mstart - first index you will use in first coordinate direction (often 0)
2669: .  nstart - first index in the second coordinate direction (often 0)
2670: -  pstart - first index in the third coordinate direction (often 0)

2672:    Output Parameter:
2673: .  a - location to put pointer to the array

2675:    Level: developer

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

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

2685: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2686:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2687:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2688: @*/
2689: PetscErrorCode  VecGetArray3dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2690: {
2691:   PetscErrorCode    ierr;
2692:   PetscInt          i,N,j;
2693:   const PetscScalar *aa;
2694:   PetscScalar       **b;

2700:   VecGetLocalSize(x,&N);
2701:   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);
2702:   VecGetArrayRead(x,&aa);

2704:   PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2705:   b    = (PetscScalar**)((*a) + m);
2706:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2707:   for (i=0; i<m; i++)
2708:     for (j=0; j<n; j++)
2709:       b[i*n+j] = (PetscScalar *)aa + i*n*p + j*p - pstart;

2711:   *a -= mstart;
2712:   return(0);
2713: }

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

2718:    Logically Collective

2720:    Input Parameters:
2721: +  x - the vector
2722: .  m - first dimension of three dimensional array
2723: .  n - second dimension of the three dimensional array
2724: .  p - third dimension of the three dimensional array
2725: .  mstart - first index you will use in first coordinate direction (often 0)
2726: .  nstart - first index in the second coordinate direction (often 0)
2727: .  pstart - first index in the third coordinate direction (often 0)
2728: -  a - location of pointer to array obtained from VecGetArray3dRead()

2730:    Level: developer

2732:    Notes:
2733:    For regular PETSc vectors this routine does not involve any copies. For
2734:    any special vectors that do not store local vector data in a contiguous
2735:    array, this routine will copy the data back into the underlying
2736:    vector data structure from the array obtained with VecGetArray().

2738:    This routine actually zeros out the a pointer.

2740: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2741:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2742:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2743: @*/
2744: PetscErrorCode  VecRestoreArray3dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2745: {
2747:   void           *dummy;

2753:   dummy = (void*)(*a + mstart);
2754:   PetscFree(dummy);
2755:   VecRestoreArrayRead(x,NULL);
2756:   return(0);
2757: }

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

2764:    Logically Collective

2766:    Input Parameter:
2767: +  x - the vector
2768: .  m - first dimension of four dimensional array
2769: .  n - second dimension of four dimensional array
2770: .  p - third dimension of four dimensional array
2771: .  q - fourth dimension of four dimensional array
2772: .  mstart - first index you will use in first coordinate direction (often 0)
2773: .  nstart - first index in the second coordinate direction (often 0)
2774: .  pstart - first index in the third coordinate direction (often 0)
2775: -  qstart - first index in the fourth coordinate direction (often 0)

2777:    Output Parameter:
2778: .  a - location to put pointer to the array

2780:    Level: beginner

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

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

2790: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2791:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2792:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2793: @*/
2794: PetscErrorCode  VecGetArray4dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2795: {
2796:   PetscErrorCode    ierr;
2797:   PetscInt          i,N,j,k;
2798:   const PetscScalar *aa;
2799:   PetscScalar       ***b,**c;

2805:   VecGetLocalSize(x,&N);
2806:   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);
2807:   VecGetArrayRead(x,&aa);

2809:   PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2810:   b    = (PetscScalar***)((*a) + m);
2811:   c    = (PetscScalar**)(b + m*n);
2812:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2813:   for (i=0; i<m; i++)
2814:     for (j=0; j<n; j++)
2815:       b[i*n+j] = c + i*n*p + j*p - pstart;
2816:   for (i=0; i<m; i++)
2817:     for (j=0; j<n; j++)
2818:       for (k=0; k<p; k++)
2819:         c[i*n*p+j*p+k] = (PetscScalar*) aa + i*n*p*q + j*p*q + k*q - qstart;
2820:   *a -= mstart;
2821:   return(0);
2822: }

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

2827:    Logically Collective

2829:    Input Parameters:
2830: +  x - the vector
2831: .  m - first dimension of four dimensional array
2832: .  n - second dimension of the four dimensional array
2833: .  p - third dimension of the four dimensional array
2834: .  q - fourth dimension of the four dimensional array
2835: .  mstart - first index you will use in first coordinate direction (often 0)
2836: .  nstart - first index in the second coordinate direction (often 0)
2837: .  pstart - first index in the third coordinate direction (often 0)
2838: .  qstart - first index in the fourth coordinate direction (often 0)
2839: -  a - location of pointer to array obtained from VecGetArray4dRead()

2841:    Level: beginner

2843:    Notes:
2844:    For regular PETSc vectors this routine does not involve any copies. For
2845:    any special vectors that do not store local vector data in a contiguous
2846:    array, this routine will copy the data back into the underlying
2847:    vector data structure from the array obtained with VecGetArray().

2849:    This routine actually zeros out the a pointer.

2851: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2852:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2853:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2854: @*/
2855: PetscErrorCode  VecRestoreArray4dRead(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2856: {
2858:   void           *dummy;

2864:   dummy = (void*)(*a + mstart);
2865:   PetscFree(dummy);
2866:   VecRestoreArrayRead(x,NULL);
2867:   return(0);
2868: }

2870: #if defined(PETSC_USE_DEBUG)

2872: /*@
2873:    VecLockGet  - Gets the current lock status of a vector

2875:    Logically Collective on Vec

2877:    Input Parameter:
2878: .  x - the vector

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

2884:    Level: beginner

2886: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockReadPush(), VecLockReadPop()
2887: @*/
2888: PetscErrorCode VecLockGet(Vec x,PetscInt *state)
2889: {
2892:   *state = x->lock;
2893:   return(0);
2894: }

2896: /*@
2897:    VecLockReadPush  - Pushes a read-only lock on a vector to prevent it from writing

2899:    Logically Collective on Vec

2901:    Input Parameter:
2902: .  x - the vector

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

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

2910:    Level: beginner

2912: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockReadPop(), VecLockGet()
2913: @*/
2914: PetscErrorCode VecLockReadPush(Vec x)
2915: {
2918:   if (x->lock < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector is already locked for exclusive write access but you want to read it");
2919:   x->lock++;
2920:   return(0);
2921: }

2923: /*@
2924:    VecLockReadPop  - Pops a read-only lock from a vector

2926:    Logically Collective on Vec

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

2931:    Level: beginner

2933: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockReadPush(), VecLockGet()
2934: @*/
2935: PetscErrorCode VecLockReadPop(Vec x)
2936: {
2939:   x->lock--;
2940:   if (x->lock < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector has been unlocked from read-only access too many times");
2941:   return(0);
2942: }

2944: /*@C
2945:    VecLockWriteSet_Private  - Lock or unlock a vector for exclusive read/write access

2947:    Logically Collective on Vec

2949:    Input Parameter:
2950: +  x   - the vector
2951: -  flg - PETSC_TRUE to lock the vector for writing; PETSC_FALSE to unlock it.

2953:    Notes:
2954:     The function is usefull in split-phase computations, which usually have a begin phase and an end phase.
2955:     One can call VecLockWriteSet_Private(x,PETSC_TRUE) in the begin phase to lock a vector for exclusive
2956:     access, and call VecLockWriteSet_Private(x,PETSC_FALSE) in the end phase to unlock the vector from exclusive
2957:     access. In this way, one is ensured no other operations can access the vector in between. The code may like


2960:        VecGetArray(x,&xdata); // begin phase
2961:        VecLockWriteSet_Private(v,PETSC_TRUE);

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

2965:        VecRestoreArray(x,&vdata); // end phase
2966:        VecLockWriteSet_Private(v,PETSC_FALSE);

2968:     The call can not be nested on the same vector, in other words, one can not call VecLockWriteSet_Private(x,PETSC_TRUE)
2969:     again before calling VecLockWriteSet_Private(v,PETSC_FALSE).

2971:    Level: beginner

2973: .seealso: VecRestoreArray(), VecGetArrayRead(), VecLockReadPush(), VecLockReadPop(), VecLockGet()
2974: @*/
2975: PetscErrorCode VecLockWriteSet_Private(Vec x,PetscBool flg)
2976: {
2979:   if (flg) {
2980:     if (x->lock > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector is already locked for read-only access but you want to write it");
2981:     else if (x->lock < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector is already locked for exclusive write access but you want to write it");
2982:     else x->lock = -1;
2983:   } else {
2984:     if (x->lock != -1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Vector is not locked for exclusive write access but you want to unlock it from that");
2985:     x->lock = 0;
2986:   }
2987:   return(0);
2988: }

2990: /*@
2991:    VecLockPush  - Pushes a read-only lock on a vector to prevent it from writing

2993:    Level: deprecated

2995: .seealso: VecLockReadPush()
2996: @*/
2997: PetscErrorCode VecLockPush(Vec x)
2998: {
3001:   VecLockReadPush(x);
3002:   return(0);
3003: }

3005: /*@
3006:    VecLockPop  - Pops a read-only lock from a vector

3008:    Level: deprecated

3010: .seealso: VecLockReadPop()
3011: @*/
3012: PetscErrorCode VecLockPop(Vec x)
3013: {
3016:   VecLockReadPop(x);
3017:   return(0);
3018: }

3020: #endif