Actual source code: rvector.c
1: #define PETSCVEC_DLL
3: /*
4: Provides the interface functions for vector operations that have PetscScalar/PetscReal in the signature
5: These are the vector functions the user calls.
6: */
7: #include private/vecimpl.h
10: if ((x)->map->N != (y)->map->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths"); \
11: if ((x)->map->n != (y)->map->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
16: /*@
17: VecMaxPointwiseDivide - Computes the maximum of the componentwise division max = max_i abs(x_i/y_i).
19: Collective on Vec
21: Input Parameters:
22: . x, y - the vectors
24: Output Parameter:
25: . max - the result
27: Level: advanced
29: Notes: any subset of the x and may be the same vector.
30: if a particular y_i is zero, it is treated as 1 in the above formula
32: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMax(), VecPointwiseMin(), VecPointwiseMaxAbs()
33: @*/
34: PetscErrorCode VecMaxPointwiseDivide(Vec x,Vec y,PetscReal *max)
35: {
47: (*x->ops->maxpointwisedivide)(x,y,max);
48: return(0);
49: }
53: /*@
54: VecDot - Computes the vector dot product.
56: Collective on Vec
58: Input Parameters:
59: . x, y - the vectors
61: Output Parameter:
62: . val - the dot product
64: Performance Issues:
65: $ per-processor memory bandwidth
66: $ interprocessor latency
67: $ work load inbalance that causes certain processes to arrive much earlier than others
69: Notes for Users of Complex Numbers:
70: For complex vectors, VecDot() computes
71: $ val = (x,y) = y^H x,
72: where y^H denotes the conjugate transpose of y.
74: Use VecTDot() for the indefinite form
75: $ val = (x,y) = y^T x,
76: where y^T denotes the transpose of y.
78: Level: intermediate
80: Concepts: inner product
81: Concepts: vector^inner product
83: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd()
84: @*/
85: PetscErrorCode VecDot(Vec x,Vec y,PetscScalar *val)
86: {
98: PetscLogEventBarrierBegin(VEC_DotBarrier,x,y,0,0,((PetscObject)x)->comm);
99: (*x->ops->dot)(x,y,val);
100: PetscLogEventBarrierEnd(VEC_DotBarrier,x,y,0,0,((PetscObject)x)->comm);
101: return(0);
102: }
106: /*@
107: VecNorm - Computes the vector norm.
109: Collective on Vec
111: Input Parameters:
112: + x - the vector
113: - type - one of NORM_1, NORM_2, NORM_INFINITY. Also available
114: NORM_1_AND_2, which computes both norms and stores them
115: in a two element array.
117: Output Parameter:
118: . val - the norm
120: Notes:
121: $ NORM_1 denotes sum_i |x_i|
122: $ NORM_2 denotes sqrt(sum_i (x_i)^2)
123: $ NORM_INFINITY denotes max_i |x_i|
125: Level: intermediate
127: Performance Issues:
128: $ per-processor memory bandwidth
129: $ interprocessor latency
130: $ work load inbalance that causes certain processes to arrive much earlier than others
132: Compile Option:
133: PETSC_HAVE_SLOW_BLAS_NORM2 will cause a C (loop unrolled) version of the norm to be used, rather
134: than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines
135: (as opposed to vendor provided) because the FORTRAN BLAS NRM2() routine is very slow.
137: Concepts: norm
138: Concepts: vector^norm
140: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNormAvailable(),
141: VecNormBegin(), VecNormEnd()
143: @*/
144: PetscErrorCode VecNorm(Vec x,NormType type,PetscReal *val)
145: {
146: PetscTruth flg;
154: /*
155: * Cached data?
156: */
157: if (type!=NORM_1_AND_2) {
158: PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,flg);
159: if (flg) return(0);
160: }
162: PetscLogEventBarrierBegin(VEC_NormBarrier,x,0,0,0,((PetscObject)x)->comm);
163: (*x->ops->norm)(x,type,val);
164: PetscLogEventBarrierEnd(VEC_NormBarrier,x,0,0,0,((PetscObject)x)->comm);
166: if (type!=NORM_1_AND_2) {
167: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[type],*val);
168: }
169: return(0);
170: }
174: /*@
175: VecNormAvailable - Returns the vector norm if it is already known.
177: Collective on Vec
179: Input Parameters:
180: + x - the vector
181: - type - one of NORM_1, NORM_2, NORM_INFINITY. Also available
182: NORM_1_AND_2, which computes both norms and stores them
183: in a two element array.
185: Output Parameter:
186: + available - PETSC_TRUE if the val returned is valid
187: - val - the norm
189: Notes:
190: $ NORM_1 denotes sum_i |x_i|
191: $ NORM_2 denotes sqrt(sum_i (x_i)^2)
192: $ NORM_INFINITY denotes max_i |x_i|
194: Level: intermediate
196: Performance Issues:
197: $ per-processor memory bandwidth
198: $ interprocessor latency
199: $ work load inbalance that causes certain processes to arrive much earlier than others
201: Compile Option:
202: PETSC_HAVE_SLOW_BLAS_NORM2 will cause a C (loop unrolled) version of the norm to be used, rather
203: than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines
204: (as opposed to vendor provided) because the FORTRAN BLAS NRM2() routine is very slow.
206: Concepts: norm
207: Concepts: vector^norm
209: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(), VecNorm()
210: VecNormBegin(), VecNormEnd()
212: @*/
213: PetscErrorCode VecNormAvailable(Vec x,NormType type,PetscTruth *available,PetscReal *val)
214: {
222: *available = PETSC_FALSE;
223: if (type!=NORM_1_AND_2) {
224: PetscObjectComposedDataGetReal((PetscObject)x,NormIds[type],*val,*available);
225: }
226: return(0);
227: }
231: /*@
232: VecNormalize - Normalizes a vector by 2-norm.
234: Collective on Vec
236: Input Parameters:
237: + x - the vector
239: Output Parameter:
240: . x - the normalized vector
241: - val - the vector norm before normalization
243: Level: intermediate
245: Concepts: vector^normalizing
246: Concepts: normalizing^vector
248: @*/
249: PetscErrorCode VecNormalize(Vec x,PetscReal *val)
250: {
252: PetscReal norm;
257: PetscLogEventBegin(VEC_Normalize,x,0,0,0);
258: VecNorm(x,NORM_2,&norm);
259: if (norm == 0.0) {
260: PetscInfo(x,"Vector of zero norm can not be normalized; Returning only the zero norm\n");
261: } else if (norm != 1.0) {
262: PetscScalar tmp = 1.0/norm;
263: VecScale(x,tmp);
264: }
265: if (val) *val = norm;
266: PetscLogEventEnd(VEC_Normalize,x,0,0,0);
267: return(0);
268: }
272: /*@C
273: VecMax - Determines the maximum vector component and its location.
275: Collective on Vec
277: Input Parameter:
278: . x - the vector
280: Output Parameters:
281: + val - the maximum component
282: - p - the location of val
284: Notes:
285: Returns the value PETSC_MIN and p = -1 if the vector is of length 0.
287: Level: intermediate
289: Concepts: maximum^of vector
290: Concepts: vector^maximum value
292: .seealso: VecNorm(), VecMin()
293: @*/
294: PetscErrorCode VecMax(Vec x,PetscInt *p,PetscReal *val)
295: {
302: PetscLogEventBegin(VEC_Max,x,0,0,0);
303: (*x->ops->max)(x,p,val);
304: PetscLogEventEnd(VEC_Max,x,0,0,0);
305: return(0);
306: }
310: /*@
311: VecMin - Determines the minimum vector component and its location.
313: Collective on Vec
315: Input Parameters:
316: . x - the vector
318: Output Parameter:
319: + val - the minimum component
320: - p - the location of val
322: Level: intermediate
324: Notes:
325: Returns the value PETSC_MAX and p = -1 if the vector is of length 0.
327: Concepts: minimum^of vector
328: Concepts: vector^minimum entry
330: .seealso: VecMax()
331: @*/
332: PetscErrorCode VecMin(Vec x,PetscInt *p,PetscReal *val)
333: {
340: PetscLogEventBegin(VEC_Min,x,0,0,0);
341: (*x->ops->min)(x,p,val);
342: PetscLogEventEnd(VEC_Min,x,0,0,0);
343: return(0);
344: }
348: /*@
349: VecTDot - Computes an indefinite vector dot product. That is, this
350: routine does NOT use the complex conjugate.
352: Collective on Vec
354: Input Parameters:
355: . x, y - the vectors
357: Output Parameter:
358: . val - the dot product
360: Notes for Users of Complex Numbers:
361: For complex vectors, VecTDot() computes the indefinite form
362: $ val = (x,y) = y^T x,
363: where y^T denotes the transpose of y.
365: Use VecDot() for the inner product
366: $ val = (x,y) = y^H x,
367: where y^H denotes the conjugate transpose of y.
369: Level: intermediate
371: Concepts: inner product^non-Hermitian
372: Concepts: vector^inner product
373: Concepts: non-Hermitian inner product
375: .seealso: VecDot(), VecMTDot()
376: @*/
377: PetscErrorCode VecTDot(Vec x,Vec y,PetscScalar *val)
378: {
390: PetscLogEventBegin(VEC_TDot,x,y,0,0);
391: (*x->ops->tdot)(x,y,val);
392: PetscLogEventEnd(VEC_TDot,x,y,0,0);
393: return(0);
394: }
398: /*@
399: VecScale - Scales a vector.
401: Collective on Vec
403: Input Parameters:
404: + x - the vector
405: - alpha - the scalar
407: Output Parameter:
408: . x - the scaled vector
410: Note:
411: For a vector with n components, VecScale() computes
412: $ x[i] = alpha * x[i], for i=1,...,n.
414: Level: intermediate
416: Concepts: vector^scaling
417: Concepts: scaling^vector
419: @*/
420: PetscErrorCode VecScale (Vec x, PetscScalar alpha)
421: {
422: PetscReal norms[4] = {0.0,0.0,0.0, 0.0};
423: PetscTruth flgs[4];
425: PetscInt i;
430: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
431: PetscLogEventBegin(VEC_Scale,x,0,0,0);
432: if (alpha != 1.0) {
433: (*x->ops->scale)(x,alpha);
435: /*
436: * Update cached data
437: */
438: for (i=0; i<4; i++) {
439: PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],norms[i],flgs[i]);
440: }
442: /* in general we consider this object touched */
443: PetscObjectStateIncrease((PetscObject)x);
445: for (i=0; i<4; i++) {
446: if (flgs[i]) {
447: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[i],PetscAbsScalar(alpha)*norms[i]);
448: }
449: }
450: }
451: PetscLogEventEnd(VEC_Scale,x,0,0,0);
452: return(0);
453: }
458: /*@
459: VecSet - Sets all components of a vector to a single scalar value.
461: Collective on Vec
463: Input Parameters:
464: + x - the vector
465: - alpha - the scalar
467: Output Parameter:
468: . x - the vector
470: Note:
471: For a vector of dimension n, VecSet() computes
472: $ x[i] = alpha, for i=1,...,n,
473: so that all vector entries then equal the identical
474: scalar value, alpha. Use the more general routine
475: VecSetValues() to set different vector entries.
477: You CANNOT call this after you have called VecSetValues() but before you call
478: VecAssemblyBegin/End().
480: Level: beginner
482: .seealso VecSetValues(), VecSetValuesBlocked(), VecSetRandom()
484: Concepts: vector^setting to constant
486: @*/
487: PetscErrorCode VecSet(Vec x,PetscScalar alpha)
488: {
489: PetscReal val;
495: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You cannot call this after you have called VecSetValues() but\n before you have called VecAssemblyBegin/End()");
496: #if defined (PETSC_USE_DEBUG)
497: {
498: PetscReal alpha_local,alpha_max;
499: alpha_local = PetscAbsScalar(alpha);
500: MPI_Allreduce(&alpha_local,&alpha_max,1,MPIU_REAL,MPI_MAX,((PetscObject)x)->comm);
501: if (alpha_local != alpha_max) SETERRQ(PETSC_ERR_ARG_WRONG,"Same value should be used across all processors");
502: }
503: #endif
505: PetscLogEventBegin(VEC_Set,x,0,0,0);
506: (*x->ops->set)(x,alpha);
507: PetscLogEventEnd(VEC_Set,x,0,0,0);
509: /*
510: * Update cached data
511: */
512: /* in general we consider this object touched */
513: PetscObjectStateIncrease((PetscObject)x);
515: /* however, norms can be simply set */
516: val = PetscAbsScalar(alpha);
517: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_1],x->map->N * val);
518: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_INFINITY],val);
519: val = sqrt((double)x->map->N) * val;
520: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_2],val);
521: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[NORM_FROBENIUS],val);
522: return(0);
523: }
528: /*@
529: VecAXPY - Computes y = alpha x + y.
531: Collective on Vec
533: Input Parameters:
534: + alpha - the scalar
535: - x, y - the vectors
537: Output Parameter:
538: . y - output vector
540: Level: intermediate
542: Notes: x and y must be different vectors
544: Concepts: vector^BLAS
545: Concepts: BLAS
547: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
548: @*/
549: PetscErrorCode VecAXPY(Vec y,PetscScalar alpha,Vec x)
550: {
561: PetscLogEventBegin(VEC_AXPY,x,y,0,0);
562: (*y->ops->axpy)(y,alpha,x);
563: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
564: PetscObjectStateIncrease((PetscObject)y);
565: return(0);
566: }
570: /*@
571: VecAXPBY - Computes y = alpha x + beta y.
573: Collective on Vec
575: Input Parameters:
576: + alpha,beta - the scalars
577: - x, y - the vectors
579: Output Parameter:
580: . y - output vector
582: Level: intermediate
584: Notes: x and y must be different vectors
586: Concepts: BLAS
587: Concepts: vector^BLAS
589: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
590: @*/
591: PetscErrorCode VecAXPBY(Vec y,PetscScalar alpha,PetscScalar beta,Vec x)
592: {
603: PetscLogEventBegin(VEC_AXPY,x,y,0,0);
604: (*y->ops->axpby)(y,alpha,beta,x);
605: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
606: PetscObjectStateIncrease((PetscObject)y);
607: return(0);
608: }
612: /*@
613: VecAXPBYPCZ - Computes z = alpha x + beta y + gamma z
615: Collective on Vec
617: Input Parameters:
618: + alpha,beta, gamma - the scalars
619: - x, y, z - the vectors
621: Output Parameter:
622: . z - output vector
624: Level: intermediate
626: Notes: x, y and z must be different vectors
628: alpha = 1 or gamma = 1 are handled as special cases
630: Concepts: BLAS
631: Concepts: vector^BLAS
633: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
634: @*/
635: PetscErrorCode VecAXPBYPCZ(Vec z,PetscScalar alpha,PetscScalar beta,PetscScalar gamma,Vec x,Vec y)
636: {
651: PetscLogEventBegin(VEC_AXPBYPCZ,x,y,z,0);
652: (*y->ops->axpbypcz)(z,alpha,beta,gamma,x,y);
653: PetscLogEventEnd(VEC_AXPBYPCZ,x,y,z,0);
654: PetscObjectStateIncrease((PetscObject)y);
655: return(0);
656: }
660: /*@
661: VecAYPX - Computes y = x + alpha y.
663: Collective on Vec
665: Input Parameters:
666: + alpha - the scalar
667: - x, y - the vectors
669: Output Parameter:
670: . y - output vector
672: Level: intermediate
674: Notes: x and y must be different vectors
676: Concepts: vector^BLAS
677: Concepts: BLAS
679: .seealso: VecAXPY(), VecWAXPY()
680: @*/
681: PetscErrorCode VecAYPX(Vec y,PetscScalar alpha,Vec x)
682: {
691: PetscLogEventBegin(VEC_AYPX,x,y,0,0);
692: (*y->ops->aypx)(y,alpha,x);
693: PetscLogEventEnd(VEC_AYPX,x,y,0,0);
694: PetscObjectStateIncrease((PetscObject)y);
695: return(0);
696: }
701: /*@
702: VecWAXPY - Computes w = alpha x + y.
704: Collective on Vec
706: Input Parameters:
707: + alpha - the scalar
708: - x, y - the vectors
710: Output Parameter:
711: . w - the result
713: Level: intermediate
715: Notes: Neither the vector x or y can be the same as vector w
717: Concepts: vector^BLAS
718: Concepts: BLAS
720: .seealso: VecAXPY(), VecAYPX(), VecAXPBY()
721: @*/
722: PetscErrorCode VecWAXPY(Vec w,PetscScalar alpha,Vec x,Vec y)
723: {
737: if (w == y) SETERRQ(PETSC_ERR_SUP,"Result vector w cannot be same as input vector y, suggest VecAXPY()");
738: if (w == x) SETERRQ(PETSC_ERR_SUP,"Result vector w cannot be same as input vector x, suggest VecAYPX()");
740: PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
741: (*w->ops->waxpy)(w,alpha,x,y);
742: PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
743: PetscObjectStateIncrease((PetscObject)w);
744: return(0);
745: }
750: /*@
751: VecSetValues - Inserts or adds values into certain locations of a vector.
753: Not Collective
755: Input Parameters:
756: + x - vector to insert in
757: . ni - number of elements to add
758: . ix - indices where to add
759: . y - array of values
760: - iora - either INSERT_VALUES or ADD_VALUES, where
761: ADD_VALUES adds values to any existing entries, and
762: INSERT_VALUES replaces existing entries with new values
764: Notes:
765: VecSetValues() sets x[ix[i]] = y[i], for i=0,...,ni-1.
767: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
768: options cannot be mixed without intervening calls to the assembly
769: routines.
771: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
772: MUST be called after all calls to VecSetValues() have been completed.
774: VecSetValues() uses 0-based indices in Fortran as well as in C.
776: If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
777: negative indices may be passed in ix. These rows are
778: simply ignored. This allows easily inserting element load matrices
779: with homogeneous Dirchlet boundary conditions that you don't want represented
780: in the vector.
782: Level: beginner
784: Concepts: vector^setting values
786: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
787: VecSetValue(), VecSetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecGetValues()
788: @*/
789: PetscErrorCode VecSetValues(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
790: {
798: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
799: (*x->ops->setvalues)(x,ni,ix,y,iora);
800: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
801: PetscObjectStateIncrease((PetscObject)x);
802: return(0);
803: }
807: /*@
808: VecGetValues - Gets values from certain locations of a vector. Currently
809: can only get values on the same processor
811: Collective on Vec
812:
813: Input Parameters:
814: + x - vector to get values from
815: . ni - number of elements to get
816: - ix - indices where to get them from (in global 1d numbering)
818: Output Parameter:
819: . y - array of values
821: Notes:
822: The user provides the allocated array y; it is NOT allocated in this routine
824: VecGetValues() gets y[i] = x[ix[i]], for i=0,...,ni-1.
826: VecAssemblyBegin() and VecAssemblyEnd() MUST be called before calling this
828: VecGetValues() uses 0-based indices in Fortran as well as in C.
830: If you call VecSetOption(x, VEC_IGNORE_NEGATIVE_INDICES,PETSC_TRUE),
831: negative indices may be passed in ix. These rows are
832: simply ignored.
834: Level: beginner
836: Concepts: vector^getting values
838: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecGetValuesLocal(),
839: VecGetValuesBlocked(), InsertMode, INSERT_VALUES, ADD_VALUES, VecSetValues()
840: @*/
841: PetscErrorCode VecGetValues(Vec x,PetscInt ni,const PetscInt ix[],PetscScalar y[])
842: {
850: (*x->ops->getvalues)(x,ni,ix,y);
851: return(0);
852: }
856: /*@
857: VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector.
859: Not Collective
861: Input Parameters:
862: + x - vector to insert in
863: . ni - number of blocks to add
864: . ix - indices where to add in block count, rather than element count
865: . y - array of values
866: - iora - either INSERT_VALUES or ADD_VALUES, where
867: ADD_VALUES adds values to any existing entries, and
868: INSERT_VALUES replaces existing entries with new values
870: Notes:
871: VecSetValuesBlocked() sets x[bs*ix[i]+j] = y[bs*i+j],
872: for j=0,...,bs, for i=0,...,ni-1. where bs was set with VecSetBlockSize().
874: Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
875: options cannot be mixed without intervening calls to the assembly
876: routines.
878: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
879: MUST be called after all calls to VecSetValuesBlocked() have been completed.
881: VecSetValuesBlocked() uses 0-based indices in Fortran as well as in C.
883: Negative indices may be passed in ix, these rows are
884: simply ignored. This allows easily inserting element load matrices
885: with homogeneous Dirchlet boundary conditions that you don't want represented
886: in the vector.
888: Level: intermediate
890: Concepts: vector^setting values blocked
892: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
893: VecSetValues()
894: @*/
895: PetscErrorCode VecSetValuesBlocked(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
896: {
904: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
905: (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
906: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
907: PetscObjectStateIncrease((PetscObject)x);
908: return(0);
909: }
914: /*@
915: VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
916: using a local ordering of the nodes.
918: Not Collective
920: Input Parameters:
921: + x - vector to insert in
922: . ni - number of elements to add
923: . ix - indices where to add
924: . y - array of values
925: - iora - either INSERT_VALUES or ADD_VALUES, where
926: ADD_VALUES adds values to any existing entries, and
927: INSERT_VALUES replaces existing entries with new values
929: Level: intermediate
931: Notes:
932: VecSetValuesLocal() sets x[ix[i]] = y[i], for i=0,...,ni-1.
934: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
935: options cannot be mixed without intervening calls to the assembly
936: routines.
938: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
939: MUST be called after all calls to VecSetValuesLocal() have been completed.
941: VecSetValuesLocal() uses 0-based indices in Fortran as well as in C.
943: Concepts: vector^setting values with local numbering
945: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
946: VecSetValuesBlockedLocal()
947: @*/
948: PetscErrorCode VecSetValuesLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
949: {
951: PetscInt lixp[128],*lix = lixp;
959: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
960: if (!x->ops->setvalueslocal) {
961: if (!x->mapping) {
962: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
963: }
964: if (ni > 128) {
965: PetscMalloc(ni*sizeof(PetscInt),&lix);
966: }
967: ISLocalToGlobalMappingApply(x->mapping,ni,(PetscInt*)ix,lix);
968: (*x->ops->setvalues)(x,ni,lix,y,iora);
969: if (ni > 128) {
970: PetscFree(lix);
971: }
972: } else {
973: (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
974: }
975: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
976: PetscObjectStateIncrease((PetscObject)x);
977: return(0);
978: }
982: /*@
983: VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
984: using a local ordering of the nodes.
986: Not Collective
988: Input Parameters:
989: + x - vector to insert in
990: . ni - number of blocks to add
991: . ix - indices where to add in block count, not element count
992: . y - array of values
993: - iora - either INSERT_VALUES or ADD_VALUES, where
994: ADD_VALUES adds values to any existing entries, and
995: INSERT_VALUES replaces existing entries with new values
997: Level: intermediate
999: Notes:
1000: VecSetValuesBlockedLocal() sets x[bs*ix[i]+j] = y[bs*i+j],
1001: for j=0,..bs-1, for i=0,...,ni-1, where bs has been set with VecSetBlockSize().
1003: Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1004: options cannot be mixed without intervening calls to the assembly
1005: routines.
1007: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1008: MUST be called after all calls to VecSetValuesBlockedLocal() have been completed.
1010: VecSetValuesBlockedLocal() uses 0-based indices in Fortran as well as in C.
1013: Concepts: vector^setting values blocked with local numbering
1015: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1016: VecSetLocalToGlobalMappingBlock()
1017: @*/
1018: PetscErrorCode VecSetValuesBlockedLocal(Vec x,PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora)
1019: {
1021: PetscInt lixp[128],*lix = lixp;
1028: if (!x->bmapping) {
1029: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMappingBlock()");
1030: }
1031: if (ni > 128) {
1032: PetscMalloc(ni*sizeof(PetscInt),&lix);
1033: }
1035: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1036: ISLocalToGlobalMappingApply(x->bmapping,ni,(PetscInt*)ix,lix);
1037: (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1038: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1039: if (ni > 128) {
1040: PetscFree(lix);
1041: }
1042: PetscObjectStateIncrease((PetscObject)x);
1043: return(0);
1044: }
1050: /*@
1051: VecMTDot - Computes indefinite vector multiple dot products.
1052: That is, it does NOT use the complex conjugate.
1054: Collective on Vec
1056: Input Parameters:
1057: + x - one vector
1058: . nv - number of vectors
1059: - y - array of vectors. Note that vectors are pointers
1061: Output Parameter:
1062: . val - array of the dot products
1064: Notes for Users of Complex Numbers:
1065: For complex vectors, VecMTDot() computes the indefinite form
1066: $ val = (x,y) = y^T x,
1067: where y^T denotes the transpose of y.
1069: Use VecMDot() for the inner product
1070: $ val = (x,y) = y^H x,
1071: where y^H denotes the conjugate transpose of y.
1073: Level: intermediate
1075: Concepts: inner product^multiple
1076: Concepts: vector^multiple inner products
1078: .seealso: VecMDot(), VecTDot()
1079: @*/
1080: PetscErrorCode VecMTDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1081: {
1094: PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1095: (*x->ops->mtdot)(x,nv,y,val);
1096: PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1097: return(0);
1098: }
1102: /*@
1103: VecMDot - Computes vector multiple dot products.
1105: Collective on Vec
1107: Input Parameters:
1108: + x - one vector
1109: . nv - number of vectors
1110: - y - array of vectors.
1112: Output Parameter:
1113: . val - array of the dot products (does not allocate the array)
1115: Notes for Users of Complex Numbers:
1116: For complex vectors, VecMDot() computes
1117: $ val = (x,y) = y^H x,
1118: where y^H denotes the conjugate transpose of y.
1120: Use VecMTDot() for the indefinite form
1121: $ val = (x,y) = y^T x,
1122: where y^T denotes the transpose of y.
1124: Level: intermediate
1126: Concepts: inner product^multiple
1127: Concepts: vector^multiple inner products
1129: .seealso: VecMTDot(), VecDot()
1130: @*/
1131: PetscErrorCode VecMDot(Vec x,PetscInt nv,const Vec y[],PetscScalar val[])
1132: {
1137: if (!nv) return(0);
1138: if (nv < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1147: PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,((PetscObject)x)->comm);
1148: (*x->ops->mdot)(x,nv,y,val);
1149: PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,((PetscObject)x)->comm);
1150: return(0);
1151: }
1155: /*@
1156: VecMAXPY - Computes y = y + sum alpha[j] x[j]
1158: Collective on Vec
1160: Input Parameters:
1161: + nv - number of scalars and x-vectors
1162: . alpha - array of scalars
1163: . y - one vector
1164: - x - array of vectors
1166: Level: intermediate
1168: Concepts: BLAS
1170: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1171: @*/
1172: PetscErrorCode VecMAXPY(Vec y,PetscInt nv,const PetscScalar alpha[],Vec x[])
1173: {
1178: if (!nv) return(0);
1179: if (nv < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of vectors (given %D) cannot be negative",nv);
1188: PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1189: (*y->ops->maxpy)(y,nv,alpha,x);
1190: PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1191: PetscObjectStateIncrease((PetscObject)y);
1192: return(0);
1193: }
1195: /*MC
1196: VecGetArray - Returns a pointer to a contiguous array that contains this
1197: processor's portion of the vector data. For the standard PETSc
1198: vectors, VecGetArray() returns a pointer to the local data array and
1199: does not use any copies. If the underlying vector data is not stored
1200: in a contiquous array this routine will copy the data to a contiquous
1201: array and return a pointer to that. You MUST call VecRestoreArray()
1202: when you no longer need access to the array.
1204: Synopsis:
1205: PetscErrorCode VecGetArray(Vec x,PetscScalar *a[])
1207: Not Collective
1209: Input Parameter:
1210: . x - the vector
1212: Output Parameter:
1213: . a - location to put pointer to the array
1215: Fortran Note:
1216: This routine is used differently from Fortran 77
1217: $ Vec x
1218: $ PetscScalar x_array(1)
1219: $ PetscOffset i_x
1220: $ PetscErrorCode ierr
1221: $ call VecGetArray(x,x_array,i_x,ierr)
1222: $
1223: $ Access first local entry in vector with
1224: $ value = x_array(i_x + 1)
1225: $
1226: $ ...... other code
1227: $ call VecRestoreArray(x,x_array,i_x,ierr)
1228: For Fortran 90 see VecGetArrayF90()
1230: See the Fortran chapter of the users manual and
1231: petsc/src/snes/examples/tutorials/ex5f.F for details.
1233: Level: beginner
1235: Concepts: vector^accessing local values
1237: .seealso: VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1238: M*/
1241: PetscErrorCode VecGetArray_Private(Vec x,PetscScalar *a[])
1242: {
1249: (*x->ops->getarray)(x,a);
1250: return(0);
1251: }
1256: /*@C
1257: VecGetArrays - Returns a pointer to the arrays in a set of vectors
1258: that were created by a call to VecDuplicateVecs(). You MUST call
1259: VecRestoreArrays() when you no longer need access to the array.
1261: Not Collective
1263: Input Parameter:
1264: + x - the vectors
1265: - n - the number of vectors
1267: Output Parameter:
1268: . a - location to put pointer to the array
1270: Fortran Note:
1271: This routine is not supported in Fortran.
1273: Level: intermediate
1275: .seealso: VecGetArray(), VecRestoreArrays()
1276: @*/
1277: PetscErrorCode VecGetArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1278: {
1280: PetscInt i;
1281: PetscScalar **q;
1287: if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %D",n);
1288: PetscMalloc(n*sizeof(PetscScalar*),&q);
1289: for (i=0; i<n; ++i) {
1290: VecGetArray(x[i],&q[i]);
1291: }
1292: *a = q;
1293: return(0);
1294: }
1298: /*@C
1299: VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1300: has been called.
1302: Not Collective
1304: Input Parameters:
1305: + x - the vector
1306: . n - the number of vectors
1307: - a - location of pointer to arrays obtained from VecGetArrays()
1309: Notes:
1310: For regular PETSc vectors this routine does not involve any copies. For
1311: any special vectors that do not store local vector data in a contiguous
1312: array, this routine will copy the data back into the underlying
1313: vector data structure from the arrays obtained with VecGetArrays().
1315: Fortran Note:
1316: This routine is not supported in Fortran.
1318: Level: intermediate
1320: .seealso: VecGetArrays(), VecRestoreArray()
1321: @*/
1322: PetscErrorCode VecRestoreArrays(const Vec x[],PetscInt n,PetscScalar **a[])
1323: {
1325: PetscInt i;
1326: PetscScalar **q = *a;
1333: for(i=0;i<n;++i) {
1334: VecRestoreArray(x[i],&q[i]);
1335: }
1336: PetscFree(q);
1337: return(0);
1338: }
1340: /*MC
1341: VecRestoreArray - Restores a vector after VecGetArray() has been called.
1343: Synopsis:
1344: PetscErrorCode VecRestoreArray(Vec x,PetscScalar *a[])
1346: Not Collective
1348: Input Parameters:
1349: + x - the vector
1350: - a - location of pointer to array obtained from VecGetArray()
1352: Level: beginner
1354: Notes:
1355: For regular PETSc vectors this routine does not involve any copies. For
1356: any special vectors that do not store local vector data in a contiguous
1357: array, this routine will copy the data back into the underlying
1358: vector data structure from the array obtained with VecGetArray().
1360: This routine actually zeros out the a pointer. This is to prevent accidental
1361: us of the array after it has been restored. If you pass null for a it will
1362: not zero the array pointer a.
1364: Fortran Note:
1365: This routine is used differently from Fortran 77
1366: $ Vec x
1367: $ PetscScalar x_array(1)
1368: $ PetscOffset i_x
1369: $ PetscErrorCode ierr
1370: $ call VecGetArray(x,x_array,i_x,ierr)
1371: $
1372: $ Access first local entry in vector with
1373: $ value = x_array(i_x + 1)
1374: $
1375: $ ...... other code
1376: $ call VecRestoreArray(x,x_array,i_x,ierr)
1378: See the Fortran chapter of the users manual and
1379: petsc/src/snes/examples/tutorials/ex5f.F for details.
1380: For Fortran 90 see VecRestoreArrayF90()
1382: .seealso: VecGetArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1383: M*/
1386: PetscErrorCode VecRestoreArray_Private(Vec x,PetscScalar *a[])
1387: {
1394: #if defined(PETSC_USE_DEBUG)
1395: CHKMEMQ;
1396: #endif
1397: if (x->ops->restorearray) {
1398: (*x->ops->restorearray)(x,a);
1399: }
1400: PetscObjectStateIncrease((PetscObject)x);
1401: return(0);
1402: }
1407: /*@
1408: VecPlaceArray - Allows one to replace the array in a vector with an
1409: array provided by the user. This is useful to avoid copying an array
1410: into a vector.
1412: Not Collective
1414: Input Parameters:
1415: + vec - the vector
1416: - array - the array
1418: Notes:
1419: You can return to the original array with a call to VecResetArray()
1421: Level: developer
1423: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()
1425: @*/
1426: PetscErrorCode VecPlaceArray(Vec vec,const PetscScalar array[])
1427: {
1434: if (vec->ops->placearray) {
1435: (*vec->ops->placearray)(vec,array);
1436: } else {
1437: SETERRQ(PETSC_ERR_SUP,"Cannot place array in this type of vector");
1438: }
1439: PetscObjectStateIncrease((PetscObject)vec);
1440: return(0);
1441: }
1446: /*@C
1447: VecReplaceArray - Allows one to replace the array in a vector with an
1448: array provided by the user. This is useful to avoid copying an array
1449: into a vector.
1451: Not Collective
1453: Input Parameters:
1454: + vec - the vector
1455: - array - the array
1457: Notes:
1458: This permanently replaces the array and frees the memory associated
1459: with the old array.
1461: The memory passed in MUST be obtained with PetscMalloc() and CANNOT be
1462: freed by the user. It will be freed when the vector is destroy.
1464: Not supported from Fortran
1466: Level: developer
1468: .seealso: VecGetArray(), VecRestoreArray(), VecPlaceArray(), VecResetArray()
1470: @*/
1471: PetscErrorCode VecReplaceArray(Vec vec,const PetscScalar array[])
1472: {
1478: if (vec->ops->replacearray) {
1479: (*vec->ops->replacearray)(vec,array);
1480: } else {
1481: SETERRQ(PETSC_ERR_SUP,"Cannot replace array in this type of vector");
1482: }
1483: PetscObjectStateIncrease((PetscObject)vec);
1484: return(0);
1485: }
1487: /*MC
1488: VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
1489: and makes them accessible via a Fortran90 pointer.
1491: Synopsis:
1492: VecDuplicateVecsF90(Vec x,int n,{Vec, pointer :: y(:)},integer ierr)
1494: Collective on Vec
1496: Input Parameters:
1497: + x - a vector to mimic
1498: - n - the number of vectors to obtain
1500: Output Parameters:
1501: + y - Fortran90 pointer to the array of vectors
1502: - ierr - error code
1504: Example of Usage:
1505: .vb
1506: Vec x
1507: Vec, pointer :: y(:)
1508: ....
1509: call VecDuplicateVecsF90(x,2,y,ierr)
1510: call VecSet(y(2),alpha,ierr)
1511: call VecSet(y(2),alpha,ierr)
1512: ....
1513: call VecDestroyVecsF90(y,2,ierr)
1514: .ve
1516: Notes:
1517: Not yet supported for all F90 compilers
1519: Use VecDestroyVecsF90() to free the space.
1521: Level: beginner
1523: .seealso: VecDestroyVecsF90(), VecDuplicateVecs()
1525: M*/
1527: /*MC
1528: VecRestoreArrayF90 - Restores a vector to a usable state after a call to
1529: VecGetArrayF90().
1531: Synopsis:
1532: VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
1534: Not collective
1536: Input Parameters:
1537: + x - vector
1538: - xx_v - the Fortran90 pointer to the array
1540: Output Parameter:
1541: . ierr - error code
1543: Example of Usage:
1544: .vb
1545: PetscScalar, pointer :: xx_v(:)
1546: ....
1547: call VecGetArrayF90(x,xx_v,ierr)
1548: a = xx_v(3)
1549: call VecRestoreArrayF90(x,xx_v,ierr)
1550: .ve
1551:
1552: Level: beginner
1554: .seealso: VecGetArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran
1556: M*/
1558: /*MC
1559: VecDestroyVecsF90 - Frees a block of vectors obtained with VecDuplicateVecsF90().
1561: Synopsis:
1562: VecDestroyVecsF90({Vec, pointer :: x(:)},integer n,integer ierr)
1564: Collective on Vec
1566: Input Parameters:
1567: + x - pointer to array of vector pointers
1568: - n - the number of vectors previously obtained
1570: Output Parameter:
1571: . ierr - error code
1573: Notes:
1574: Not yet supported for all F90 compilers
1576: Level: beginner
1578: .seealso: VecDestroyVecs(), VecDuplicateVecsF90()
1580: M*/
1582: /*MC
1583: VecGetArrayF90 - Accesses a vector array from Fortran90. For default PETSc
1584: vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
1585: this routine is implementation dependent. You MUST call VecRestoreArrayF90()
1586: when you no longer need access to the array.
1588: Synopsis:
1589: VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
1591: Not Collective
1593: Input Parameter:
1594: . x - vector
1596: Output Parameters:
1597: + xx_v - the Fortran90 pointer to the array
1598: - ierr - error code
1600: Example of Usage:
1601: .vb
1602: PetscScalar, pointer :: xx_v(:)
1603: ....
1604: call VecGetArrayF90(x,xx_v,ierr)
1605: a = xx_v(3)
1606: call VecRestoreArrayF90(x,xx_v,ierr)
1607: .ve
1609: Level: beginner
1611: .seealso: VecRestoreArrayF90(), VecGetArray(), VecRestoreArray(), UsingFortran
1613: M*/
1618: /*@C
1619: VecGetArray2d - Returns a pointer to a 2d contiguous array that contains this
1620: processor's portion of the vector data. You MUST call VecRestoreArray2d()
1621: when you no longer need access to the array.
1623: Not Collective
1625: Input Parameter:
1626: + x - the vector
1627: . m - first dimension of two dimensional array
1628: . n - second dimension of two dimensional array
1629: . mstart - first index you will use in first coordinate direction (often 0)
1630: - nstart - first index in the second coordinate direction (often 0)
1632: Output Parameter:
1633: . a - location to put pointer to the array
1635: Level: developer
1637: Notes:
1638: For a vector obtained from DACreateLocalVector() mstart and nstart are likely
1639: obtained from the corner indices obtained from DAGetGhostCorners() while for
1640: DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
1641: the arguments from DAGet[Ghost]Corners() are reversed in the call to VecGetArray2d().
1642:
1643: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
1645: Concepts: vector^accessing local values as 2d array
1647: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1648: VecRestoreArray2d(), DAVecGetArray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1649: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1650: @*/
1651: PetscErrorCode VecGetArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1652: {
1654: PetscInt i,N;
1655: PetscScalar *aa;
1661: VecGetLocalSize(x,&N);
1662: if (m*n != N) SETERRQ3(PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 2d array dimensions %D by %D",N,m,n);
1663: VecGetArray(x,&aa);
1665: PetscMalloc(m*sizeof(PetscScalar*),a);
1666: for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
1667: *a -= mstart;
1668: return(0);
1669: }
1673: /*@C
1674: VecRestoreArray2d - Restores a vector after VecGetArray2d() has been called.
1676: Not Collective
1678: Input Parameters:
1679: + x - the vector
1680: . m - first dimension of two dimensional array
1681: . n - second dimension of the two dimensional array
1682: . mstart - first index you will use in first coordinate direction (often 0)
1683: . nstart - first index in the second coordinate direction (often 0)
1684: - a - location of pointer to array obtained from VecGetArray2d()
1686: Level: developer
1688: Notes:
1689: For regular PETSc vectors this routine does not involve any copies. For
1690: any special vectors that do not store local vector data in a contiguous
1691: array, this routine will copy the data back into the underlying
1692: vector data structure from the array obtained with VecGetArray().
1694: This routine actually zeros out the a pointer.
1696: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
1697: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
1698: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1699: @*/
1700: PetscErrorCode VecRestoreArray2d(Vec x,PetscInt m,PetscInt n,PetscInt mstart,PetscInt nstart,PetscScalar **a[])
1701: {
1703: void *dummy;
1709: dummy = (void*)(*a + mstart);
1710: PetscFree(dummy);
1711: VecRestoreArray(x,PETSC_NULL);
1712: return(0);
1713: }
1717: /*@C
1718: VecGetArray1d - Returns a pointer to a 1d contiguous array that contains this
1719: processor's portion of the vector data. You MUST call VecRestoreArray1d()
1720: when you no longer need access to the array.
1722: Not Collective
1724: Input Parameter:
1725: + x - the vector
1726: . m - first dimension of two dimensional array
1727: - mstart - first index you will use in first coordinate direction (often 0)
1729: Output Parameter:
1730: . a - location to put pointer to the array
1732: Level: developer
1734: Notes:
1735: For a vector obtained from DACreateLocalVector() mstart are likely
1736: obtained from the corner indices obtained from DAGetGhostCorners() while for
1737: DACreateGlobalVector() they are the corner indices from DAGetCorners().
1738:
1739: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
1741: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1742: VecRestoreArray2d(), DAVecGetArray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1743: VecGetArray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1744: @*/
1745: PetscErrorCode VecGetArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
1746: {
1748: PetscInt N;
1754: VecGetLocalSize(x,&N);
1755: if (m != N) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Local array size %D does not match 1d array dimensions %D",N,m);
1756: VecGetArray(x,a);
1757: *a -= mstart;
1758: return(0);
1759: }
1763: /*@C
1764: VecRestoreArray1d - Restores a vector after VecGetArray1d() has been called.
1766: Not Collective
1768: Input Parameters:
1769: + x - the vector
1770: . m - first dimension of two dimensional array
1771: . mstart - first index you will use in first coordinate direction (often 0)
1772: - a - location of pointer to array obtained from VecGetArray21()
1774: Level: developer
1776: Notes:
1777: For regular PETSc vectors this routine does not involve any copies. For
1778: any special vectors that do not store local vector data in a contiguous
1779: array, this routine will copy the data back into the underlying
1780: vector data structure from the array obtained with VecGetArray1d().
1782: This routine actually zeros out the a pointer.
1784: Concepts: vector^accessing local values as 1d array
1786: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
1787: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
1788: VecGetArray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
1789: @*/
1790: PetscErrorCode VecRestoreArray1d(Vec x,PetscInt m,PetscInt mstart,PetscScalar *a[])
1791: {
1797: VecRestoreArray(x,PETSC_NULL);
1798: return(0);
1799: }
1804: /*@C
1805: VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this
1806: processor's portion of the vector data. You MUST call VecRestoreArray3d()
1807: when you no longer need access to the array.
1809: Not Collective
1811: Input Parameter:
1812: + x - the vector
1813: . m - first dimension of three dimensional array
1814: . n - second dimension of three dimensional array
1815: . p - third dimension of three dimensional array
1816: . mstart - first index you will use in first coordinate direction (often 0)
1817: . nstart - first index in the second coordinate direction (often 0)
1818: - pstart - first index in the third coordinate direction (often 0)
1820: Output Parameter:
1821: . a - location to put pointer to the array
1823: Level: developer
1825: Notes:
1826: For a vector obtained from DACreateLocalVector() mstart, nstart, and pstart are likely
1827: obtained from the corner indices obtained from DAGetGhostCorners() while for
1828: DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
1829: the arguments from DAGet[Ghost]Corners() are reversed in the call to VecGetArray3d().
1830:
1831: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
1833: Concepts: vector^accessing local values as 3d array
1835: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1836: VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1837: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1838: @*/
1839: PetscErrorCode VecGetArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
1840: {
1842: PetscInt i,N,j;
1843: PetscScalar *aa,**b;
1849: VecGetLocalSize(x,&N);
1850: if (m*n*p != N) SETERRQ4(PETSC_ERR_ARG_INCOMP,"Local array size %D does not match 3d array dimensions %D by %D by %D",N,m,n,p);
1851: VecGetArray(x,&aa);
1853: PetscMalloc(m*sizeof(PetscScalar**)+m*n*sizeof(PetscScalar*),a);
1854: b = (PetscScalar **)((*a) + m);
1855: for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
1856: for (i=0; i<m; i++) {
1857: for (j=0; j<n; j++) {
1858: b[i*n+j] = aa + i*n*p + j*p - pstart;
1859: }
1860: }
1861: *a -= mstart;
1862: return(0);
1863: }
1867: /*@C
1868: VecRestoreArray3d - Restores a vector after VecGetArray3d() has been called.
1870: Not Collective
1872: Input Parameters:
1873: + x - the vector
1874: . m - first dimension of three dimensional array
1875: . n - second dimension of the three dimensional array
1876: . p - third dimension of the three dimensional array
1877: . mstart - first index you will use in first coordinate direction (often 0)
1878: . nstart - first index in the second coordinate direction (often 0)
1879: . pstart - first index in the third coordinate direction (often 0)
1880: - a - location of pointer to array obtained from VecGetArray3d()
1882: Level: developer
1884: Notes:
1885: For regular PETSc vectors this routine does not involve any copies. For
1886: any special vectors that do not store local vector data in a contiguous
1887: array, this routine will copy the data back into the underlying
1888: vector data structure from the array obtained with VecGetArray().
1890: This routine actually zeros out the a pointer.
1892: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
1893: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
1894: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
1895: @*/
1896: PetscErrorCode VecRestoreArray3d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscScalar ***a[])
1897: {
1899: void *dummy;
1905: dummy = (void*)(*a + mstart);
1906: PetscFree(dummy);
1907: VecRestoreArray(x,PETSC_NULL);
1908: return(0);
1909: }
1913: /*@C
1914: VecGetArray4d - Returns a pointer to a 4d contiguous array that contains this
1915: processor's portion of the vector data. You MUST call VecRestoreArray4d()
1916: when you no longer need access to the array.
1918: Not Collective
1920: Input Parameter:
1921: + x - the vector
1922: . m - first dimension of four dimensional array
1923: . n - second dimension of four dimensional array
1924: . p - third dimension of four dimensional array
1925: . q - fourth dimension of four dimensional array
1926: . mstart - first index you will use in first coordinate direction (often 0)
1927: . nstart - first index in the second coordinate direction (often 0)
1928: . pstart - first index in the third coordinate direction (often 0)
1929: - qstart - first index in the fourth coordinate direction (often 0)
1931: Output Parameter:
1932: . a - location to put pointer to the array
1934: Level: beginner
1936: Notes:
1937: For a vector obtained from DACreateLocalVector() mstart, nstart, and pstart are likely
1938: obtained from the corner indices obtained from DAGetGhostCorners() while for
1939: DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
1940: the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray3d().
1941:
1942: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
1944: Concepts: vector^accessing local values as 3d array
1946: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
1947: VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
1948: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
1949: @*/
1950: PetscErrorCode VecGetArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
1951: {
1953: PetscInt i,N,j,k;
1954: PetscScalar *aa,***b,**c;
1960: VecGetLocalSize(x,&N);
1961: if (m*n*p*q != N) SETERRQ5(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);
1962: VecGetArray(x,&aa);
1964: PetscMalloc(m*sizeof(PetscScalar***)+m*n*sizeof(PetscScalar**)+m*n*p*sizeof(PetscScalar*),a);
1965: b = (PetscScalar ***)((*a) + m);
1966: c = (PetscScalar **)(b + m*n);
1967: for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
1968: for (i=0; i<m; i++) {
1969: for (j=0; j<n; j++) {
1970: b[i*n+j] = c + i*n*p + j*p - pstart;
1971: }
1972: }
1973: for (i=0; i<m; i++) {
1974: for (j=0; j<n; j++) {
1975: for (k=0; k<p; k++) {
1976: c[i*n*p+j*p+k] = aa + i*n*p*q + j*p*q + k*q - qstart;
1977: }
1978: }
1979: }
1980: *a -= mstart;
1981: return(0);
1982: }
1986: /*@C
1987: VecRestoreArray4d - Restores a vector after VecGetArray3d() has been called.
1989: Not Collective
1991: Input Parameters:
1992: + x - the vector
1993: . m - first dimension of four dimensional array
1994: . n - second dimension of the four dimensional array
1995: . p - third dimension of the four dimensional array
1996: . q - fourth dimension of the four dimensional array
1997: . mstart - first index you will use in first coordinate direction (often 0)
1998: . nstart - first index in the second coordinate direction (often 0)
1999: . pstart - first index in the third coordinate direction (often 0)
2000: . qstart - first index in the fourth coordinate direction (often 0)
2001: - a - location of pointer to array obtained from VecGetArray4d()
2003: Level: beginner
2005: Notes:
2006: For regular PETSc vectors this routine does not involve any copies. For
2007: any special vectors that do not store local vector data in a contiguous
2008: array, this routine will copy the data back into the underlying
2009: vector data structure from the array obtained with VecGetArray().
2011: This routine actually zeros out the a pointer.
2013: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2014: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2015: VecGetArray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2016: @*/
2017: PetscErrorCode VecRestoreArray4d(Vec x,PetscInt m,PetscInt n,PetscInt p,PetscInt q,PetscInt mstart,PetscInt nstart,PetscInt pstart,PetscInt qstart,PetscScalar ****a[])
2018: {
2020: void *dummy;
2026: dummy = (void*)(*a + mstart);
2027: PetscFree(dummy);
2028: VecRestoreArray(x,PETSC_NULL);
2029: return(0);
2030: }