Actual source code: rvector.c
petsc-3.5.4 2015-05-23
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 (x->map->N == 0) {
582: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],0.0l);
583: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],0.0);
584: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],0.0);
585: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],0.0);
586: } else if (val > PETSC_MAX_REAL/x->map->N) {
587: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
588: } else {
589: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],x->map->N * val);
590: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
591: val = PetscSqrtReal((PetscReal)x->map->N) * val;
592: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],val);
593: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],val);
594: }
595: return(0);
596: }
601: /*@
602: VecAXPY - Computes y = alpha x + y.
604: Logically Collective on Vec
606: Input Parameters:
607: + alpha - the scalar
608: - x, y - the vectors
610: Output Parameter:
611: . y - output vector
613: Level: intermediate
615: Notes: x and y MUST be different vectors
617: Concepts: vector^BLAS
618: Concepts: BLAS
620: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
621: @*/
622: PetscErrorCode VecAXPY(Vec y,PetscScalar alpha,Vec x)
623: {
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->axpy)(y,alpha,x);
638: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
639: PetscObjectStateIncrease((PetscObject)y);
640: return(0);
641: }
645: /*@
646: VecAXPBY - Computes y = alpha x + beta y.
648: Logically Collective on Vec
650: Input Parameters:
651: + alpha,beta - the scalars
652: - x, y - the vectors
654: Output Parameter:
655: . y - output vector
657: Level: intermediate
659: Notes: x and y MUST be different vectors
661: Concepts: BLAS
662: Concepts: vector^BLAS
664: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
665: @*/
666: PetscErrorCode VecAXPBY(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
667: {
677: if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y cannot be the same vector");
681: PetscLogEventBegin(VEC_AXPY,x,y,0,0);
682: (*y->ops->axpby)(y,alpha,beta,x);
683: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
684: PetscObjectStateIncrease((PetscObject)y);
685: return(0);
686: }
690: /*@
691: VecAXPBYPCZ - Computes z = alpha x + beta y + gamma z
693: Logically Collective on Vec
695: Input Parameters:
696: + alpha,beta, gamma - the scalars
697: - x, y, z - the vectors
699: Output Parameter:
700: . z - output vector
702: Level: intermediate
704: Notes: x, y and z must be different vectors
706: Developer Note: alpha = 1 or gamma = 1 or gamma = 0.0 are handled as special cases
708: Concepts: BLAS
709: Concepts: vector^BLAS
711: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
712: @*/
713: PetscErrorCode VecAXPBYPCZ(Vec z,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec x,Vec y)
714: {
728: if (x == y || x == z) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");
729: if (y == z) SETERRQ(PetscObjectComm((PetscObject)y),PETSC_ERR_ARG_IDN,"x, y, and z must be different vectors");
734: PetscLogEventBegin(VEC_AXPBYPCZ,x,y,z,0);
735: (*y->ops->axpbypcz)(z,alpha,beta,gamma,x,y);
736: PetscLogEventEnd(VEC_AXPBYPCZ,x,y,z,0);
737: PetscObjectStateIncrease((PetscObject)z);
738: return(0);
739: }
743: /*@
744: VecAYPX - Computes y = x + alpha y.
746: Logically Collective on Vec
748: Input Parameters:
749: + alpha - the scalar
750: - x, y - the vectors
752: Output Parameter:
753: . y - output vector
755: Level: intermediate
757: Notes: x and y MUST be different vectors
759: Concepts: vector^BLAS
760: Concepts: BLAS
762: .seealso: VecAXPY(), VecWAXPY()
763: @*/
764: PetscErrorCode VecAYPX(Vec y,PetscScalar alpha,Vec x)
765: {
773: if (x == y) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
776: PetscLogEventBegin(VEC_AYPX,x,y,0,0);
777: (*y->ops->aypx)(y,alpha,x);
778: PetscLogEventEnd(VEC_AYPX,x,y,0,0);
779: PetscObjectStateIncrease((PetscObject)y);
780: return(0);
781: }
786: /*@
787: VecWAXPY - Computes w = alpha x + y.
789: Logically Collective on Vec
791: Input Parameters:
792: + alpha - the scalar
793: - x, y - the vectors
795: Output Parameter:
796: . w - the result
798: Level: intermediate
800: Notes: w cannot be either x or y, but x and y can be the same
802: Concepts: vector^BLAS
803: Concepts: BLAS
805: .seealso: VecAXPY(), VecAYPX(), VecAXPBY()
806: @*/
807: PetscErrorCode VecWAXPY(Vec w,PetscScalar alpha,Vec x,Vec y)
808: {
822: if (w == y) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector y, suggest VecAXPY()");
823: if (w == x) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Result vector w cannot be same as input vector x, suggest VecAYPX()");
826: PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
827: (*w->ops->waxpy)(w,alpha,x,y);
828: PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
829: PetscObjectStateIncrease((PetscObject)w);
830: return(0);
831: }
836: /*@
837: VecSetValues - Inserts or adds values into certain locations of a vector.
839: Not Collective
841: Input Parameters:
842: + x - vector to insert in
843: . ni - number of elements to add
844: . ix - indices where to add
845: . y - array of values
846: - iora - either INSERT_VALUES or ADD_VALUES, where
847: ADD_VALUES adds values to any existing entries, and
848: INSERT_VALUES replaces existing entries with new values
850: Notes:
851: VecSetValues() sets x[ix[i]] = y[i], for i=0,...,ni-1.
853: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
854: options cannot be mixed without intervening calls to the assembly
855: routines.
857: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
858: MUST be called after all calls to VecSetValues() have been completed.
860: VecSetValues() uses 0-based indices in Fortran as well as in C.
862: If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
863: negative indices may be passed in ix. These rows are
864: simply ignored. This allows easily inserting element load matrices
865: with homogeneous Dirchlet boundary conditions that you don't want represented
866: in the vector.
868: Level: beginner
870: Concepts: vector^setting values
872: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
873: VecSetValue(), VecSetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecGetValues()
874: @*/
875: PetscErrorCode VecSetValues(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
876: {
881: if (!ni) return(0);
885: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
886: (*x->ops->setvalues)(x,ni,ix,y,iora);
887: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
888: PetscObjectStateIncrease((PetscObject)x);
889: return(0);
890: }
894: /*@
895: VecGetValues - Gets values from certain locations of a vector. Currently
896: can only get values on the same processor
898: Not Collective
900: Input Parameters:
901: + x - vector to get values from
902: . ni - number of elements to get
903: - ix - indices where to get them from (in global 1d numbering)
905: Output Parameter:
906: . y - array of values
908: Notes:
909: The user provides the allocated array y; it is NOT allocated in this routine
911: VecGetValues() gets y[i] = x[ix[i]], for i=0,...,ni-1.
913: VecAssemblyBegin() and VecAssemblyEnd() MUST be called before calling this
915: VecGetValues() uses 0-based indices in Fortran as well as in C.
917: If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
918: negative indices may be passed in ix. These rows are
919: simply ignored.
921: Level: beginner
923: Concepts: vector^getting values
925: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecGetValuesLocal(),
926: VecGetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecSetValues()
927: @*/
928: PetscErrorCode VecGetValues(Vec x,PetscInt ni,const PetscInt ix[],PetscScalar y[])
929: {
934: if (!ni) return(0);
938: (*x->ops->getvalues)(x,ni,ix,y);
939: return(0);
940: }
944: /*@
945: VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector.
947: Not Collective
949: Input Parameters:
950: + x - vector to insert in
951: . ni - number of blocks to add
952: . ix - indices where to add in block count, rather than element count
953: . y - array of values
954: - iora - either INSERT_VALUES or ADD_VALUES, where
955: ADD_VALUES adds values to any existing entries, and
956: INSERT_VALUES replaces existing entries with new values
958: Notes:
959: VecSetValuesBlocked() sets x[bs*ix[i]+j] = y[bs*i+j],
960: for j=0,...,bs, for i=0,...,ni-1. where bs was set with VecSetBlockSize().
962: Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
963: options cannot be mixed without intervening calls to the assembly
964: routines.
966: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
967: MUST be called after all calls to VecSetValuesBlocked() have been completed.
969: VecSetValuesBlocked() uses 0-based indices in Fortran as well as in C.
971: Negative indices may be passed in ix, these rows are
972: simply ignored. This allows easily inserting element load matrices
973: with homogeneous Dirchlet boundary conditions that you don't want represented
974: in the vector.
976: Level: intermediate
978: Concepts: vector^setting values blocked
980: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
981: VecSetValues()
982: @*/
983: PetscErrorCode VecSetValuesBlocked(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
984: {
992: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
993: (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
994: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
995: PetscObjectStateIncrease((PetscObject)x);
996: return(0);
997: }
1002: /*@
1003: VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
1004: using a local ordering of the nodes.
1006: Not Collective
1008: Input Parameters:
1009: + x - vector to insert in
1010: . ni - number of elements to add
1011: . ix - indices where to add
1012: . y - array of values
1013: - iora - either INSERT_VALUES or ADD_VALUES, where
1014: ADD_VALUES adds values to any existing entries, and
1015: INSERT_VALUES replaces existing entries with new values
1017: Level: intermediate
1019: Notes:
1020: VecSetValuesLocal() sets x[ix[i]] = y[i], for i=0,...,ni-1.
1022: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
1023: options cannot be mixed without intervening calls to the assembly
1024: routines.
1026: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1027: MUST be called after all calls to VecSetValuesLocal() have been completed.
1029: VecSetValuesLocal() uses 0-based indices in Fortran as well as in C.
1031: Concepts: vector^setting values with local numbering
1033: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
1034: VecSetValuesBlockedLocal()
1035: @*/
1036: PetscErrorCode VecSetValuesLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1037: {
1039: PetscInt lixp[128],*lix = lixp;
1043: if (!ni) return(0);
1048: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1049: if (!x->ops->setvalueslocal) {
1050: if (!x->map->mapping) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
1051: if (ni > 128) {
1052: PetscMalloc1(ni,&lix);
1053: }
1054: ISLocalToGlobalMappingApply(x->map->mapping,ni,(PetscInt*)ix,lix);
1055: (*x->ops->setvalues)(x,ni,lix,y,iora);
1056: if (ni > 128) {
1057: PetscFree(lix);
1058: }
1059: } else {
1060: (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
1061: }
1062: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1063: PetscObjectStateIncrease((PetscObject)x);
1064: return(0);
1065: }
1069: /*@
1070: VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1071: using a local ordering of the nodes.
1073: Not Collective
1075: Input Parameters:
1076: + x - vector to insert in
1077: . ni - number of blocks to add
1078: . ix - indices where to add in block count, not element count
1079: . y - array of values
1080: - iora - either INSERT_VALUES or ADD_VALUES, where
1081: ADD_VALUES adds values to any existing entries, and
1082: INSERT_VALUES replaces existing entries with new values
1084: Level: intermediate
1086: Notes:
1087: VecSetValuesBlockedLocal() sets x[bs*ix[i]+j] = y[bs*i+j],
1088: for j=0,..bs-1, for i=0,...,ni-1, where bs has been set with VecSetBlockSize().
1090: Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1091: options cannot be mixed without intervening calls to the assembly
1092: routines.
1094: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1095: MUST be called after all calls to VecSetValuesBlockedLocal() have been completed.
1097: VecSetValuesBlockedLocal() uses 0-based indices in Fortran as well as in C.
1100: Concepts: vector^setting values blocked with local numbering
1102: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1103: VecSetLocalToGlobalMapping()
1104: @*/
1105: PetscErrorCode VecSetValuesBlockedLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1106: {
1108: PetscInt lixp[128],*lix = lixp;
1115: if (ni > 128) {
1116: PetscMalloc1(ni,&lix);
1117: }
1119: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1120: ISLocalToGlobalMappingApplyBlock(x->map->mapping,ni,(PetscInt*)ix,lix);
1121: (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1122: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1123: if (ni > 128) {
1124: PetscFree(lix);
1125: }
1126: PetscObjectStateIncrease((PetscObject)x);
1127: return(0);
1128: }
1132: /*@
1133: VecMTDot - Computes indefinite vector multiple dot products.
1134: That is, it does NOT use the complex conjugate.
1136: Collective on Vec
1138: Input Parameters:
1139: + x - one vector
1140: . nv - number of vectors
1141: - y - array of vectors. Note that vectors are pointers
1143: Output Parameter:
1144: . val - array of the dot products
1146: Notes for Users of Complex Numbers:
1147: For complex vectors, VecMTDot() computes the indefinite form
1148: $ val = (x,y) = y^T x,
1149: where y^T denotes the transpose of y.
1151: Use VecMDot() for the inner product
1152: $ val = (x,y) = y^H x,
1153: where y^H denotes the conjugate transpose of y.
1155: Level: intermediate
1157: Concepts: inner product^multiple
1158: Concepts: vector^multiple inner products
1160: .seealso: VecMDot(), VecTDot()
1161: @*/
1162: PetscErrorCode VecMTDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1163: {
1176: PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1177: (*x->ops->mtdot)(x,nv,y,val);
1178: PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1179: return(0);
1180: }
1184: /*@
1185: VecMDot - Computes vector multiple dot products.
1187: Collective on Vec
1189: Input Parameters:
1190: + x - one vector
1191: . nv - number of vectors
1192: - y - array of vectors.
1194: Output Parameter:
1195: . val - array of the dot products (does not allocate the array)
1197: Notes for Users of Complex Numbers:
1198: For complex vectors, VecMDot() computes
1199: $ val = (x,y) = y^H x,
1200: where y^H denotes the conjugate transpose of y.
1202: Use VecMTDot() for the indefinite form
1203: $ val = (x,y) = y^T x,
1204: where y^T denotes the transpose of y.
1206: Level: intermediate
1208: Concepts: inner product^multiple
1209: Concepts: vector^multiple inner products
1211: .seealso: VecMTDot(), VecDot()
1212: @*/
1213: PetscErrorCode VecMDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1214: {
1219: if (!nv) return(0);
1220: if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1229: PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,PetscObjectComm((PetscObject)x));
1230: (*x->ops->mdot)(x,nv,y,val);
1231: PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,PetscObjectComm((PetscObject)x));
1232: return(0);
1233: }
1237: /*@
1238: VecMAXPY - Computes y = y + sum alpha[j] x[j]
1240: Logically Collective on Vec
1242: Input Parameters:
1243: + nv - number of scalars and x-vectors
1244: . alpha - array of scalars
1245: . y - one vector
1246: - x - array of vectors
1248: Level: intermediate
1250: Notes: y cannot be any of the x vectors
1252: Concepts: BLAS
1254: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1255: @*/
1256: PetscErrorCode VecMAXPY(Vec y,PetscInt nv,const PetscScalar alpha[],Vec x[])
1257: {
1259: PetscInt i;
1263: if (!nv) return(0);
1264: if (nv < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1274: PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1275: (*y->ops->maxpy)(y,nv,alpha,x);
1276: PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1277: PetscObjectStateIncrease((PetscObject)y);
1278: return(0);
1279: }
1283: /*@
1284: VecGetSubVector - Gets a vector representing part of another vector
1286: Collective on IS (and Vec if nonlocal entries are needed)
1288: Input Arguments:
1289: + X - vector from which to extract a subvector
1290: - is - index set representing portion of X to extract
1292: Output Arguments:
1293: . Y - subvector corresponding to is
1295: Level: advanced
1297: Notes:
1298: The subvector Y should be returned with VecRestoreSubVector().
1300: This function may return a subvector without making a copy, therefore it is not safe to use the original vector while
1301: modifying the subvector. Other non-overlapping subvectors can still be obtained from X using this function.
1303: .seealso: MatGetSubMatrix()
1304: @*/
1305: PetscErrorCode VecGetSubVector(Vec X,IS is,Vec *Y)
1306: {
1307: PetscErrorCode ierr;
1308: Vec Z;
1314: if (X->ops->getsubvector) {
1315: (*X->ops->getsubvector)(X,is,&Z);
1316: } else { /* Default implementation currently does no caching */
1317: PetscInt gstart,gend,start;
1318: PetscBool contiguous,gcontiguous;
1319: VecGetOwnershipRange(X,&gstart,&gend);
1320: ISContiguousLocal(is,gstart,gend,&start,&contiguous);
1321: MPI_Allreduce(&contiguous,&gcontiguous,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)is));
1322: if (gcontiguous) { /* We can do a no-copy implementation */
1323: PetscInt n,N,bs;
1324: PetscScalar *x;
1325: PetscMPIInt size;
1326: ISGetLocalSize(is,&n);
1327: VecGetArray(X,&x);
1328: VecGetBlockSize(X,&bs);
1329: if (n%bs) bs = 1;
1330: MPI_Comm_size(PetscObjectComm((PetscObject)X),&size);
1331: if (size == 1) {
1332: VecCreateSeqWithArray(PetscObjectComm((PetscObject)X),bs,n,x+start,&Z);
1333: } else {
1334: ISGetSize(is,&N);
1335: VecCreateMPIWithArray(PetscObjectComm((PetscObject)X),bs,n,N,x+start,&Z);
1336: }
1337: VecRestoreArray(X,&x);
1338: } else { /* Have to create a scatter and do a copy */
1339: VecScatter scatter;
1340: PetscInt n,N;
1341: ISGetLocalSize(is,&n);
1342: ISGetSize(is,&N);
1343: VecCreate(PetscObjectComm((PetscObject)is),&Z);
1344: VecSetSizes(Z,n,N);
1345: VecSetType(Z,((PetscObject)X)->type_name);
1346: VecScatterCreate(X,is,Z,NULL,&scatter);
1347: VecScatterBegin(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1348: VecScatterEnd(scatter,X,Z,INSERT_VALUES,SCATTER_FORWARD);
1349: PetscObjectCompose((PetscObject)Z,"VecGetSubVector_Scatter",(PetscObject)scatter);
1350: VecScatterDestroy(&scatter);
1351: }
1352: }
1353: /* Record the state when the subvector was gotten so we know whether its values need to be put back */
1354: if (VecGetSubVectorSavedStateId < 0) {PetscObjectComposedDataRegister(&VecGetSubVectorSavedStateId);}
1355: PetscObjectComposedDataSetInt((PetscObject)Z,VecGetSubVectorSavedStateId,1);
1356: *Y = Z;
1357: return(0);
1358: }
1362: /*@
1363: VecRestoreSubVector - Restores a subvector extracted using VecGetSubVector()
1365: Collective on IS (and Vec if nonlocal entries need to be written)
1367: Input Arguments:
1368: + X - vector from which subvector was obtained
1369: . is - index set representing the subset of X
1370: - Y - subvector being restored
1372: Level: advanced
1374: .seealso: VecGetSubVector()
1375: @*/
1376: PetscErrorCode VecRestoreSubVector(Vec X,IS is,Vec *Y)
1377: {
1385: if (X->ops->restoresubvector) {
1386: (*X->ops->restoresubvector)(X,is,Y);
1387: } else {
1388: PETSC_UNUSED PetscObjectState dummystate = 0;
1389: PetscBool valid;
1390: PetscObjectComposedDataGetInt((PetscObject)*Y,VecGetSubVectorSavedStateId,dummystate,valid);
1391: if (!valid) {
1392: VecScatter scatter;
1394: PetscObjectQuery((PetscObject)*Y,"VecGetSubVector_Scatter",(PetscObject*)&scatter);
1395: if (scatter) {
1396: VecScatterBegin(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1397: VecScatterEnd(scatter,*Y,X,INSERT_VALUES,SCATTER_REVERSE);
1398: }
1399: }
1400: VecDestroy(Y);
1401: }
1402: return(0);
1403: }
1407: /*@C
1408: VecGetArray - Returns a pointer to a contiguous array that contains this
1409: processor's portion of the vector data. For the standard PETSc
1410: vectors, VecGetArray() returns a pointer to the local data array and
1411: does not use any copies. If the underlying vector data is not stored
1412: in a contiquous array this routine will copy the data to a contiquous
1413: array and return a pointer to that. You MUST call VecRestoreArray()
1414: when you no longer need access to the array.
1416: Logically Collective on Vec
1418: Input Parameter:
1419: . x - the vector
1421: Output Parameter:
1422: . a - location to put pointer to the array
1424: Fortran Note:
1425: This routine is used differently from Fortran 77
1426: $ Vec x
1427: $ PetscScalar x_array(1)
1428: $ PetscOffset i_x
1429: $ PetscErrorCode ierr
1430: $ call VecGetArray(x,x_array,i_x,ierr)
1431: $
1432: $ Access first local entry in vector with
1433: $ value = x_array(i_x + 1)
1434: $
1435: $ ...... other code
1436: $ call VecRestoreArray(x,x_array,i_x,ierr)
1437: For Fortran 90 see VecGetArrayF90()
1439: See the Fortran chapter of the users manual and
1440: petsc/src/snes/examples/tutorials/ex5f.F for details.
1442: Level: beginner
1444: Concepts: vector^accessing local values
1446: .seealso: VecRestoreArray(), VecGetArrayRead(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1447: @*/
1448: PetscErrorCode VecGetArray(Vec x,PetscScalar **a)
1449: {
1454: if (x->petscnative) {
1455: #if defined(PETSC_HAVE_CUSP)
1456: if (x->valid_GPU_array == PETSC_CUSP_GPU) {
1457: VecCUSPCopyFromGPU(x);
1458: }
1459: #endif
1460: #if defined(PETSC_HAVE_VIENNACL)
1461: if (x->valid_GPU_array == PETSC_VIENNACL_GPU) {
1462: VecViennaCLCopyFromGPU(x);
1463: }
1464: #endif
1465: *a = *((PetscScalar **)x->data);
1466: } else {
1467: (*x->ops->getarray)(x,a);
1468: }
1469: return(0);
1470: }
1474: /*@C
1475: VecGetArrayRead - Get read-only pointer to contiguous array containing this processor's portion of the vector data.
1477: Not Collective
1479: Input Parameters:
1480: . x - the vector
1482: Output Parameter:
1483: . a - the array
1485: Level: beginner
1487: Notes:
1488: The array must be returned using a matching call to VecRestoreArrayRead().
1490: Unlike VecGetArray(), this routine is not collective and preserves cached information like vector norms.
1492: Standard PETSc vectors use contiguous storage so that this routine does not perform a copy. Other vector
1493: implementations may require a copy, but must such implementations should cache the contiguous representation so that
1494: only one copy is performed when this routine is called multiple times in sequence.
1496: .seealso: VecGetArray(), VecRestoreArray()
1497: @*/
1498: PetscErrorCode VecGetArrayRead(Vec x,const PetscScalar **a)
1499: {
1504: if (x->petscnative) {
1505: #if defined(PETSC_HAVE_CUSP)
1506: if (x->valid_GPU_array == PETSC_CUSP_GPU) {
1507: VecCUSPCopyFromGPU(x);
1508: }
1509: #endif
1510: #if defined(PETSC_HAVE_VIENNACL)
1511: if (x->valid_GPU_array == PETSC_VIENNACL_GPU) {
1512: VecViennaCLCopyFromGPU(x);
1513: }
1514: #endif
1515: *a = *((PetscScalar **)x->data);
1516: } else if (x->ops->getarrayread) {
1517: (*x->ops->getarrayread)(x,a);
1518: } else {
1519: (*x->ops->getarray)(x,(PetscScalar**)a);
1520: }
1521: return(0);
1522: }
1526: /*@C
1527: VecGetArrays - Returns a pointer to the arrays in a set of vectors
1528: that were created by a call to VecDuplicateVecs(). You MUST call
1529: VecRestoreArrays() when you no longer need access to the array.
1531: Logically Collective on Vec
1533: Input Parameter:
1534: + x - the vectors
1535: - n - the number of vectors
1537: Output Parameter:
1538: . a - location to put pointer to the array
1540: Fortran Note:
1541: This routine is not supported in Fortran.
1543: Level: intermediate
1545: .seealso: VecGetArray(), VecRestoreArrays()
1546: @*/
1547: PetscErrorCode VecGetArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1548: {
1550: PetscInt i;
1551: PetscScalar **q;
1557: if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %D",n);
1558: PetscMalloc1(n,&q);
1559: for (i=0; i<n; ++i) {
1560: VecGetArray(x[i],&q[i]);
1561: }
1562: *a = q;
1563: return(0);
1564: }
1568: /*@C
1569: VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1570: has been called.
1572: Logically Collective on Vec
1574: Input Parameters:
1575: + x - the vector
1576: . n - the number of vectors
1577: - a - location of pointer to arrays obtained from VecGetArrays()
1579: Notes:
1580: For regular PETSc vectors this routine does not involve any copies. For
1581: any special vectors that do not store local vector data in a contiguous
1582: array, this routine will copy the data back into the underlying
1583: vector data structure from the arrays obtained with VecGetArrays().
1585: Fortran Note:
1586: This routine is not supported in Fortran.
1588: Level: intermediate
1590: .seealso: VecGetArrays(), VecRestoreArray()
1591: @*/
1592: PetscErrorCode VecRestoreArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1593: {
1595: PetscInt i;
1596: PetscScalar **q = *a;
1603: for (i=0; i<n; ++i) {
1604: VecRestoreArray(x[i],&q[i]);
1605: }
1606: PetscFree(q);
1607: return(0);
1608: }
1612: /*@C
1613: VecRestoreArray - Restores a vector after VecGetArray() has been called.
1615: Logically Collective on Vec
1617: Input Parameters:
1618: + x - the vector
1619: - a - location of pointer to array obtained from VecGetArray()
1621: Level: beginner
1623: Notes:
1624: For regular PETSc vectors this routine does not involve any copies. For
1625: any special vectors that do not store local vector data in a contiguous
1626: array, this routine will copy the data back into the underlying
1627: vector data structure from the array obtained with VecGetArray().
1629: This routine actually zeros out the a pointer. This is to prevent accidental
1630: us of the array after it has been restored. If you pass null for a it will
1631: not zero the array pointer a.
1633: Fortran Note:
1634: This routine is used differently from Fortran 77
1635: $ Vec x
1636: $ PetscScalar x_array(1)
1637: $ PetscOffset i_x
1638: $ PetscErrorCode ierr
1639: $ call VecGetArray(x,x_array,i_x,ierr)
1640: $
1641: $ Access first local entry in vector with
1642: $ value = x_array(i_x + 1)
1643: $
1644: $ ...... other code
1645: $ call VecRestoreArray(x,x_array,i_x,ierr)
1647: See the Fortran chapter of the users manual and
1648: petsc/src/snes/examples/tutorials/ex5f.F for details.
1649: For Fortran 90 see VecRestoreArrayF90()
1651: .seealso: VecGetArray(), VecRestoreArrayRead(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1652: @*/
1653: PetscErrorCode VecRestoreArray(Vec x,PetscScalar **a)
1654: {
1659: if (x->petscnative) {
1660: #if defined(PETSC_HAVE_CUSP)
1661: x->valid_GPU_array = PETSC_CUSP_CPU;
1662: #endif
1663: #if defined(PETSC_HAVE_VIENNACL)
1664: x->valid_GPU_array = PETSC_VIENNACL_CPU;
1665: #endif
1666: } else {
1667: (*x->ops->restorearray)(x,a);
1668: }
1669: if (a) *a = NULL;
1670: PetscObjectStateIncrease((PetscObject)x);
1671: return(0);
1672: }
1676: /*@C
1677: VecRestoreArrayRead - Restore array obtained with VecGetArrayRead()
1679: Not Collective
1681: Input Parameters:
1682: + vec - the vector
1683: - array - the array
1685: Level: beginner
1687: .seealso: VecGetArray(), VecRestoreArray()
1688: @*/
1689: PetscErrorCode VecRestoreArrayRead(Vec x,const PetscScalar **a)
1690: {
1695: if (x->petscnative) {
1696: #if defined(PETSC_HAVE_VIENNACL)
1697: x->valid_GPU_array = PETSC_VIENNACL_CPU;
1698: #endif
1699: } else if (x->ops->restorearrayread) {
1700: (*x->ops->restorearrayread)(x,a);
1701: } else {
1702: (*x->ops->restorearray)(x,(PetscScalar**)a);
1703: }
1704: if (a) *a = NULL;
1705: return(0);
1706: }
1710: /*@
1711: VecPlaceArray - Allows one to replace the array in a vector with an
1712: array provided by the user. This is useful to avoid copying an array
1713: into a vector.
1715: Not Collective
1717: Input Parameters:
1718: + vec - the vector
1719: - array - the array
1721: Notes:
1722: You can return to the original array with a call to VecResetArray()
1724: Level: developer
1726: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()
1728: @*/
1729: PetscErrorCode VecPlaceArray(Vec vec,const PetscScalar array[])
1730: {
1737: if (vec->ops->placearray) {
1738: (*vec->ops->placearray)(vec,array);
1739: } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot place array in this type of vector");
1740: PetscObjectStateIncrease((PetscObject)vec);
1741: return(0);
1742: }
1746: /*@C
1747: VecReplaceArray - Allows one to replace the array in a vector with an
1748: array provided by the user. This is useful to avoid copying an array
1749: into a vector.
1751: Not Collective
1753: Input Parameters:
1754: + vec - the vector
1755: - array - the array
1757: Notes:
1758: This permanently replaces the array and frees the memory associated
1759: with the old array.
1761: The memory passed in MUST be obtained with PetscMalloc() and CANNOT be
1762: freed by the user. It will be freed when the vector is destroy.
1764: Not supported from Fortran
1766: Level: developer
1768: .seealso: VecGetArray(), VecRestoreArray(), VecPlaceArray(), VecResetArray()
1770: @*/
1771: PetscErrorCode VecReplaceArray(Vec vec,const PetscScalar array[])
1772: {
1778: if (vec->ops->replacearray) {
1779: (*vec->ops->replacearray)(vec,array);
1780: } else SETERRQ(PetscObjectComm((PetscObject)vec),PETSC_ERR_SUP,"Cannot replace array in this type of vector");
1781: PetscObjectStateIncrease((PetscObject)vec);
1782: return(0);
1783: }
1785: /*MC
1786: VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
1787: and makes them accessible via a Fortran90 pointer.
1789: Synopsis:
1790: VecDuplicateVecsF90(Vec x,PetscInt n,{Vec, pointer :: y(:)},integer ierr)
1792: Collective on Vec
1794: Input Parameters:
1795: + x - a vector to mimic
1796: - n - the number of vectors to obtain
1798: Output Parameters:
1799: + y - Fortran90 pointer to the array of vectors
1800: - ierr - error code
1802: Example of Usage:
1803: .vb
1804: Vec x
1805: Vec, pointer :: y(:)
1806: ....
1807: call VecDuplicateVecsF90(x,2,y,ierr)
1808: call VecSet(y(2),alpha,ierr)
1809: call VecSet(y(2),alpha,ierr)
1810: ....
1811: call VecDestroyVecsF90(2,y,ierr)
1812: .ve
1814: Notes:
1815: Not yet supported for all F90 compilers
1817: Use VecDestroyVecsF90() to free the space.
1819: Level: beginner
1821: .seealso: VecDestroyVecsF90(), VecDuplicateVecs()
1823: M*/
1825: /*MC
1826: VecRestoreArrayF90 - Restores a vector to a usable state after a call to
1827: VecGetArrayF90().
1829: Synopsis:
1830: VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
1832: Logically Collective on Vec
1834: Input Parameters:
1835: + x - vector
1836: - xx_v - the Fortran90 pointer to the array
1838: Output Parameter:
1839: . ierr - error code
1841: Example of Usage:
1842: .vb
1843: PetscScalar, pointer :: xx_v(:)
1844: ....
1845: call VecGetArrayF90(x,xx_v,ierr)
1846: a = xx_v(3)
1847: call VecRestoreArrayF90(x,xx_v,ierr)
1848: .ve
1850: Level: beginner
1852: .seealso: VecGetArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran
1854: M*/
1856: /*MC
1857: VecDestroyVecsF90 - Frees a block of vectors obtained with VecDuplicateVecsF90().
1859: Synopsis:
1860: VecDestroyVecsF90(PetscInt n,{Vec, pointer :: x(:)},PetscErrorCode ierr)
1862: Collective on Vec
1864: Input Parameters:
1865: + n - the number of vectors previously obtained
1866: - x - pointer to array of vector pointers
1868: Output Parameter:
1869: . ierr - error code
1871: Notes:
1872: Not yet supported for all F90 compilers
1874: Level: beginner
1876: .seealso: VecDestroyVecs(), VecDuplicateVecsF90()
1878: M*/
1880: /*MC
1881: VecGetArrayF90 - Accesses a vector array from Fortran90. For default PETSc
1882: vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
1883: this routine is implementation dependent. You MUST call VecRestoreArrayF90()
1884: when you no longer need access to the array.
1886: Synopsis:
1887: VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
1889: Logically Collective on Vec
1891: Input Parameter:
1892: . x - vector
1894: Output Parameters:
1895: + xx_v - the Fortran90 pointer to the array
1896: - ierr - error code
1898: Example of Usage:
1899: .vb
1900: PetscScalar, pointer :: xx_v(:)
1901: ....
1902: call VecGetArrayF90(x,xx_v,ierr)
1903: a = xx_v(3)
1904: call VecRestoreArrayF90(x,xx_v,ierr)
1905: .ve
1907: Level: beginner
1909: .seealso: VecRestoreArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran
1911: M*/
1916: /*@C
1917: VecGetArray2d - Returns a pointer to a 2d contiguous array that contains this
1918: processor's portion of the vector data. You MUST call VecRestoreArray2d()
1919: when you no longer need access to the array.
1921: Logically Collective
1923: Input Parameter:
1924: + x - the vector
1925: . m - first dimension of two dimensional array
1926: . n - second dimension of two dimensional array
1927: . mstart - first index you will use in first coordinate direction (often 0)
1928: - nstart - first index in the second coordinate direction (often 0)
1930: Output Parameter:
1931: . a - location to put pointer to the array
1933: Level: developer
1935: Notes:
1936: For a vector obtained from DMCreateLocalVector() mstart and nstart are likely
1937: obtained from the corner indices obtained from DMDAGetGhostCorners() while for
1938: DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
1939: the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().
1941: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
1943: Concepts: vector^accessing local values as 2d array
1945: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1946: VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1947: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1948: @*/
1949: PetscErrorCode VecGetArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1950: {
1952: PetscInt i,N;
1953: PetscScalar *aa;
1959: VecGetLocalSize(x,&N);
1960: 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);
1961: VecGetArray(x,&aa);
1963: PetscMalloc1(m,a);
1964: for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
1965: *a -= mstart;
1966: return(0);
1967: }
1971: /*@C
1972: VecRestoreArray2d - Restores a vector after VecGetArray2d() has been called.
1974: Logically Collective
1976: Input Parameters:
1977: + x - the vector
1978: . m - first dimension of two dimensional array
1979: . n - second dimension of the two dimensional array
1980: . mstart - first index you will use in first coordinate direction (often 0)
1981: . nstart - first index in the second coordinate direction (often 0)
1982: - a - location of pointer to array obtained from VecGetArray2d()
1984: Level: developer
1986: Notes:
1987: For regular PETSc vectors this routine does not involve any copies. For
1988: any special vectors that do not store local vector data in a contiguous
1989: array, this routine will copy the data back into the underlying
1990: vector data structure from the array obtained with VecGetArray().
1992: This routine actually zeros out the a pointer.
1994: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
1995: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
1996: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1997: @*/
1998: PetscErrorCode VecRestoreArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1999: {
2001: void *dummy;
2007: dummy = (void*)(*a + mstart);
2008: PetscFree(dummy);
2009: VecRestoreArray(x,NULL);
2010: return(0);
2011: }
2015: /*@C
2016: VecGetArray1d - Returns a pointer to a 1d contiguous array that contains this
2017: processor's portion of the vector data. You MUST call VecRestoreArray1d()
2018: when you no longer need access to the array.
2020: Logically Collective
2022: Input Parameter:
2023: + x - the vector
2024: . m - first dimension of two dimensional array
2025: - mstart - first index you will use in first coordinate direction (often 0)
2027: Output Parameter:
2028: . a - location to put pointer to the array
2030: Level: developer
2032: Notes:
2033: For a vector obtained from DMCreateLocalVector() mstart are likely
2034: obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2035: DMCreateGlobalVector() they are the corner indices from DMDAGetCorners().
2037: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2039: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2040: VecRestoreArray2d(), DMDAVecGetArray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2041: VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2042: @*/
2043: PetscErrorCode VecGetArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2044: {
2046: PetscInt N;
2052: VecGetLocalSize(x,&N);
2053: if (m != N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
2054: VecGetArray(x,a);
2055: *a -= mstart;
2056: return(0);
2057: }
2061: /*@C
2062: VecRestoreArray1d - Restores a vector after VecGetArray1d() has been called.
2064: Logically Collective
2066: Input Parameters:
2067: + x - the vector
2068: . m - first dimension of two dimensional array
2069: . mstart - first index you will use in first coordinate direction (often 0)
2070: - a - location of pointer to array obtained from VecGetArray21()
2072: Level: developer
2074: Notes:
2075: For regular PETSc vectors this routine does not involve any copies. For
2076: any special vectors that do not store local vector data in a contiguous
2077: array, this routine will copy the data back into the underlying
2078: vector data structure from the array obtained with VecGetArray1d().
2080: This routine actually zeros out the a pointer.
2082: Concepts: vector^accessing local values as 1d array
2084: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2085: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2086: VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2087: @*/
2088: PetscErrorCode VecRestoreArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
2089: {
2095: VecRestoreArray(x,NULL);
2096: return(0);
2097: }
2102: /*@C
2103: VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this
2104: processor's portion of the vector data. You MUST call VecRestoreArray3d()
2105: when you no longer need access to the array.
2107: Logically Collective
2109: Input Parameter:
2110: + x - the vector
2111: . m - first dimension of three dimensional array
2112: . n - second dimension of three dimensional array
2113: . p - third dimension of three dimensional array
2114: . mstart - first index you will use in first coordinate direction (often 0)
2115: . nstart - first index in the second coordinate direction (often 0)
2116: - pstart - first index in the third coordinate direction (often 0)
2118: Output Parameter:
2119: . a - location to put pointer to the array
2121: Level: developer
2123: Notes:
2124: For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2125: obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2126: DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2127: the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().
2129: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2131: Concepts: vector^accessing local values as 3d array
2133: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2134: VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2135: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2136: @*/
2137: PetscErrorCode VecGetArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2138: {
2140: PetscInt i,N,j;
2141: PetscScalar *aa,**b;
2147: VecGetLocalSize(x,&N);
2148: 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);
2149: VecGetArray(x,&aa);
2151: PetscMalloc1(m*sizeof(PetscScalar**)+m*n,a);
2152: b = (PetscScalar**)((*a) + m);
2153: for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2154: for (i=0; i<m; i++)
2155: for (j=0; j<n; j++)
2156: b[i*n+j] = aa + i*n*p + j*p - pstart;
2158: *a -= mstart;
2159: return(0);
2160: }
2164: /*@C
2165: VecRestoreArray3d - Restores a vector after VecGetArray3d() has been called.
2167: Logically Collective
2169: Input Parameters:
2170: + x - the vector
2171: . m - first dimension of three dimensional array
2172: . n - second dimension of the three dimensional array
2173: . p - third dimension of the three dimensional array
2174: . mstart - first index you will use in first coordinate direction (often 0)
2175: . nstart - first index in the second coordinate direction (often 0)
2176: . pstart - first index in the third coordinate direction (often 0)
2177: - a - location of pointer to array obtained from VecGetArray3d()
2179: Level: developer
2181: Notes:
2182: For regular PETSc vectors this routine does not involve any copies. For
2183: any special vectors that do not store local vector data in a contiguous
2184: array, this routine will copy the data back into the underlying
2185: vector data structure from the array obtained with VecGetArray().
2187: This routine actually zeros out the a pointer.
2189: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2190: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2191: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2192: @*/
2193: PetscErrorCode VecRestoreArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
2194: {
2196: void *dummy;
2202: dummy = (void*)(*a + mstart);
2203: PetscFree(dummy);
2204: VecRestoreArray(x,NULL);
2205: return(0);
2206: }
2210: /*@C
2211: VecGetArray4d - Returns a pointer to a 4d contiguous array that contains this
2212: processor's portion of the vector data. You MUST call VecRestoreArray4d()
2213: when you no longer need access to the array.
2215: Logically Collective
2217: Input Parameter:
2218: + x - the vector
2219: . m - first dimension of four dimensional array
2220: . n - second dimension of four dimensional array
2221: . p - third dimension of four dimensional array
2222: . q - fourth dimension of four dimensional array
2223: . mstart - first index you will use in first coordinate direction (often 0)
2224: . nstart - first index in the second coordinate direction (often 0)
2225: . pstart - first index in the third coordinate direction (often 0)
2226: - qstart - first index in the fourth coordinate direction (often 0)
2228: Output Parameter:
2229: . a - location to put pointer to the array
2231: Level: beginner
2233: Notes:
2234: For a vector obtained from DMCreateLocalVector() mstart, nstart, and pstart are likely
2235: obtained from the corner indices obtained from DMDAGetGhostCorners() while for
2236: DMCreateGlobalVector() they are the corner indices from DMDAGetCorners(). In both cases
2237: the arguments from DMDAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().
2239: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2241: Concepts: vector^accessing local values as 3d array
2243: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2244: VecRestoreArray2d(), DMDAVecGetarray(), DMDAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2245: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2246: @*/
2247: PetscErrorCode VecGetArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2248: {
2250: PetscInt i,N,j,k;
2251: PetscScalar *aa,***b,**c;
2257: VecGetLocalSize(x,&N);
2258: 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);
2259: VecGetArray(x,&aa);
2261: PetscMalloc1(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p,a);
2262: b = (PetscScalar***)((*a) + m);
2263: c = (PetscScalar**)(b + m*n);
2264: for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2265: for (i=0; i<m; i++)
2266: for (j=0; j<n; j++)
2267: b[i*n+j] = c + i*n*p + j*p - pstart;
2268: for (i=0; i<m; i++)
2269: for (j=0; j<n; j++)
2270: for (k=0; k<p; k++)
2271: c[i*n*p+j*p+k] = aa + i*n*p*q + j*p*q + k*q - qstart;
2272: *a -= mstart;
2273: return(0);
2274: }
2278: /*@C
2279: VecRestoreArray4d - Restores a vector after VecGetArray3d() has been called.
2281: Logically Collective
2283: Input Parameters:
2284: + x - the vector
2285: . m - first dimension of four dimensional array
2286: . n - second dimension of the four dimensional array
2287: . p - third dimension of the four dimensional array
2288: . q - fourth dimension of the four dimensional array
2289: . mstart - first index you will use in first coordinate direction (often 0)
2290: . nstart - first index in the second coordinate direction (often 0)
2291: . pstart - first index in the third coordinate direction (often 0)
2292: . qstart - first index in the fourth coordinate direction (often 0)
2293: - a - location of pointer to array obtained from VecGetArray4d()
2295: Level: beginner
2297: Notes:
2298: For regular PETSc vectors this routine does not involve any copies. For
2299: any special vectors that do not store local vector data in a contiguous
2300: array, this routine will copy the data back into the underlying
2301: vector data structure from the array obtained with VecGetArray().
2303: This routine actually zeros out the a pointer.
2305: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2306: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DMDAVecGetArray(), DMDAVecRestoreArray()
2307: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2308: @*/
2309: PetscErrorCode VecRestoreArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2310: {
2312: void *dummy;
2318: dummy = (void*)(*a + mstart);
2319: PetscFree(dummy);
2320: VecRestoreArray(x,NULL);
2321: return(0);
2322: }