Actual source code: rvector.c

petsc-3.5.0 2014-06-30
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>       /*I  "petscvec.h"   I*/
  7: static PetscInt VecGetSubVectorSavedStateId = -1;

 10:   if ((x)->map->N != (y)->map->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths %d != %d", (x)->map->N, (y)->map->N); \
 11:   if ((x)->map->n != (y)->map->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths %d != %d", (x)->map->n, (y)->map->n);

 15: PETSC_EXTERN PetscErrorCode VecValidValues(Vec vec,PetscInt argnum,PetscBool begin)
 16: {
 17: #if defined(PETSC_USE_DEBUG)
 18:   PetscErrorCode    ierr;
 19:   PetscInt          n,i;
 20:   const PetscScalar *x;

 23: #if defined(PETSC_HAVE_CUSP)
 24:   if ((vec->petscnative || vec->ops->getarray) && vec->valid_GPU_array != PETSC_CUSP_GPU) {
 25: #else
 26:   if (vec->petscnative || vec->ops->getarray) {
 27: #endif
 28:     VecGetLocalSize(vec,&n);
 29:     VecGetArrayRead(vec,&x);
 30:     for (i=0; i<n; i++) {
 31:       if (begin) {
 32:         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);
 33:       } else {
 34:         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);
 35:       }
 36:     }
 37:     VecRestoreArrayRead(vec,&x);
 38:   }
 39:   return(0);
 40: #else
 41:   return 0;
 42: #endif
 43: }

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

 50:    Logically Collective on Vec

 52:    Input Parameters:
 53: .  x, y  - the vectors

 55:    Output Parameter:
 56: .  max - the result

 58:    Level: advanced

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

 63: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMax(), VecPointwiseMin(), VecPointwiseMaxAbs()
 64: @*/
 65: PetscErrorCode  VecMaxPointwiseDivide(Vec x,Vec y,PetscReal *max)
 66: {


 78:   (*x->ops->maxpointwisedivide)(x,y,max);
 79:   return(0);
 80: }

 84: /*@
 85:    VecDot - Computes the vector dot product.

 87:    Collective on Vec

 89:    Input Parameters:
 90: .  x, y - the vectors

 92:    Output Parameter:
 93: .  val - the dot product

 95:    Performance Issues:
 96: $    per-processor memory bandwidth
 97: $    interprocessor latency
 98: $    work load inbalance that causes certain processes to arrive much earlier than others

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

107:    Use VecTDot() for the indefinite form
108: $     val = (x,y) = y^T x,
109:    where y^T denotes the transpose of y.

111:    Level: intermediate

113:    Concepts: inner product
114:    Concepts: vector^inner product

116: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDotRealPart()
117: @*/
118: PetscErrorCode  VecDot(Vec x,Vec y,PetscScalar *val)
119: {


131:   PetscLogEventBarrierBegin(VEC_DotBarrier,x,y,0,0,PetscObjectComm((PetscObject)x));
132:   (*x->ops->dot)(x,y,val);
133:   PetscLogEventBarrierEnd(VEC_DotBarrier,x,y,0,0,PetscObjectComm((PetscObject)x));
134:   return(0);
135: }

139: /*@
140:    VecDotRealPart - Computes the real part of the vector dot product.

142:    Collective on Vec

144:    Input Parameters:
145: .  x, y - the vectors

147:    Output Parameter:
148: .  val - the real part of the dot product;

150:    Performance Issues:
151: $    per-processor memory bandwidth
152: $    interprocessor latency
153: $    work load inbalance that causes certain processes to arrive much earlier than others

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

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

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

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

165:    Level: intermediate

167:    Concepts: inner product
168:    Concepts: vector^inner product

170: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecDot(), VecDotNorm2()
171: @*/
172: PetscErrorCode  VecDotRealPart(Vec x,Vec y,PetscReal *val)
173: {
175:   PetscScalar    fdot;

178:   VecDot(x,y,&fdot);
179:   *val = PetscRealPart(fdot);
180:   return(0);
181: }

185: /*@
186:    VecNorm  - Computes the vector norm.

188:    Collective on Vec

190:    Input Parameters:
191: +  x - the vector
192: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
193:           NORM_1_AND_2, which computes both norms and stores them
194:           in a two element array.

196:    Output Parameter:
197: .  val - the norm

199:    Notes:
200: $     NORM_1 denotes sum_i |x_i|
201: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
202: $     NORM_INFINITY denotes max_i |x_i|

204:    Level: intermediate

206:    Performance Issues:
207: $    per-processor memory bandwidth
208: $    interprocessor latency
209: $    work load inbalance that causes certain processes to arrive much earlier than others

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

216:    Concepts: norm
217:    Concepts: vector^norm

219: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNormAvailable(),
220:           VecNormBegin(), VecNormEnd()

222: @*/
223: PetscErrorCode  VecNorm(Vec x,NormType type,PetscReal *val)
224: {
225:   PetscBool      flg;

232:   if (((PetscObject)x)->precision != sizeof(PetscReal)) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_SUP,"Wrong precision of input argument");

234:   /*
235:    * Cached data?
236:    */
237:   if (type!=NORM_1_AND_2) {
238:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,flg);
239:     if (flg) return(0);
240:   }
241:   PetscLogEventBarrierBegin(VEC_NormBarrier,x,0,0,0,PetscObjectComm((PetscObject)x));
242:   (*x->ops->norm)(x,type,val);
243:   PetscLogEventBarrierEnd(VEC_NormBarrier,x,0,0,0,PetscObjectComm((PetscObject)x));

245:   if (type!=NORM_1_AND_2) {
246:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[type],*val);
247:   }
248:   return(0);
249: }

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

256:    Not Collective

258:    Input Parameters:
259: +  x - the vector
260: -  type - one of NORM_1, NORM_2, NORM_INFINITY.  Also available
261:           NORM_1_AND_2, which computes both norms and stores them
262:           in a two element array.

264:    Output Parameter:
265: +  available - PETSC_TRUE if the val returned is valid
266: -  val - the norm

268:    Notes:
269: $     NORM_1 denotes sum_i |x_i|
270: $     NORM_2 denotes sqrt(sum_i (x_i)^2)
271: $     NORM_INFINITY denotes max_i |x_i|

273:    Level: intermediate

275:    Performance Issues:
276: $    per-processor memory bandwidth
277: $    interprocessor latency
278: $    work load inbalance that causes certain processes to arrive much earlier than others

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

285:    Concepts: norm
286:    Concepts: vector^norm

288: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNorm()
289:           VecNormBegin(), VecNormEnd()

291: @*/
292: PetscErrorCode  VecNormAvailable(Vec x,NormType type,PetscBool  *available,PetscReal *val)
293: {


301:   *available = PETSC_FALSE;
302:   if (type!=NORM_1_AND_2) {
303:     PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,*available);
304:   }
305:   return(0);
306: }

310: /*@
311:    VecNormalize - Normalizes a vector by 2-norm.

313:    Collective on Vec

315:    Input Parameters:
316: +  x - the vector

318:    Output Parameter:
319: .  x - the normalized vector
320: -  val - the vector norm before normalization

322:    Level: intermediate

324:    Concepts: vector^normalizing
325:    Concepts: normalizing^vector

327: @*/
328: PetscErrorCode  VecNormalize(Vec x,PetscReal *val)
329: {
331:   PetscReal      norm;

336:   PetscLogEventBegin(VEC_Normalize,x,0,0,0);
337:   VecNorm(x,NORM_2,&norm);
338:   if (norm == 0.0) {
339:     PetscInfo(x,"Vector of zero norm can not be normalized; Returning only the zero norm\n");
340:   } else if (norm != 1.0) {
341:     PetscScalar tmp = 1.0/norm;
342:     VecScale(x,tmp);
343:   }
344:   if (val) *val = norm;
345:   PetscLogEventEnd(VEC_Normalize,x,0,0,0);
346:   return(0);
347: }

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

354:    Collective on Vec

356:    Input Parameter:
357: .  x - the vector

359:    Output Parameters:
360: +  val - the maximum component
361: -  p - the location of val (pass NULL if you don't want this)

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

366:    Returns the smallest index with the maximum value
367:    Level: intermediate

369:    Concepts: maximum^of vector
370:    Concepts: vector^maximum value

372: .seealso: VecNorm(), VecMin()
373: @*/
374: PetscErrorCode  VecMax(Vec x,PetscInt *p,PetscReal *val)
375: {

382:   PetscLogEventBegin(VEC_Max,x,0,0,0);
383:   (*x->ops->max)(x,p,val);
384:   PetscLogEventEnd(VEC_Max,x,0,0,0);
385:   return(0);
386: }

390: /*@
391:    VecMin - Determines the minimum vector component and its location.

393:    Collective on Vec

395:    Input Parameters:
396: .  x - the vector

398:    Output Parameter:
399: +  val - the minimum component
400: -  p - the location of val (pass NULL if you don't want this location)

402:    Level: intermediate

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

407:    This returns the smallest index with the minumum value

409:    Concepts: minimum^of vector
410:    Concepts: vector^minimum entry

412: .seealso: VecMax()
413: @*/
414: PetscErrorCode  VecMin(Vec x,PetscInt *p,PetscReal *val)
415: {

422:   PetscLogEventBegin(VEC_Min,x,0,0,0);
423:   (*x->ops->min)(x,p,val);
424:   PetscLogEventEnd(VEC_Min,x,0,0,0);
425:   return(0);
426: }

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

434:    Collective on Vec

436:    Input Parameters:
437: .  x, y - the vectors

439:    Output Parameter:
440: .  val - the dot product

442:    Notes for Users of Complex Numbers:
443:    For complex vectors, VecTDot() computes the indefinite form
444: $     val = (x,y) = y^T x,
445:    where y^T denotes the transpose of y.

447:    Use VecDot() for the inner product
448: $     val = (x,y) = y^H x,
449:    where y^H denotes the conjugate transpose of y.

451:    Level: intermediate

453:    Concepts: inner product^non-Hermitian
454:    Concepts: vector^inner product
455:    Concepts: non-Hermitian inner product

457: .seealso: VecDot(), VecMTDot()
458: @*/
459: PetscErrorCode  VecTDot(Vec x,Vec y,PetscScalar *val)
460: {


472:   PetscLogEventBegin(VEC_TDot,x,y,0,0);
473:   (*x->ops->tdot)(x,y,val);
474:   PetscLogEventEnd(VEC_TDot,x,y,0,0);
475:   return(0);
476: }

480: /*@
481:    VecScale - Scales a vector.

483:    Not collective on Vec

485:    Input Parameters:
486: +  x - the vector
487: -  alpha - the scalar

489:    Output Parameter:
490: .  x - the scaled vector

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

496:    Level: intermediate

498:    Concepts: vector^scaling
499:    Concepts: scaling^vector

501: @*/
502: PetscErrorCode  VecScale(Vec x, PetscScalar alpha)
503: {
504:   PetscReal      norms[4] = {0.0,0.0,0.0, 0.0};
505:   PetscBool      flgs[4];
507:   PetscInt       i;

512:   if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
513:   PetscLogEventBegin(VEC_Scale,x,0,0,0);
514:   if (alpha != (PetscScalar)1.0) {
515:     /* get current stashed norms */
516:     for (i=0; i<4; i++) {
517:       PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],norms[i],flgs[i]);
518:     }
519:     (*x->ops->scale)(x,alpha);
520:     PetscObjectStateIncrease((PetscObject)x);
521:     /* put the scaled stashed norms back into the Vec */
522:     for (i=0; i<4; i++) {
523:       if (flgs[i]) {
524:         PetscObjectComposedDataSetReal((PetscObject)x,NormIds[i],PetscAbsScalar(alpha)*norms[i]);
525:       }
526:     }
527:   }
528:   PetscLogEventEnd(VEC_Scale,x,0,0,0);
529:   return(0);
530: }

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

537:    Logically Collective on Vec

539:    Input Parameters:
540: +  x  - the vector
541: -  alpha - the scalar

543:    Output Parameter:
544: .  x  - the vector

546:    Note:
547:    For a vector of dimension n, VecSet() computes
548: $     x[i] = alpha, for i=1,...,n,
549:    so that all vector entries then equal the identical
550:    scalar value, alpha.  Use the more general routine
551:    VecSetValues() to set different vector entries.

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

556:    Level: beginner

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

560:    Concepts: vector^setting to constant

562: @*/
563: PetscErrorCode  VecSet(Vec x,PetscScalar alpha)
564: {
565:   PetscReal      val;

571:   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()");

574:   PetscLogEventBegin(VEC_Set,x,0,0,0);
575:   (*x->ops->set)(x,alpha);
576:   PetscLogEventEnd(VEC_Set,x,0,0,0);
577:   PetscObjectStateIncrease((PetscObject)x);

579:   /*  norms can be simply set (if |alpha|*N not too large) */
580:   val  = PetscAbsScalar(alpha);
581:   if (val > PETSC_MAX_REAL/x->map->N) {
582:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
583:   } else {
584:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],x->map->N * val);
585:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
586:     val  = PetscSqrtReal((PetscReal)x->map->N) * val;
587:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],val);
588:     PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],val);
589:   }
590:   return(0);
591: }


596: /*@
597:    VecAXPY - Computes y = alpha x + y.

599:    Logically Collective on Vec

601:    Input Parameters:
602: +  alpha - the scalar
603: -  x, y  - the vectors

605:    Output Parameter:
606: .  y - output vector

608:    Level: intermediate

610:    Notes: x and y MUST be different vectors

612:    Concepts: vector^BLAS
613:    Concepts: BLAS

615: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
616: @*/
617: PetscErrorCode  VecAXPY(Vec y,PetscScalar alpha,Vec x)
618: {

628:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");

631:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
632:   (*y->ops->axpy)(y,alpha,x);
633:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
634:   PetscObjectStateIncrease((PetscObject)y);
635:   return(0);
636: }

640: /*@
641:    VecAXPBY - Computes y = alpha x + beta y.

643:    Logically Collective on Vec

645:    Input Parameters:
646: +  alpha,beta - the scalars
647: -  x, y  - the vectors

649:    Output Parameter:
650: .  y - output vector

652:    Level: intermediate

654:    Notes: x and y MUST be different vectors

656:    Concepts: BLAS
657:    Concepts: vector^BLAS

659: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
660: @*/
661: PetscErrorCode  VecAXPBY(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
662: {

672:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");

676:   PetscLogEventBegin(VEC_AXPY,x,y,0,0);
677:   (*y->ops->axpby)(y,alpha,beta,x);
678:   PetscLogEventEnd(VEC_AXPY,x,y,0,0);
679:   PetscObjectStateIncrease((PetscObject)y);
680:   return(0);
681: }

685: /*@
686:    VecAXPBYPCZ - Computes z = alpha x + beta y + gamma z

688:    Logically Collective on Vec

690:    Input Parameters:
691: +  alpha,beta, gamma - the scalars
692: -  x, y, z  - the vectors

694:    Output Parameter:
695: .  z - output vector

697:    Level: intermediate

699:    Notes: x, y and z must be different vectors

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

703:    Concepts: BLAS
704:    Concepts: vector^BLAS

706: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
707: @*/
708: PetscErrorCode  VecAXPBYPCZ(Vec z,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec x,Vec y)
709: {

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

729:   PetscLogEventBegin(VEC_AXPBYPCZ,x,y,z,0);
730:   (*y->ops->axpbypcz)(z,alpha,beta,gamma,x,y);
731:   PetscLogEventEnd(VEC_AXPBYPCZ,x,y,z,0);
732:   PetscObjectStateIncrease((PetscObject)z);
733:   return(0);
734: }

738: /*@
739:    VecAYPX - Computes y = x + alpha y.

741:    Logically Collective on Vec

743:    Input Parameters:
744: +  alpha - the scalar
745: -  x, y  - the vectors

747:    Output Parameter:
748: .  y - output vector

750:    Level: intermediate

752:    Notes: x and y MUST be different vectors

754:    Concepts: vector^BLAS
755:    Concepts: BLAS

757: .seealso: VecAXPY(), VecWAXPY()
758: @*/
759: PetscErrorCode  VecAYPX(Vec y,PetscScalar alpha,Vec x)
760: {

768:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y must be different vectors");

771:   PetscLogEventBegin(VEC_AYPX,x,y,0,0);
772:    (*y->ops->aypx)(y,alpha,x);
773:   PetscLogEventEnd(VEC_AYPX,x,y,0,0);
774:   PetscObjectStateIncrease((PetscObject)y);
775:   return(0);
776: }


781: /*@
782:    VecWAXPY - Computes w = alpha x + y.

784:    Logically Collective on Vec

786:    Input Parameters:
787: +  alpha - the scalar
788: -  x, y  - the vectors

790:    Output Parameter:
791: .  w - the result

793:    Level: intermediate

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

797:    Concepts: vector^BLAS
798:    Concepts: BLAS

800: .seealso: VecAXPY(), VecAYPX(), VecAXPBY()
801: @*/
802: PetscErrorCode  VecWAXPY(Vec w,PetscScalar alpha,Vec x,Vec y)
803: {

817:   if (w == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector y, suggest VecAXPY()");
818:   if (w == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector x, suggest VecAYPX()");

821:   PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
822:    (*w->ops->waxpy)(w,alpha,x,y);
823:   PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
824:   PetscObjectStateIncrease((PetscObject)w);
825:   return(0);
826: }


831: /*@
832:    VecSetValues - Inserts or adds values into certain locations of a vector.

834:    Not Collective

836:    Input Parameters:
837: +  x - vector to insert in
838: .  ni - number of elements to add
839: .  ix - indices where to add
840: .  y - array of values
841: -  iora - either INSERT_VALUES or ADD_VALUES, where
842:    ADD_VALUES adds values to any existing entries, and
843:    INSERT_VALUES replaces existing entries with new values

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

848:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
849:    options cannot be mixed without intervening calls to the assembly
850:    routines.

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

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

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

863:    Level: beginner

865:    Concepts: vector^setting values

867: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
868:            VecSetValue(), VecSetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecGetValues()
869: @*/
870: PetscErrorCode  VecSetValues(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
871: {

879:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
880:   (*x->ops->setvalues)(x,ni,ix,y,iora);
881:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
882:   PetscObjectStateIncrease((PetscObject)x);
883:   return(0);
884: }

888: /*@
889:    VecGetValues - Gets values from certain locations of a vector. Currently
890:           can only get values on the same processor

892:     Not Collective

894:    Input Parameters:
895: +  x - vector to get values from
896: .  ni - number of elements to get
897: -  ix - indices where to get them from (in global 1d numbering)

899:    Output Parameter:
900: .   y - array of values

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

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

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

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

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

915:    Level: beginner

917:    Concepts: vector^getting values

919: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecGetValuesLocal(),
920:            VecGetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecSetValues()
921: @*/
922: PetscErrorCode  VecGetValues(Vec x,PetscInt ni,const PetscInt ix[],PetscScalar y[])
923: {

931:   (*x->ops->getvalues)(x,ni,ix,y);
932:   return(0);
933: }

937: /*@
938:    VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector.

940:    Not Collective

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

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

955:    Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
956:    options cannot be mixed without intervening calls to the assembly
957:    routines.

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

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

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

969:    Level: intermediate

971:    Concepts: vector^setting values blocked

973: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
974:            VecSetValues()
975: @*/
976: PetscErrorCode  VecSetValuesBlocked(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
977: {

985:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
986:   (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
987:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
988:   PetscObjectStateIncrease((PetscObject)x);
989:   return(0);
990: }


995: /*@
996:    VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
997:    using a local ordering of the nodes.

999:    Not Collective

1001:    Input Parameters:
1002: +  x - vector to insert in
1003: .  ni - number of elements to add
1004: .  ix - indices where to add
1005: .  y - array of values
1006: -  iora - either INSERT_VALUES or ADD_VALUES, where
1007:    ADD_VALUES adds values to any existing entries, and
1008:    INSERT_VALUES replaces existing entries with new values

1010:    Level: intermediate

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

1015:    Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
1016:    options cannot be mixed without intervening calls to the assembly
1017:    routines.

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

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

1024:    Concepts: vector^setting values with local numbering

1026: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
1027:            VecSetValuesBlockedLocal()
1028: @*/
1029: PetscErrorCode  VecSetValuesLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1030: {
1032:   PetscInt       lixp[128],*lix = lixp;


1040:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1041:   if (!x->ops->setvalueslocal) {
1042:     if (!x->map->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
1043:     if (ni > 128) {
1044:       PetscMalloc1(ni,&lix);
1045:     }
1046:     ISLocalToGlobalMappingApply(x->map->mapping,ni,(PetscInt*)ix,lix);
1047:     (*x->ops->setvalues)(x,ni,lix,y,iora);
1048:     if (ni > 128) {
1049:       PetscFree(lix);
1050:     }
1051:   } else {
1052:     (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
1053:   }
1054:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1055:   PetscObjectStateIncrease((PetscObject)x);
1056:   return(0);
1057: }

1061: /*@
1062:    VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1063:    using a local ordering of the nodes.

1065:    Not Collective

1067:    Input Parameters:
1068: +  x - vector to insert in
1069: .  ni - number of blocks to add
1070: .  ix - indices where to add in block count, not element count
1071: .  y - array of values
1072: -  iora - either INSERT_VALUES or ADD_VALUES, where
1073:    ADD_VALUES adds values to any existing entries, and
1074:    INSERT_VALUES replaces existing entries with new values

1076:    Level: intermediate

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

1082:    Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1083:    options cannot be mixed without intervening calls to the assembly
1084:    routines.

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

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


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

1094: .seealso:  VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1095:            VecSetLocalToGlobalMapping()
1096: @*/
1097: PetscErrorCode  VecSetValuesBlockedLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1098: {
1100:   PetscInt       lixp[128],*lix = lixp;

1107:   if (ni > 128) {
1108:     PetscMalloc1(ni,&lix);
1109:   }

1111:   PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1112:   ISLocalToGlobalMappingApplyBlock(x->map->mapping,ni,(PetscInt*)ix,lix);
1113:   (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1114:   PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1115:   if (ni > 128) {
1116:     PetscFree(lix);
1117:   }
1118:   PetscObjectStateIncrease((PetscObject)x);
1119:   return(0);
1120: }

1124: /*@
1125:    VecMTDot - Computes indefinite vector multiple dot products.
1126:    That is, it does NOT use the complex conjugate.

1128:    Collective on Vec

1130:    Input Parameters:
1131: +  x - one vector
1132: .  nv - number of vectors
1133: -  y - array of vectors.  Note that vectors are pointers

1135:    Output Parameter:
1136: .  val - array of the dot products

1138:    Notes for Users of Complex Numbers:
1139:    For complex vectors, VecMTDot() computes the indefinite form
1140: $      val = (x,y) = y^T x,
1141:    where y^T denotes the transpose of y.

1143:    Use VecMDot() for the inner product
1144: $      val = (x,y) = y^H x,
1145:    where y^H denotes the conjugate transpose of y.

1147:    Level: intermediate

1149:    Concepts: inner product^multiple
1150:    Concepts: vector^multiple inner products

1152: .seealso: VecMDot(), VecTDot()
1153: @*/
1154: PetscErrorCode  VecMTDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1155: {


1168:   PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1169:   (*x->ops->mtdot)(x,nv,y,val);
1170:   PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1171:   return(0);
1172: }

1176: /*@
1177:    VecMDot - Computes vector multiple dot products.

1179:    Collective on Vec

1181:    Input Parameters:
1182: +  x - one vector
1183: .  nv - number of vectors
1184: -  y - array of vectors.

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

1189:    Notes for Users of Complex Numbers:
1190:    For complex vectors, VecMDot() computes
1191: $     val = (x,y) = y^H x,
1192:    where y^H denotes the conjugate transpose of y.

1194:    Use VecMTDot() for the indefinite form
1195: $     val = (x,y) = y^T x,
1196:    where y^T denotes the transpose of y.

1198:    Level: intermediate

1200:    Concepts: inner product^multiple
1201:    Concepts: vector^multiple inner products

1203: .seealso: VecMTDot(), VecDot()
1204: @*/
1205: PetscErrorCode  VecMDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1206: {

1211:   if (!nv) return(0);
1212:   if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);

1221:   PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,PetscObjectComm((PetscObject)x));
1222:   (*x->ops->mdot)(x,nv,y,val);
1223:   PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,PetscObjectComm((PetscObject)x));
1224:   return(0);
1225: }

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

1232:    Logically Collective on Vec

1234:    Input Parameters:
1235: +  nv - number of scalars and x-vectors
1236: .  alpha - array of scalars
1237: .  y - one vector
1238: -  x - array of vectors

1240:    Level: intermediate

1242:    Notes: y cannot be any of the x vectors

1244:    Concepts: BLAS

1246: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1247: @*/
1248: PetscErrorCode  VecMAXPY(Vec y,PetscInt nv,const PetscScalar alpha[],Vec x[])
1249: {
1251:   PetscInt       i;

1255:   if (!nv) return(0);
1256:   if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);

1266:   PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1267:   (*y->ops->maxpy)(y,nv,alpha,x);
1268:   PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1269:   PetscObjectStateIncrease((PetscObject)y);
1270:   return(0);
1271: }

1275: /*@
1276:    VecGetSubVector - Gets a vector representing part of another vector

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

1280:    Input Arguments:
1281: + X - vector from which to extract a subvector
1282: - is - index set representing portion of X to extract

1284:    Output Arguments:
1285: . Y - subvector corresponding to is

1287:    Level: advanced

1289:    Notes:
1290:    The subvector Y should be returned with VecRestoreSubVector().

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

1295: .seealso: MatGetSubMatrix()
1296: @*/
1297: PetscErrorCode  VecGetSubVector(Vec X,IS is,Vec *Y)
1298: {
1299:   PetscErrorCode   ierr;
1300:   Vec              Z;
1301:   PetscObjectState state;

1307:   if (X->ops->getsubvector) {
1308:     (*X->ops->getsubvector)(X,is,&Z);
1309:   } else {                      /* Default implementation currently does no caching */
1310:     PetscInt  gstart,gend,start;
1311:     PetscBool contiguous,gcontiguous;
1312:     VecGetOwnershipRange(X,&gstart,&gend);
1313:     ISContiguousLocal(is,gstart,gend,&start,&contiguous);
1314:     MPI_Allreduce(&contiguous,&gcontiguous,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)is));
1315:     if (gcontiguous) {          /* We can do a no-copy implementation */
1316:       PetscInt    n,N,bs;
1317:       PetscScalar *x;
1318:       PetscMPIInt size;
1319:       ISGetLocalSize(is,&n);
1320:       VecGetArray(X,&x);
1321:       VecGetBlockSize(X,&bs);
1322:       if (n%bs) bs = 1;
1323:       MPI_Comm_size(PetscObjectComm((PetscObject)X),&size);
1324:       if (size == 1) {
1325:         VecCreateSeqWithArray(PetscObjectComm((PetscObject)X),bs,n,x+start,&Z);
1326:       } else {
1327:         ISGetSize(is,&N);
1328:         VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),bs,n,N,x+start,&Z);
1329:       }
1330:       VecRestoreArray(X,&x);
1331:     } else {                    /* Have to create a scatter and do a copy */
1332:       VecScatter scatter;
1333:       PetscInt   n,N;
1334:       ISGetLocalSize(is,&n);
1335:       ISGetSize(is,&N);
1336:       VecCreate(PetscObjectComm((PetscObject)is),&Z);
1337:       VecSetSizes(Z,n,N);
1338:       VecSetType(Z,((PetscObject)X)->type_name);
1339:       VecScatterCreate(X,is,Z,NULL,&scatter);
1340:       VecScatterBegin(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1341:       VecScatterEnd(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1342:       VecScatterDestroy(&scatter);
1343:     }
1344:   }
1345:   /* Record the state when the subvector was gotten so we know whether its values need to be put back */
1346:   if (VecGetSubVectorSavedStateId < 0) {PetscObjectComposedDataRegister(&VecGetSubVectorSavedStateId);}
1347:   PetscObjectStateGet((PetscObject)Z,&state);
1348:   PetscObjectComposedDataSetInt((PetscObject)Z,VecGetSubVectorSavedStateId,state);
1349:   *Y   = Z;
1350:   return(0);
1351: }

1355: /*@
1356:    VecRestoreSubVector - Restores a subvector extracted using VecGetSubVector()

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

1360:    Input Arguments:
1361: + X - vector from which subvector was obtained
1362: . is - index set representing the subset of X
1363: - Y - subvector being restored

1365:    Level: advanced

1367: .seealso: VecGetSubVector()
1368: @*/
1369: PetscErrorCode  VecRestoreSubVector(Vec X,IS is,Vec *Y)
1370: {

1378:   if (X->ops->restoresubvector) {
1379:     (*X->ops->restoresubvector)(X,is,Y);
1380:   } else {
1381:     PetscObjectState savedstate=0,newstate;
1382:     PetscBool valid;
1383:     PetscObjectComposedDataGetInt((PetscObject)*Y,VecGetSubVectorSavedStateId,savedstate,valid);
1384:     PetscObjectStateGet((PetscObject)*Y,&newstate);
1385:     if (valid && savedstate < newstate) {
1386:       /* We might need to copy entries back, first check whether we have no-copy view */
1387:       PetscInt  gstart,gend,start;
1388:       PetscBool contiguous,gcontiguous;
1389:       VecGetOwnershipRange(X,&gstart,&gend);
1390:       ISContiguousLocal(is,gstart,gend,&start,&contiguous);
1391:       MPI_Allreduce(&contiguous,&gcontiguous,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)is));
1392:       if (!gcontiguous) SETERRQ(PetscObjectComm((PetscObject)is),PETSC_ERR_SUP,"Unhandled case, values have been changed and need to be copied back into X");
1393:     }
1394:     VecDestroy(Y);
1395:   }
1396:   return(0);
1397: }

1401: /*@C
1402:    VecGetArray - Returns a pointer to a contiguous array that contains this
1403:    processor's portion of the vector data. For the standard PETSc
1404:    vectors, VecGetArray() returns a pointer to the local data array and
1405:    does not use any copies. If the underlying vector data is not stored
1406:    in a contiquous array this routine will copy the data to a contiquous
1407:    array and return a pointer to that. You MUST call VecRestoreArray()
1408:    when you no longer need access to the array.

1410:    Logically Collective on Vec

1412:    Input Parameter:
1413: .  x - the vector

1415:    Output Parameter:
1416: .  a - location to put pointer to the array

1418:    Fortran Note:
1419:    This routine is used differently from Fortran 77
1420: $    Vec         x
1421: $    PetscScalar x_array(1)
1422: $    PetscOffset i_x
1423: $    PetscErrorCode ierr
1424: $       call VecGetArray(x,x_array,i_x,ierr)
1425: $
1426: $   Access first local entry in vector with
1427: $      value = x_array(i_x + 1)
1428: $
1429: $      ...... other code
1430: $       call VecRestoreArray(x,x_array,i_x,ierr)
1431:    For Fortran 90 see VecGetArrayF90()

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

1436:    Level: beginner

1438:    Concepts: vector^accessing local values

1440: .seealso: VecRestoreArray(), VecGetArrayRead(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1441: @*/
1442: PetscErrorCode VecGetArray(Vec x,PetscScalar **a)
1443: {

1448:   if (x->petscnative) {
1449: #if defined(PETSC_HAVE_CUSP)
1450:     if (x->valid_GPU_array == PETSC_CUSP_GPU) {
1451:       VecCUSPCopyFromGPU(x);
1452:     }
1453: #endif
1454: #if defined(PETSC_HAVE_VIENNACL)
1455:     if (x->valid_GPU_array == PETSC_VIENNACL_GPU) {
1456:       VecViennaCLCopyFromGPU(x);
1457:     }
1458: #endif
1459:     *a = *((PetscScalar **)x->data);
1460:   } else {
1461:     (*x->ops->getarray)(x,a);
1462:   }
1463:   return(0);
1464: }

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

1471:    Not Collective

1473:    Input Parameters:
1474: .  x - the vector

1476:    Output Parameter:
1477: .  a - the array

1479:    Level: beginner

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

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

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

1490: .seealso: VecGetArray(), VecRestoreArray()
1491: @*/
1492: PetscErrorCode VecGetArrayRead(Vec x,const PetscScalar **a)
1493: {

1498:   if (x->petscnative) {
1499: #if defined(PETSC_HAVE_CUSP)
1500:     if (x->valid_GPU_array == PETSC_CUSP_GPU) {
1501:       VecCUSPCopyFromGPU(x);
1502:     }
1503: #endif
1504: #if defined(PETSC_HAVE_VIENNACL)
1505:     if (x->valid_GPU_array == PETSC_VIENNACL_GPU) {
1506:       VecViennaCLCopyFromGPU(x);
1507:     }
1508: #endif
1509:     *a = *((PetscScalar **)x->data);
1510:   } else if (x->ops->getarrayread) {
1511:     (*x->ops->getarrayread)(x,a);
1512:   } else {
1513:     (*x->ops->getarray)(x,(PetscScalar**)a);
1514:   }
1515:   return(0);
1516: }

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

1525:    Logically Collective on Vec

1527:    Input Parameter:
1528: +  x - the vectors
1529: -  n - the number of vectors

1531:    Output Parameter:
1532: .  a - location to put pointer to the array

1534:    Fortran Note:
1535:    This routine is not supported in Fortran.

1537:    Level: intermediate

1539: .seealso: VecGetArray(), VecRestoreArrays()
1540: @*/
1541: PetscErrorCode  VecGetArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1542: {
1544:   PetscInt       i;
1545:   PetscScalar    **q;

1551:   if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %D",n);
1552:   PetscMalloc1(n,&q);
1553:   for (i=0; i<n; ++i) {
1554:     VecGetArray(x[i],&q[i]);
1555:   }
1556:   *a = q;
1557:   return(0);
1558: }

1562: /*@C
1563:    VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1564:    has been called.

1566:    Logically Collective on Vec

1568:    Input Parameters:
1569: +  x - the vector
1570: .  n - the number of vectors
1571: -  a - location of pointer to arrays obtained from VecGetArrays()

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

1579:    Fortran Note:
1580:    This routine is not supported in Fortran.

1582:    Level: intermediate

1584: .seealso: VecGetArrays(), VecRestoreArray()
1585: @*/
1586: PetscErrorCode  VecRestoreArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1587: {
1589:   PetscInt       i;
1590:   PetscScalar    **q = *a;


1597:   for (i=0; i<n; ++i) {
1598:     VecRestoreArray(x[i],&q[i]);
1599:   }
1600:   PetscFree(q);
1601:   return(0);
1602: }

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

1609:    Logically Collective on Vec

1611:    Input Parameters:
1612: +  x - the vector
1613: -  a - location of pointer to array obtained from VecGetArray()

1615:    Level: beginner

1617:    Notes:
1618:    For regular PETSc vectors this routine does not involve any copies. For
1619:    any special vectors that do not store local vector data in a contiguous
1620:    array, this routine will copy the data back into the underlying
1621:    vector data structure from the array obtained with VecGetArray().

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

1627:    Fortran Note:
1628:    This routine is used differently from Fortran 77
1629: $    Vec         x
1630: $    PetscScalar x_array(1)
1631: $    PetscOffset i_x
1632: $    PetscErrorCode ierr
1633: $       call VecGetArray(x,x_array,i_x,ierr)
1634: $
1635: $   Access first local entry in vector with
1636: $      value = x_array(i_x + 1)
1637: $
1638: $      ...... other code
1639: $       call VecRestoreArray(x,x_array,i_x,ierr)

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

1645: .seealso: VecGetArray(), VecRestoreArrayRead(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1646: @*/
1647: PetscErrorCode VecRestoreArray(Vec x,PetscScalar **a)
1648: {

1653:   if (x->petscnative) {
1654: #if defined(PETSC_HAVE_CUSP)
1655:     x->valid_GPU_array = PETSC_CUSP_CPU;
1656: #endif
1657: #if defined(PETSC_HAVE_VIENNACL)
1658:     x->valid_GPU_array = PETSC_VIENNACL_CPU;
1659: #endif
1660:   } else {
1661:     (*x->ops->restorearray)(x,a);
1662:   }
1663:   if (a) *a = NULL;
1664:   PetscObjectStateIncrease((PetscObject)x);
1665:   return(0);
1666: }

1670: /*@C
1671:    VecRestoreArrayRead - Restore array obtained with VecGetArrayRead()

1673:    Not Collective

1675:    Input Parameters:
1676: +  vec - the vector
1677: -  array - the array

1679:    Level: beginner

1681: .seealso: VecGetArray(), VecRestoreArray()
1682: @*/
1683: PetscErrorCode VecRestoreArrayRead(Vec x,const PetscScalar **a)
1684: {

1689:   if (x->petscnative) {
1690: #if defined(PETSC_HAVE_VIENNACL)
1691:     x->valid_GPU_array = PETSC_VIENNACL_CPU;
1692: #endif
1693:   } else if (x->ops->restorearrayread) {
1694:     (*x->ops->restorearrayread)(x,a);
1695:   } else {
1696:     (*x->ops->restorearray)(x,(PetscScalar**)a);
1697:   }
1698:   if (a) *a = NULL;
1699:   return(0);
1700: }

1704: /*@
1705:    VecPlaceArray - Allows one to replace the array in a vector with an
1706:    array provided by the user. This is useful to avoid copying an array
1707:    into a vector.

1709:    Not Collective

1711:    Input Parameters:
1712: +  vec - the vector
1713: -  array - the array

1715:    Notes:
1716:    You can return to the original array with a call to VecResetArray()

1718:    Level: developer

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

1722: @*/
1723: PetscErrorCode  VecPlaceArray(Vec vec,const PetscScalar array[])
1724: {

1731:   if (vec->ops->placearray) {
1732:     (*vec->ops->placearray)(vec,array);
1733:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot place array in this type of vector");
1734:   PetscObjectStateIncrease((PetscObject)vec);
1735:   return(0);
1736: }

1740: /*@C
1741:    VecReplaceArray - Allows one to replace the array in a vector with an
1742:    array provided by the user. This is useful to avoid copying an array
1743:    into a vector.

1745:    Not Collective

1747:    Input Parameters:
1748: +  vec - the vector
1749: -  array - the array

1751:    Notes:
1752:    This permanently replaces the array and frees the memory associated
1753:    with the old array.

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

1758:    Not supported from Fortran

1760:    Level: developer

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

1764: @*/
1765: PetscErrorCode  VecReplaceArray(Vec vec,const PetscScalar array[])
1766: {

1772:   if (vec->ops->replacearray) {
1773:     (*vec->ops->replacearray)(vec,array);
1774:   } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot replace array in this type of vector");
1775:   PetscObjectStateIncrease((PetscObject)vec);
1776:   return(0);
1777: }

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

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

1786:     Collective on Vec

1788:     Input Parameters:
1789: +   x - a vector to mimic
1790: -   n - the number of vectors to obtain

1792:     Output Parameters:
1793: +   y - Fortran90 pointer to the array of vectors
1794: -   ierr - error code

1796:     Example of Usage:
1797: .vb
1798:     Vec x
1799:     Vec, pointer :: y(:)
1800:     ....
1801:     call VecDuplicateVecsF90(x,2,y,ierr)
1802:     call VecSet(y(2),alpha,ierr)
1803:     call VecSet(y(2),alpha,ierr)
1804:     ....
1805:     call VecDestroyVecsF90(2,y,ierr)
1806: .ve

1808:     Notes:
1809:     Not yet supported for all F90 compilers

1811:     Use VecDestroyVecsF90() to free the space.

1813:     Level: beginner

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

1817: M*/

1819: /*MC
1820:     VecRestoreArrayF90 - Restores a vector to a usable state after a call to
1821:     VecGetArrayF90().

1823:     Synopsis:
1824:     VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

1826:     Logically Collective on Vec

1828:     Input Parameters:
1829: +   x - vector
1830: -   xx_v - the Fortran90 pointer to the array

1832:     Output Parameter:
1833: .   ierr - error code

1835:     Example of Usage:
1836: .vb
1837:     PetscScalar, pointer :: xx_v(:)
1838:     ....
1839:     call VecGetArrayF90(x,xx_v,ierr)
1840:     a = xx_v(3)
1841:     call VecRestoreArrayF90(x,xx_v,ierr)
1842: .ve

1844:     Level: beginner

1846: .seealso:  VecGetArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran

1848: M*/

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

1853:     Synopsis:
1854:     VecDestroyVecsF90(PetscInt n,{Vec, pointer :: x(:)},PetscErrorCode ierr)

1856:     Collective on Vec

1858:     Input Parameters:
1859: +   n - the number of vectors previously obtained
1860: -   x - pointer to array of vector pointers

1862:     Output Parameter:
1863: .   ierr - error code

1865:     Notes:
1866:     Not yet supported for all F90 compilers

1868:     Level: beginner

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

1872: M*/

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

1880:     Synopsis:
1881:     VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)

1883:     Logically Collective on Vec

1885:     Input Parameter:
1886: .   x - vector

1888:     Output Parameters:
1889: +   xx_v - the Fortran90 pointer to the array
1890: -   ierr - error code

1892:     Example of Usage:
1893: .vb
1894:     PetscScalar, pointer :: xx_v(:)
1895:     ....
1896:     call VecGetArrayF90(x,xx_v,ierr)
1897:     a = xx_v(3)
1898:     call VecRestoreArrayF90(x,xx_v,ierr)
1899: .ve

1901:     Level: beginner

1903: .seealso:  VecRestoreArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran

1905: M*/


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

1915:    Logically Collective

1917:    Input Parameter:
1918: +  x - the vector
1919: .  m - first dimension of two dimensional array
1920: .  n - second dimension of two dimensional array
1921: .  mstart - first index you will use in first coordinate direction (often 0)
1922: -  nstart - first index in the second coordinate direction (often 0)

1924:    Output Parameter:
1925: .  a - location to put pointer to the array

1927:    Level: developer

1929:   Notes:
1930:    For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
1931:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
1932:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
1933:    the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().

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

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

1939: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1940:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1941:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1942: @*/
1943: PetscErrorCode  VecGetArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1944: {
1946:   PetscInt       i,N;
1947:   PetscScalar    *aa;

1953:   VecGetLocalSize(x,&N);
1954:   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);
1955:   VecGetArray(x,&aa);

1957:   PetscMalloc1(m,a);
1958:   for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
1959:   *a -= mstart;
1960:   return(0);
1961: }

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

1968:    Logically Collective

1970:    Input Parameters:
1971: +  x - the vector
1972: .  m - first dimension of two dimensional array
1973: .  n - second dimension of the two dimensional array
1974: .  mstart - first index you will use in first coordinate direction (often 0)
1975: .  nstart - first index in the second coordinate direction (often 0)
1976: -  a - location of pointer to array obtained from VecGetArray2d()

1978:    Level: developer

1980:    Notes:
1981:    For regular PETSc vectors this routine does not involve any copies. For
1982:    any special vectors that do not store local vector data in a contiguous
1983:    array, this routine will copy the data back into the underlying
1984:    vector data structure from the array obtained with VecGetArray().

1986:    This routine actually zeros out the a pointer.

1988: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
1989:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
1990:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1991: @*/
1992: PetscErrorCode  VecRestoreArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1993: {
1995:   void           *dummy;

2001:   dummy = (void*)(*a + mstart);
2002:   PetscFree(dummy);
2003:   VecRestoreArray(x,NULL);
2004:   return(0);
2005: }

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

2014:    Logically Collective

2016:    Input Parameter:
2017: +  x - the vector
2018: .  m - first dimension of two dimensional array
2019: -  mstart - first index you will use in first coordinate direction (often 0)

2021:    Output Parameter:
2022: .  a - location to put pointer to the array

2024:    Level: developer

2026:   Notes:
2027:    For a vector obtained from DMCreateLocalVector() mstart are likely
2028:    obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2029:    DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().

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

2033: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2034:           VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2035:           VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2036: @*/
2037: PetscErrorCode  VecGetArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2038: {
2040:   PetscInt       N;

2046:   VecGetLocalSize(x,&N);
2047:   if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2048:   VecGetArray(x,a);
2049:   *a  -= mstart;
2050:   return(0);
2051: }

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

2058:    Logically Collective

2060:    Input Parameters:
2061: +  x - the vector
2062: .  m - first dimension of two dimensional array
2063: .  mstart - first index you will use in first coordinate direction (often 0)
2064: -  a - location of pointer to array obtained from VecGetArray21()

2066:    Level: developer

2068:    Notes:
2069:    For regular PETSc vectors this routine does not involve any copies. For
2070:    any special vectors that do not store local vector data in a contiguous
2071:    array, this routine will copy the data back into the underlying
2072:    vector data structure from the array obtained with VecGetArray1d().

2074:    This routine actually zeros out the a pointer.

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

2078: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2079:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2080:           VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2081: @*/
2082: PetscErrorCode  VecRestoreArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2083: {

2089:   VecRestoreArray(x,NULL);
2090:   return(0);
2091: }


2096: /*@C
2097:    VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this
2098:    processor's portion of the vector data.  You MUST call VecRestoreArray3d()
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 three dimensional array
2106: .  n - second dimension of three dimensional array
2107: .  p - third dimension of three dimensional array
2108: .  mstart - first index you will use in first coordinate direction (often 0)
2109: .  nstart - first index in the second coordinate direction (often 0)
2110: -  pstart - first index in the third coordinate direction (often 0)

2112:    Output Parameter:
2113: .  a - location to put pointer to the array

2115:    Level: developer

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

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

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

2127: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2128:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2129:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2130: @*/
2131: PetscErrorCode  VecGetArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2132: {
2134:   PetscInt       i,N,j;
2135:   PetscScalar    *aa,**b;

2141:   VecGetLocalSize(x,&N);
2142:   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);
2143:   VecGetArray(x,&aa);

2145:   PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2146:   b    = (PetscScalar**)((*a) + m);
2147:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2148:   for (i=0; i<m; i++)
2149:     for (j=0; j<n; j++)
2150:       b[i*n+j] = aa + i*n*p + j*p - pstart;

2152:   *a -= mstart;
2153:   return(0);
2154: }

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

2161:    Logically Collective

2163:    Input Parameters:
2164: +  x - the vector
2165: .  m - first dimension of three dimensional array
2166: .  n - second dimension of the three dimensional array
2167: .  p - third dimension of the three dimensional array
2168: .  mstart - first index you will use in first coordinate direction (often 0)
2169: .  nstart - first index in the second coordinate direction (often 0)
2170: .  pstart - first index in the third coordinate direction (often 0)
2171: -  a - location of pointer to array obtained from VecGetArray3d()

2173:    Level: developer

2175:    Notes:
2176:    For regular PETSc vectors this routine does not involve any copies. For
2177:    any special vectors that do not store local vector data in a contiguous
2178:    array, this routine will copy the data back into the underlying
2179:    vector data structure from the array obtained with VecGetArray().

2181:    This routine actually zeros out the a pointer.

2183: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2184:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2185:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2186: @*/
2187: PetscErrorCode  VecRestoreArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2188: {
2190:   void           *dummy;

2196:   dummy = (void*)(*a + mstart);
2197:   PetscFree(dummy);
2198:   VecRestoreArray(x,NULL);
2199:   return(0);
2200: }

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

2209:    Logically Collective

2211:    Input Parameter:
2212: +  x - the vector
2213: .  m - first dimension of four dimensional array
2214: .  n - second dimension of four dimensional array
2215: .  p - third dimension of four dimensional array
2216: .  q - fourth dimension of four dimensional array
2217: .  mstart - first index you will use in first coordinate direction (often 0)
2218: .  nstart - first index in the second coordinate direction (often 0)
2219: .  pstart - first index in the third coordinate direction (often 0)
2220: -  qstart - first index in the fourth coordinate direction (often 0)

2222:    Output Parameter:
2223: .  a - location to put pointer to the array

2225:    Level: beginner

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

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

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

2237: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2238:           VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2239:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2240: @*/
2241: PetscErrorCode  VecGetArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2242: {
2244:   PetscInt       i,N,j,k;
2245:   PetscScalar    *aa,***b,**c;

2251:   VecGetLocalSize(x,&N);
2252:   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);
2253:   VecGetArray(x,&aa);

2255:   PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2256:   b    = (PetscScalar***)((*a) + m);
2257:   c    = (PetscScalar**)(b + m*n);
2258:   for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2259:   for (i=0; i<m; i++)
2260:     for (j=0; j<n; j++)
2261:       b[i*n+j] = c + i*n*p + j*p - pstart;
2262:   for (i=0; i<m; i++)
2263:     for (j=0; j<n; j++)
2264:       for (k=0; k<p; k++)
2265:         c[i*n*p+j*p+k] = aa + i*n*p*q + j*p*q + k*q - qstart;
2266:   *a -= mstart;
2267:   return(0);
2268: }

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

2275:    Logically Collective

2277:    Input Parameters:
2278: +  x - the vector
2279: .  m - first dimension of four dimensional array
2280: .  n - second dimension of the four dimensional array
2281: .  p - third dimension of the four dimensional array
2282: .  q - fourth dimension of the four dimensional array
2283: .  mstart - first index you will use in first coordinate direction (often 0)
2284: .  nstart - first index in the second coordinate direction (often 0)
2285: .  pstart - first index in the third coordinate direction (often 0)
2286: .  qstart - first index in the fourth coordinate direction (often 0)
2287: -  a - location of pointer to array obtained from VecGetArray4d()

2289:    Level: beginner

2291:    Notes:
2292:    For regular PETSc vectors this routine does not involve any copies. For
2293:    any special vectors that do not store local vector data in a contiguous
2294:    array, this routine will copy the data back into the underlying
2295:    vector data structure from the array obtained with VecGetArray().

2297:    This routine actually zeros out the a pointer.

2299: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2300:           VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2301:           VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2302: @*/
2303: PetscErrorCode  VecRestoreArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2304: {
2306:   void           *dummy;

2312:   dummy = (void*)(*a + mstart);
2313:   PetscFree(dummy);
2314:   VecRestoreArray(x,NULL);
2315:   return(0);
2316: }