Actual source code: vector.c
1: /*$Id: vector.c,v 1.229 2001/03/28 19:41:06 balay Exp $*/
2: /*
3: Provides the interface functions for all vector operations.
4: These are the vector functions the user calls.
5: */
6: #include "src/vec/vecimpl.h" /*I "petscvec.h" I*/
8: /*@
9: VecSetBlockSize - Sets the blocksize for future calls to VecSetValuesBlocked()
10: and VecSetValuesBlockedLocal().
12: Collective on Vec
14: Input Parameter:
15: + v - the vector
16: - bs - the blocksize
18: Notes:
19: All vectors obtained by VecDuplicate() inherit the same blocksize.
21: Level: advanced
23: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMappingBlocked(), VecGetBlockSize()
25: Concepts: block size^vectors
26: @*/
27: int VecSetBlockSize(Vec v,int bs)
28: {
31: if (bs <= 0) bs = 1;
32: if (bs == v->bs) return(0);
33: if (v->bs != -1) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot reset blocksize. Current size %d new %d",v->bs,bs);
34: if (v->N % bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Vector length not divisible by blocksize %d %d",v->N,bs);
35: if (v->n % bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Local vector length not divisible by blocksize %d %d",v->n,bs);
36:
37: v->bs = bs;
38: v->bstash.bs = bs; /* use the same blocksize for the vec's block-stash */
39: return(0);
40: }
42: /*@
43: VecGetBlockSize - Gets the blocksize for the vector, i.e. what is used for VecSetValuesBlocked()
44: and VecSetValuesBlockedLocal().
46: Collective on Vec
48: Input Parameter:
49: . v - the vector
51: Output Parameter:
52: . bs - the blocksize
54: Notes:
55: All vectors obtained by VecDuplicate() inherit the same blocksize.
57: Level: advanced
59: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMappingBlocked(), VecSetBlockSize()
61: Concepts: vector^block size
62: Concepts: block^vector
64: @*/
65: int VecGetBlockSize(Vec v,int *bs)
66: {
69: *bs = v->bs;
70: return(0);
71: }
73: /*@
74: VecValid - Checks whether a vector object is valid.
76: Not Collective
78: Input Parameter:
79: . v - the object to check
81: Output Parameter:
82: flg - flag indicating vector status, either
83: PETSC_TRUE if vector is valid, or PETSC_FALSE otherwise.
85: Level: developer
87: @*/
88: int VecValid(Vec v,PetscTruth *flg)
89: {
92: if (!v) *flg = PETSC_FALSE;
93: else if (v->cookie != VEC_COOKIE) *flg = PETSC_FALSE;
94: else *flg = PETSC_TRUE;
95: return(0);
96: }
98: /*@
99: VecDot - Computes the vector dot product.
101: Collective on Vec
103: Input Parameters:
104: . x, y - the vectors
106: Output Parameter:
107: . alpha - the dot product
109: Performance Issues:
110: + per-processor memory bandwidth
111: . interprocessor latency
112: - work load inbalance that causes certain processes to arrive much earlier than
113: others
115: Notes for Users of Complex Numbers:
116: For complex vectors, VecDot() computes
117: $ val = (x,y) = y^H x,
118: where y^H denotes the conjugate transpose of y.
120: Use VecTDot() for the indefinite form
121: $ val = (x,y) = y^T x,
122: where y^T denotes the transpose of y.
124: Level: intermediate
126: Concepts: inner product
127: Concepts: vector^inner product
129: .seealso: VecMDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd()
130: @*/
131: int VecDot(Vec x,Vec y,Scalar *val)
132: {
143: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
144: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
146: PetscLogEventBarrierBegin(VEC_DotBarrier,x,y,0,0,x->comm);
147: (*x->ops->dot)(x,y,val);
148: PetscLogEventBarrierEnd(VEC_DotBarrier,x,y,0,0,x->comm);
149: /*
150: The next block is for incremental debugging
151: */
152: if (PetscCompare) {
153: int flag;
154: MPI_Comm_compare(PETSC_COMM_WORLD,x->comm,&flag);
155: if (flag != MPI_UNEQUAL) {
156: PetscCompareScalar(*val);
157: }
158: }
159: return(0);
160: }
162: /*@
163: VecNorm - Computes the vector norm.
165: Collective on Vec
167: Input Parameters:
168: + x - the vector
169: - type - one of NORM_1, NORM_2, NORM_INFINITY. Also available
170: NORM_1_AND_2, which computes both norms and stores them
171: in a two element array.
173: Output Parameter:
174: . val - the norm
176: Notes:
177: $ NORM_1 denotes sum_i |x_i|
178: $ NORM_2 denotes sqrt(sum_i (x_i)^2)
179: $ NORM_INFINITY denotes max_i |x_i|
181: Level: intermediate
183: Performance Issues:
184: + per-processor memory bandwidth
185: . interprocessor latency
186: - work load inbalance that causes certain processes to arrive much earlier than
187: others
189: Compile Option:
190: PETSC_HAVE_SLOW_NRM2 will cause a C (loop unrolled) version of the norm to be used, rather
191: than the BLAS. This should probably only be used when one is using the FORTRAN BLAS routines
192: (as opposed to vendor provided) because the FORTRAN BLAS NRM2() routine is very slow.
194: Concepts: norm
195: Concepts: vector^norm
197: .seealso: VecDot(), VecTDot(), VecNorm(), VecDotBegin(), VecDotEnd(),
198: VecNormBegin(), VecNormEnd()
200: @*/
201: int VecNorm(Vec x,NormType type,PetscReal *val)
202: {
208: PetscLogEventBarrierBegin(VEC_NormBarrier,x,0,0,0,x->comm);
209: (*x->ops->norm)(x,type,val);
210: PetscLogEventBarrierEnd(VEC_NormBarrier,x,0,0,0,x->comm);
211: /*
212: The next block is for incremental debugging
213: */
214: if (PetscCompare) {
215: int flag;
216: MPI_Comm_compare(PETSC_COMM_WORLD,x->comm,&flag);
217: if (flag != MPI_UNEQUAL) {
218: PetscCompareDouble(*val);
219: }
220: }
221: return(0);
222: }
224: /*@C
225: VecMax - Determines the maximum vector component and its location.
227: Collective on Vec
229: Input Parameter:
230: . x - the vector
232: Output Parameters:
233: + val - the maximum component
234: - p - the location of val
236: Notes:
237: Returns the value PETSC_MIN and p = -1 if the vector is of length 0.
239: Level: intermediate
241: Concepts: maximum^of vector
242: Concepts: vector^maximum value
244: .seealso: VecNorm(), VecMin()
245: @*/
246: int VecMax(Vec x,int *p,PetscReal *val)
247: {
254: PetscLogEventBegin(VEC_Max,x,0,0,0);
255: (*x->ops->max)(x,p,val);
256: PetscLogEventEnd(VEC_Max,x,0,0,0);
257: return(0);
258: }
260: /*@
261: VecMin - Determines the minimum vector component and its location.
263: Collective on Vec
265: Input Parameters:
266: . x - the vector
268: Output Parameter:
269: + val - the minimum component
270: - p - the location of val
272: Level: intermediate
274: Notes:
275: Returns the value PETSC_MAX and p = -1 if the vector is of length 0.
277: Concepts: minimum^of vector
278: Concepts: vector^minimum entry
280: .seealso: VecMax()
281: @*/
282: int VecMin(Vec x,int *p,PetscReal *val)
283: {
290: PetscLogEventBegin(VEC_Min,x,0,0,0);
291: (*x->ops->min)(x,p,val);
292: PetscLogEventEnd(VEC_Min,x,0,0,0);
293: return(0);
294: }
296: /*@
297: VecTDot - Computes an indefinite vector dot product. That is, this
298: routine does NOT use the complex conjugate.
300: Collective on Vec
302: Input Parameters:
303: . x, y - the vectors
305: Output Parameter:
306: . val - the dot product
308: Notes for Users of Complex Numbers:
309: For complex vectors, VecTDot() computes the indefinite form
310: $ val = (x,y) = y^T x,
311: where y^T denotes the transpose of y.
313: Use VecDot() for the inner product
314: $ val = (x,y) = y^H x,
315: where y^H denotes the conjugate transpose of y.
317: Level: intermediate
319: Concepts: inner product^non-Hermitian
320: Concepts: vector^inner product
321: Concepts: non-Hermitian inner product
323: .seealso: VecDot(), VecMTDot()
324: @*/
325: int VecTDot(Vec x,Vec y,Scalar *val)
326: {
337: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
338: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
340: PetscLogEventBegin(VEC_TDot,x,y,0,0);
341: (*x->ops->tdot)(x,y,val);
342: PetscLogEventEnd(VEC_TDot,x,y,0,0);
343: return(0);
344: }
346: /*@
347: VecScale - Scales a vector.
349: Collective on Vec
351: Input Parameters:
352: + x - the vector
353: - alpha - the scalar
355: Output Parameter:
356: . x - the scaled vector
358: Note:
359: For a vector with n components, VecScale() computes
360: $ x[i] = alpha * x[i], for i=1,...,n.
362: Level: intermediate
364: Concepts: vector^scaling
365: Concepts: scaling^vector
367: @*/
368: int VecScale(const Scalar *alpha,Vec x)
369: {
376: PetscLogEventBegin(VEC_Scale,x,0,0,0);
377: (*x->ops->scale)(alpha,x);
378: PetscLogEventEnd(VEC_Scale,x,0,0,0);
379: return(0);
380: }
382: /*@
383: VecCopy - Copies a vector.
385: Collective on Vec
387: Input Parameter:
388: . x - the vector
390: Output Parameter:
391: . y - the copy
393: Notes:
394: For default parallel PETSc vectors, both x and y must be distributed in
395: the same manner; local copies are done.
397: Level: beginner
399: .seealso: VecDuplicate()
400: @*/
401: int VecCopy(Vec x,Vec y)
402: {
411: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
412: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
414: PetscLogEventBegin(VEC_Copy,x,y,0,0);
415: (*x->ops->copy)(x,y);
416: PetscLogEventEnd(VEC_Copy,x,y,0,0);
417: return(0);
418: }
420: /*@
421: VecSet - Sets all components of a vector to a single scalar value.
423: Collective on Vec
425: Input Parameters:
426: + alpha - the scalar
427: - x - the vector
429: Output Parameter:
430: . x - the vector
432: Note:
433: For a vector of dimension n, VecSet() computes
434: $ x[i] = alpha, for i=1,...,n,
435: so that all vector entries then equal the identical
436: scalar value, alpha. Use the more general routine
437: VecSetValues() to set different vector entries.
439: Level: beginner
441: .seealso VecSetValues(), VecSetValuesBlocked(), VecSetRandom()
443: Concepts: vector^setting to constant
445: @*/
446: int VecSet(const Scalar *alpha,Vec x)
447: {
455: PetscLogEventBegin(VEC_Set,x,0,0,0);
456: (*x->ops->set)(alpha,x);
457: PetscLogEventEnd(VEC_Set,x,0,0,0);
458: return(0);
459: }
461: /*@C
462: VecSetRandom - Sets all components of a vector to random numbers.
464: Collective on Vec
466: Input Parameters:
467: + rctx - the random number context, formed by PetscRandomCreate(), or PETSC_NULL and
468: it will create one internally.
469: - x - the vector
471: Output Parameter:
472: . x - the vector
474: Example of Usage:
475: .vb
476: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
477: VecSetRandom(rctx,x);
478: PetscRandomDestroy(rctx);
479: .ve
481: Level: intermediate
483: Concepts: vector^setting to random
484: Concepts: random^vector
486: .seealso: VecSet(), VecSetValues(), PetscRandomCreate(), PetscRandomDestroy()
487: @*/
488: int VecSetRandom(PetscRandom rctx,Vec x)
489: {
490: int ierr;
491: PetscRandom rand = 0;
498: if (!rctx) {
499: MPI_Comm comm;
500: PetscObjectGetComm((PetscObject)x,&comm);
501: PetscRandomCreate(comm,RANDOM_DEFAULT,&rand);
502: rctx = rand;
503: }
505: PetscLogEventBegin(VEC_SetRandom,x,rctx,0,0);
506: (*x->ops->setrandom)(rctx,x);
507: PetscLogEventEnd(VEC_SetRandom,x,rctx,0,0);
508:
509: if (rand) {
510: PetscRandomDestroy(rand);
511: }
512: return(0);
513: }
515: /*@
516: VecAXPY - Computes y = alpha x + y.
518: Collective on Vec
520: Input Parameters:
521: + alpha - the scalar
522: - x, y - the vectors
524: Output Parameter:
525: . y - output vector
527: Level: intermediate
529: Concepts: vector^BLAS
530: Concepts: BLAS
532: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY()
533: @*/
534: int VecAXPY(const Scalar *alpha,Vec x,Vec y)
535: {
546: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
547: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
549: PetscLogEventBegin(VEC_AXPY,x,y,0,0);
550: (*x->ops->axpy)(alpha,x,y);
551: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
552: return(0);
553: }
555: /*@
556: VecAXPBY - Computes y = alpha x + beta y.
558: Collective on Vec
560: Input Parameters:
561: + alpha,beta - the scalars
562: - x, y - the vectors
564: Output Parameter:
565: . y - output vector
567: Level: intermediate
569: Concepts: BLAS
570: Concepts: vector^BLAS
572: .seealso: VecAYPX(), VecMAXPY(), VecWAXPY(), VecAXPY()
573: @*/
574: int VecAXPBY(const Scalar *alpha,const Scalar *beta,Vec x,Vec y)
575: {
587: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
588: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
590: PetscLogEventBegin(VEC_AXPY,x,y,0,0);
591: (*x->ops->axpby)(alpha,beta,x,y);
592: PetscLogEventEnd(VEC_AXPY,x,y,0,0);
593: return(0);
594: }
596: /*@
597: VecAYPX - Computes y = x + alpha y.
599: Collective on Vec
601: Input Parameters:
602: + alpha - the scalar
603: - x, y - the vectors
605: Output Parameter:
606: . y - output vector
608: Level: intermediate
610: Concepts: vector^BLAS
611: Concepts: BLAS
613: .seealso: VecAXPY(), VecWAXPY()
614: @*/
615: int VecAYPX(const Scalar *alpha,Vec x,Vec y)
616: {
627: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
628: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
630: PetscLogEventBegin(VEC_AYPX,x,y,0,0);
631: (*x->ops->aypx)(alpha,x,y);
632: PetscLogEventEnd(VEC_AYPX,x,y,0,0);
633: return(0);
634: }
636: /*@
637: VecSwap - Swaps the vectors x and y.
639: Collective on Vec
641: Input Parameters:
642: . x, y - the vectors
644: Level: advanced
646: Concepts: vector^swapping values
648: @*/
649: int VecSwap(Vec x,Vec y)
650: {
660: if (x->N != y->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
661: if (x->n != y->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
663: PetscLogEventBegin(VEC_Swap,x,y,0,0);
664: (*x->ops->swap)(x,y);
665: PetscLogEventEnd(VEC_Swap,x,y,0,0);
666: return(0);
667: }
669: /*@
670: VecWAXPY - Computes w = alpha x + y.
672: Collective on Vec
674: Input Parameters:
675: + alpha - the scalar
676: - x, y - the vectors
678: Output Parameter:
679: . w - the result
681: Level: intermediate
683: Concepts: vector^BLAS
684: Concepts: BLAS
686: .seealso: VecAXPY(), VecAYPX()
687: @*/
688: int VecWAXPY(const Scalar *alpha,Vec x,Vec y,Vec w)
689: {
702: if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
703: if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
705: PetscLogEventBegin(VEC_WAXPY,x,y,w,0);
706: (*x->ops->waxpy)(alpha,x,y,w);
707: PetscLogEventEnd(VEC_WAXPY,x,y,w,0);
708: return(0);
709: }
711: /*@
712: VecPointwiseMult - Computes the componentwise multiplication w = x*y.
714: Collective on Vec
716: Input Parameters:
717: . x, y - the vectors
719: Output Parameter:
720: . w - the result
722: Level: advanced
724: Notes: any subset of the x, y, and w may be the same vector.
726: Concepts: vector^pointwise multiply
728: .seealso: VecPointwiseDivide()
729: @*/
730: int VecPointwiseMult(Vec x,Vec y,Vec w)
731: {
743: if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
744: if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
746: PetscLogEventBegin(VEC_PMult,x,y,w,0);
747: (*x->ops->pointwisemult)(x,y,w);
748: PetscLogEventEnd(VEC_PMult,x,y,w,0);
749: return(0);
750: }
752: /*@
753: VecPointwiseDivide - Computes the componentwise division w = x/y.
755: Collective on Vec
757: Input Parameters:
758: . x, y - the vectors
760: Output Parameter:
761: . w - the result
763: Level: advanced
765: Notes: any subset of the x, y, and w may be the same vector.
767: Concepts: vector^pointwise divide
769: .seealso: VecPointwiseMult()
770: @*/
771: int VecPointwiseDivide(Vec x,Vec y,Vec w)
772: {
784: if (x->N != y->N || x->N != w->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
785: if (x->n != y->n || x->n != w->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
787: (*x->ops->pointwisedivide)(x,y,w);
788: return(0);
789: }
791: /*@C
792: VecDuplicate - Creates a new vector of the same type as an existing vector.
794: Collective on Vec
796: Input Parameters:
797: . v - a vector to mimic
799: Output Parameter:
800: . newv - location to put new vector
802: Notes:
803: VecDuplicate() does not copy the vector, but rather allocates storage
804: for the new vector. Use VecCopy() to copy a vector.
806: Use VecDestroy() to free the space. Use VecDuplicateVecs() to get several
807: vectors.
809: Level: beginner
811: .seealso: VecDestroy(), VecDuplicateVecs(), VecCreate(), VecCopy()
812: @*/
813: int VecDuplicate(Vec x,Vec *newv)
814: {
821: (*x->ops->duplicate)(x,newv);
822: return(0);
823: }
825: /*@C
826: VecDestroy - Destroys a vector.
828: Collective on Vec
830: Input Parameters:
831: . v - the vector
833: Level: beginner
835: .seealso: VecDuplicate()
836: @*/
837: int VecDestroy(Vec v)
838: {
843: if (--v->refct > 0) return(0);
844: /* destroy the internal part */
845: if (v->ops->destroy) {
846: (*v->ops->destroy)(v);
847: }
848: /* destroy the external/common part */
849: if (v->mapping) {
850: ISLocalToGlobalMappingDestroy(v->mapping);
851: }
852: if (v->bmapping) {
853: ISLocalToGlobalMappingDestroy(v->bmapping);
854: }
855: if (v->map) {
856: MapDestroy(v->map);
857: }
858: PetscLogObjectDestroy(v);
859: PetscHeaderDestroy(v);
860: return(0);
861: }
863: /*@C
864: VecDuplicateVecs - Creates several vectors of the same type as an existing vector.
866: Collective on Vec
868: Input Parameters:
869: + m - the number of vectors to obtain
870: - v - a vector to mimic
872: Output Parameter:
873: . V - location to put pointer to array of vectors
875: Notes:
876: Use VecDestroyVecs() to free the space. Use VecDuplicate() to form a single
877: vector.
879: Fortran Note:
880: The Fortran interface is slightly different from that given below, it
881: requires one to pass in V a Vec (integer) array of size at least m.
882: See the Fortran chapter of the users manual and petsc/src/vec/examples for details.
884: Level: intermediate
886: .seealso: VecDestroyVecs(), VecDuplicate(), VecCreate(), VecDuplicateVecsF90()
887: @*/
888: int VecDuplicateVecs(Vec v,int m,Vec *V[])
889: {
896: (*v->ops->duplicatevecs)(v, m,V);
897: return(0);
898: }
900: /*@C
901: VecDestroyVecs - Frees a block of vectors obtained with VecDuplicateVecs().
903: Collective on Vec
905: Input Parameters:
906: + vv - pointer to array of vector pointers
907: - m - the number of vectors previously obtained
909: Fortran Note:
910: The Fortran interface is slightly different from that given below.
911: See the Fortran chapter of the users manual and
912: petsc/src/vec/examples for details.
914: Level: intermediate
916: .seealso: VecDuplicateVecs(), VecDestroyVecsF90()
917: @*/
918: int VecDestroyVecs(const Vec vv[],int m)
919: {
923: if (!vv) SETERRQ(PETSC_ERR_ARG_BADPTR,"Null vectors");
926: (*(*vv)->ops->destroyvecs)(vv,m);
927: return(0);
928: }
930: /*@
931: VecSetValues - Inserts or adds values into certain locations of a vector.
933: Input Parameters:
934: Not Collective
936: + x - vector to insert in
937: . ni - number of elements to add
938: . ix - indices where to add
939: . y - array of values
940: - iora - either INSERT_VALUES or ADD_VALUES, where
941: ADD_VALUES adds values to any existing entries, and
942: INSERT_VALUES replaces existing entries with new values
944: Notes:
945: VecSetValues() sets x[ix[i]] = y[i], for i=0,...,ni-1.
947: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
948: options cannot be mixed without intervening calls to the assembly
949: routines.
951: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
952: MUST be called after all calls to VecSetValues() have been completed.
954: VecSetValues() uses 0-based indices in Fortran as well as in C.
956: Negative indices may be passed in ix, these rows are
957: simply ignored. This allows easily inserting element load matrices
958: with homogeneous Dirchlet boundary conditions that you don't want represented
959: in the vector.
961: Level: beginner
963: Concepts: vector^setting values
965: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesLocal(),
966: VecSetValue(), VecSetValuesBlocked()
967: @*/
968: int VecSetValues(Vec x,int ni,const int ix[],const Scalar y[],InsertMode iora)
969: {
977: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
978: (*x->ops->setvalues)(x,ni,ix,y,iora);
979: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
980: return(0);
981: }
983: /*@
984: VecSetValuesBlocked - Inserts or adds blocks of values into certain locations of a vector.
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, rather than 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: Notes:
998: VecSetValuesBlocked() sets x[ix[bs*i]+j] = y[bs*i+j],
999: for j=0,...,bs, for i=0,...,ni-1. where bs was set with VecSetBlockSize().
1001: Calls to VecSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES
1002: options cannot be mixed without intervening calls to the assembly
1003: routines.
1005: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1006: MUST be called after all calls to VecSetValuesBlocked() have been completed.
1008: VecSetValuesBlocked() uses 0-based indices in Fortran as well as in C.
1010: Negative indices may be passed in ix, these rows are
1011: simply ignored. This allows easily inserting element load matrices
1012: with homogeneous Dirchlet boundary conditions that you don't want represented
1013: in the vector.
1015: Level: intermediate
1017: Concepts: vector^setting values blocked
1019: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(),
1020: VecSetValues()
1021: @*/
1022: int VecSetValuesBlocked(Vec x,int ni,const int ix[],const Scalar y[],InsertMode iora)
1023: {
1031: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1032: (*x->ops->setvaluesblocked)(x,ni,ix,y,iora);
1033: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1034: return(0);
1035: }
1037: /*MC
1038: VecSetValue - Set a single entry into a vector.
1040: Synopsis:
1041: void VecSetValue(Vec v,int row,Scalar value, InsertMode mode);
1043: Not Collective
1045: Input Parameters:
1046: + v - the vector
1047: . row - the row location of the entry
1048: . value - the value to insert
1049: - mode - either INSERT_VALUES or ADD_VALUES
1051: Notes:
1052: For efficiency one should use VecSetValues() and set several or
1053: many values simultaneously if possible.
1055: Note that VecSetValue() does NOT return an error code (since this
1056: is checked internally).
1058: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1059: MUST be called after all calls to VecSetValues() have been completed.
1061: VecSetValues() uses 0-based indices in Fortran as well as in C.
1063: Level: beginner
1065: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal()
1066: M*/
1068: /*@
1069: VecSetLocalToGlobalMapping - Sets a local numbering to global numbering used
1070: by the routine VecSetValuesLocal() to allow users to insert vector entries
1071: using a local (per-processor) numbering.
1073: Collective on Vec
1075: Input Parameters:
1076: + x - vector
1077: - mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()
1079: Notes:
1080: All vectors obtained with VecDuplicate() from this vector inherit the same mapping.
1082: Level: intermediate
1084: Concepts: vector^setting values with local numbering
1086: seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1087: VecSetLocalToGlobalMappingBlocked(), VecSetValuesBlockedLocal()
1088: @*/
1089: int VecSetLocalToGlobalMapping(Vec x,ISLocalToGlobalMapping mapping)
1090: {
1096: if (x->mapping) {
1097: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1098: }
1100: if (x->ops->setlocaltoglobalmapping) {
1101: (*x->ops->setlocaltoglobalmapping)(x,mapping);
1102: } else {
1103: x->mapping = mapping;
1104: PetscObjectReference((PetscObject)mapping);
1105: }
1106: return(0);
1107: }
1109: /*@
1110: VecSetLocalToGlobalMappingBlock - Sets a local numbering to global numbering used
1111: by the routine VecSetValuesBlockedLocal() to allow users to insert vector entries
1112: using a local (per-processor) numbering.
1114: Collective on Vec
1116: Input Parameters:
1117: + x - vector
1118: - mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()
1120: Notes:
1121: All vectors obtained with VecDuplicate() from this vector inherit the same mapping.
1123: Level: intermediate
1125: Concepts: vector^setting values blocked with local numbering
1127: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
1128: VecSetLocalToGlobalMapping(), VecSetValuesBlockedLocal()
1129: @*/
1130: int VecSetLocalToGlobalMappingBlock(Vec x,ISLocalToGlobalMapping mapping)
1131: {
1137: if (x->bmapping) {
1138: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Mapping already set for vector");
1139: }
1140: x->bmapping = mapping;
1141: PetscObjectReference((PetscObject)mapping);
1142: return(0);
1143: }
1145: /*@
1146: VecSetValuesLocal - Inserts or adds values into certain locations of a vector,
1147: using a local ordering of the nodes.
1149: Not Collective
1151: Input Parameters:
1152: + x - vector to insert in
1153: . ni - number of elements to add
1154: . ix - indices where to add
1155: . y - array of values
1156: - iora - either INSERT_VALUES or ADD_VALUES, where
1157: ADD_VALUES adds values to any existing entries, and
1158: INSERT_VALUES replaces existing entries with new values
1160: Level: intermediate
1162: Notes:
1163: VecSetValuesLocal() sets x[ix[i]] = y[i], for i=0,...,ni-1.
1165: Calls to VecSetValues() with the INSERT_VALUES and ADD_VALUES
1166: options cannot be mixed without intervening calls to the assembly
1167: routines.
1169: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1170: MUST be called after all calls to VecSetValuesLocal() have been completed.
1172: VecSetValuesLocal() uses 0-based indices in Fortran as well as in C.
1174: Concepts: vector^setting values with local numbering
1176: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetLocalToGlobalMapping(),
1177: VecSetValuesBlockedLocal()
1178: @*/
1179: int VecSetValuesLocal(Vec x,int ni,const int ix[],const Scalar y[],InsertMode iora)
1180: {
1181: int ierr,lixp[128],*lix = lixp;
1189: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1190: if (!x->ops->setvalueslocal) {
1191: if (!x->mapping) {
1192: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMapping()");
1193: }
1194: if (ni > 128) {
1195: PetscMalloc(ni*sizeof(int),&lix);
1196: }
1197: ISLocalToGlobalMappingApply(x->mapping,ni,(int*)ix,lix);
1198: (*x->ops->setvalues)(x,ni,lix,y,iora);
1199: if (ni > 128) {
1200: PetscFree(lix);
1201: }
1202: } else {
1203: (*x->ops->setvalueslocal)(x,ni,ix,y,iora);
1204: }
1205: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1206: return(0);
1207: }
1209: /*@
1210: VecSetValuesBlockedLocal - Inserts or adds values into certain locations of a vector,
1211: using a local ordering of the nodes.
1213: Not Collective
1215: Input Parameters:
1216: + x - vector to insert in
1217: . ni - number of blocks to add
1218: . ix - indices where to add in block count, not element count
1219: . y - array of values
1220: - iora - either INSERT_VALUES or ADD_VALUES, where
1221: ADD_VALUES adds values to any existing entries, and
1222: INSERT_VALUES replaces existing entries with new values
1224: Level: intermediate
1226: Notes:
1227: VecSetValuesBlockedLocal() sets x[bs*ix[i]+j] = y[bs*i+j],
1228: for j=0,..bs-1, for i=0,...,ni-1, where bs has been set with VecSetBlockSize().
1230: Calls to VecSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES
1231: options cannot be mixed without intervening calls to the assembly
1232: routines.
1234: These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd()
1235: MUST be called after all calls to VecSetValuesBlockedLocal() have been completed.
1237: VecSetValuesBlockedLocal() uses 0-based indices in Fortran as well as in C.
1240: Concepts: vector^setting values blocked with local numbering
1242: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesBlocked(),
1243: VecSetLocalToGlobalMappingBlocked()
1244: @*/
1245: int VecSetValuesBlockedLocal(Vec x,int ni,const int ix[],const Scalar y[],InsertMode iora)
1246: {
1247: int ierr,lixp[128],*lix = lixp;
1254: if (!x->bmapping) {
1255: SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Local to global never set with VecSetLocalToGlobalMappingBlocked()");
1256: }
1257: if (ni > 128) {
1258: PetscMalloc(ni*sizeof(int),&lix);
1259: }
1261: PetscLogEventBegin(VEC_SetValues,x,0,0,0);
1262: ISLocalToGlobalMappingApply(x->bmapping,ni,(int*)ix,lix);
1263: (*x->ops->setvaluesblocked)(x,ni,lix,y,iora);
1264: PetscLogEventEnd(VEC_SetValues,x,0,0,0);
1265: if (ni > 128) {
1266: PetscFree(lix);
1267: }
1268: return(0);
1269: }
1271: /*@
1272: VecAssemblyBegin - Begins assembling the vector. This routine should
1273: be called after completing all calls to VecSetValues().
1275: Collective on Vec
1277: Input Parameter:
1278: . vec - the vector
1280: Level: beginner
1282: Concepts: assembly^vectors
1284: .seealso: VecAssemblyEnd(), VecSetValues()
1285: @*/
1286: int VecAssemblyBegin(Vec vec)
1287: {
1288: int ierr;
1289: PetscTruth flg;
1295: PetscOptionsHasName(vec->prefix,"-vec_view_stash",&flg);
1296: if (flg) {
1297: VecStashView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1298: }
1300: PetscLogEventBegin(VEC_AssemblyBegin,vec,0,0,0);
1301: if (vec->ops->assemblybegin) {
1302: (*vec->ops->assemblybegin)(vec);
1303: }
1304: PetscLogEventEnd(VEC_AssemblyBegin,vec,0,0,0);
1305: return(0);
1306: }
1308: /*@
1309: VecAssemblyEnd - Completes assembling the vector. This routine should
1310: be called after VecAssemblyBegin().
1312: Collective on Vec
1314: Input Parameter:
1315: . vec - the vector
1317: Options Database Keys:
1318: . -vec_view - Prints vector in ASCII format
1319: . -vec_view_matlab - Prints vector in Matlab format
1320: . -vec_view_draw - Activates vector viewing using drawing tools
1321: . -display <name> - Sets display name (default is host)
1322: . -draw_pause <sec> - Sets number of seconds to pause after display
1323: . -vec_view_socket - Activates vector viewing using a socket
1324: - -vec_view_ams - Activates vector viewing using the ALICE Memory Snooper (AMS)
1325:
1326: Level: beginner
1328: .seealso: VecAssemblyBegin(), VecSetValues()
1329: @*/
1330: int VecAssemblyEnd(Vec vec)
1331: {
1332: int ierr;
1333: PetscTruth flg;
1337: PetscLogEventBegin(VEC_AssemblyEnd,vec,0,0,0);
1339: if (vec->ops->assemblyend) {
1340: (*vec->ops->assemblyend)(vec);
1341: }
1342: PetscLogEventEnd(VEC_AssemblyEnd,vec,0,0,0);
1343: PetscOptionsHasName(vec->prefix,"-vec_view",&flg);
1344: if (flg) {
1345: VecView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1346: }
1347: PetscOptionsHasName(PETSC_NULL,"-vec_view_matlab",&flg);
1348: if (flg) {
1349: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_(vec->comm),PETSC_VIEWER_ASCII_MATLAB);
1350: VecView(vec,PETSC_VIEWER_STDOUT_(vec->comm));
1351: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_(vec->comm));
1352: }
1353: PetscOptionsHasName(PETSC_NULL,"-vec_view_draw",&flg);
1354: if (flg) {
1355: VecView(vec,PETSC_VIEWER_DRAW_(vec->comm));
1356: PetscViewerFlush(PETSC_VIEWER_DRAW_(vec->comm));
1357: }
1358: PetscOptionsHasName(PETSC_NULL,"-vec_view_draw_lg",&flg);
1359: if (flg) {
1360: PetscViewerSetFormat(PETSC_VIEWER_DRAW_(vec->comm),PETSC_VIEWER_DRAW_LG);
1361: VecView(vec,PETSC_VIEWER_DRAW_(vec->comm));
1362: PetscViewerFlush(PETSC_VIEWER_DRAW_(vec->comm));
1363: }
1364: PetscOptionsHasName(PETSC_NULL,"-vec_view_socket",&flg);
1365: if (flg) {
1366: VecView(vec,PETSC_VIEWER_SOCKET_(vec->comm));
1367: PetscViewerFlush(PETSC_VIEWER_SOCKET_(vec->comm));
1368: }
1369: PetscOptionsHasName(PETSC_NULL,"-vec_view_binary",&flg);
1370: if (flg) {
1371: VecView(vec,PETSC_VIEWER_BINARY_(vec->comm));
1372: PetscViewerFlush(PETSC_VIEWER_BINARY_(vec->comm));
1373: }
1374: #if defined(PETSC_HAVE_AMS)
1375: PetscOptionsHasName(PETSC_NULL,"-vec_view_ams",&flg);
1376: if (flg) {
1377: VecView(vec,PETSC_VIEWER_AMS_(vec->comm));
1378: }
1379: #endif
1380: return(0);
1381: }
1384: /*@C
1385: VecMTDot - Computes indefinite vector multiple dot products.
1386: That is, it does NOT use the complex conjugate.
1388: Collective on Vec
1390: Input Parameters:
1391: + nv - number of vectors
1392: . x - one vector
1393: - y - array of vectors. Note that vectors are pointers
1395: Output Parameter:
1396: . val - array of the dot products
1398: Notes for Users of Complex Numbers:
1399: For complex vectors, VecMTDot() computes the indefinite form
1400: $ val = (x,y) = y^T x,
1401: where y^T denotes the transpose of y.
1403: Use VecMDot() for the inner product
1404: $ val = (x,y) = y^H x,
1405: where y^H denotes the conjugate transpose of y.
1407: Level: intermediate
1409: Concepts: inner product^multiple
1410: Concepts: vector^multiple inner products
1412: .seealso: VecMDot(), VecTDot()
1413: @*/
1414: int VecMTDot(int nv,Vec x,const Vec y[],Scalar *val)
1415: {
1426: if (x->N != (*y)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1427: if (x->n != (*y)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
1429: PetscLogEventBegin(VEC_MTDot,x,*y,0,0);
1430: (*x->ops->mtdot)(nv,x,y,val);
1431: PetscLogEventEnd(VEC_MTDot,x,*y,0,0);
1432: return(0);
1433: }
1435: /*@C
1436: VecMDot - Computes vector multiple dot products.
1438: Collective on Vec
1440: Input Parameters:
1441: + nv - number of vectors
1442: . x - one vector
1443: - y - array of vectors.
1445: Output Parameter:
1446: . val - array of the dot products
1448: Notes for Users of Complex Numbers:
1449: For complex vectors, VecMDot() computes
1450: $ val = (x,y) = y^H x,
1451: where y^H denotes the conjugate transpose of y.
1453: Use VecMTDot() for the indefinite form
1454: $ val = (x,y) = y^T x,
1455: where y^T denotes the transpose of y.
1457: Level: intermediate
1459: Concepts: inner product^multiple
1460: Concepts: vector^multiple inner products
1462: .seealso: VecMTDot(), VecDot()
1463: @*/
1464: int VecMDot(int nv,Vec x,const Vec y[],Scalar *val)
1465: {
1476: if (x->N != (*y)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1477: if (x->n != (*y)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
1479: PetscLogEventBarrierBegin(VEC_MDotBarrier,x,*y,0,0,x->comm);
1480: (*x->ops->mdot)(nv,x,y,val);
1481: PetscLogEventBarrierEnd(VEC_MDotBarrier,x,*y,0,0,x->comm);
1482: return(0);
1483: }
1485: /*@C
1486: VecMAXPY - Computes y = y + sum alpha[j] x[j]
1488: Collective on Vec
1490: Input Parameters:
1491: + nv - number of scalars and x-vectors
1492: . alpha - array of scalars
1493: . y - one vector
1494: - x - array of vectors
1496: Level: intermediate
1498: Concepts: BLAS
1500: .seealso: VecAXPY(), VecWAXPY(), VecAYPX()
1501: @*/
1502: int VecMAXPY(int nv,const Scalar *alpha,Vec y,Vec *x)
1503: {
1514: if (y->N != (*x)->N) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths");
1515: if (y->n != (*x)->n) SETERRQ(PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths");
1517: PetscLogEventBegin(VEC_MAXPY,*x,y,0,0);
1518: (*y->ops->maxpy)(nv,alpha,y,x);
1519: PetscLogEventEnd(VEC_MAXPY,*x,y,0,0);
1520: return(0);
1521: }
1523: /*@C
1524: VecGetArray - Returns a pointer to a contiguous array that contains this
1525: processor's portion of the vector data. For the standard PETSc
1526: vectors, VecGetArray() returns a pointer to the local data array and
1527: does not use any copies. If the underlying vector data is not stored
1528: in a contiquous array this routine will copy the data to a contiquous
1529: array and return a pointer to that. You MUST call VecRestoreArray()
1530: when you no longer need access to the array.
1532: Not Collective
1534: Input Parameter:
1535: . x - the vector
1537: Output Parameter:
1538: . a - location to put pointer to the array
1540: Fortran Note:
1541: This routine is used differently from Fortran 77
1542: $ Vec x
1543: $ Scalar x_array(1)
1544: $ PetscOffset i_x
1545: $ int ierr
1546: $ call VecGetArray(x,x_array,i_x,ierr)
1547: $
1548: $ Access first local entry in vector with
1549: $ value = x_array(i_x + 1)
1550: $
1551: $ ...... other code
1552: $ call VecRestoreArray(x,x_array,i_x,ierr)
1553: For Fortran 90 see VecGetArrayF90()
1555: See the Fortran chapter of the users manual and
1556: petsc/src/snes/examples/tutorials/ex5f.F for details.
1558: Level: beginner
1560: Concepts: vector^accessing local values
1562: .seealso: VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(), VecGetArray2d()
1563: @*/
1564: int VecGetArray(Vec x,Scalar *a[])
1565: {
1572: (*x->ops->getarray)(x,a);
1573: return(0);
1574: }
1577: /*@C
1578: VecGetArrays - Returns a pointer to the arrays in a set of vectors
1579: that were created by a call to VecDuplicateVecs(). You MUST call
1580: VecRestoreArrays() when you no longer need access to the array.
1582: Not Collective
1584: Input Parameter:
1585: + x - the vectors
1586: - n - the number of vectors
1588: Output Parameter:
1589: . a - location to put pointer to the array
1591: Fortran Note:
1592: This routine is not supported in Fortran.
1594: Level: intermediate
1596: .seealso: VecGetArray(), VecRestoreArrays()
1597: @*/
1598: int VecGetArrays(const Vec x[],int n,Scalar **a[])
1599: {
1600: int i,ierr;
1601: Scalar **q;
1606: if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must get at least one array n = %d",n);
1607: PetscMalloc(n*sizeof(Scalar*),&q);
1608: for (i=0; i<n; ++i) {
1609: VecGetArray(x[i],&q[i]);
1610: }
1611: *a = q;
1612: return(0);
1613: }
1615: /*@C
1616: VecRestoreArrays - Restores a group of vectors after VecGetArrays()
1617: has been called.
1619: Not Collective
1621: Input Parameters:
1622: + x - the vector
1623: . n - the number of vectors
1624: - a - location of pointer to arrays obtained from VecGetArrays()
1626: Notes:
1627: For regular PETSc vectors this routine does not involve any copies. For
1628: any special vectors that do not store local vector data in a contiguous
1629: array, this routine will copy the data back into the underlying
1630: vector data structure from the arrays obtained with VecGetArrays().
1632: Fortran Note:
1633: This routine is not supported in Fortran.
1635: Level: intermediate
1637: .seealso: VecGetArrays(), VecRestoreArray()
1638: @*/
1639: int VecRestoreArrays(const Vec x[],int n,Scalar **a[])
1640: {
1641: int i,ierr;
1642: Scalar **q = *a;
1648: for(i=0;i<n;++i) {
1649: VecRestoreArray(x[i],&q[i]);
1650: }
1651: PetscFree(q);
1652: return(0);
1653: }
1655: /*@C
1656: VecRestoreArray - Restores a vector after VecGetArray() has been called.
1658: Not Collective
1660: Input Parameters:
1661: + x - the vector
1662: - a - location of pointer to array obtained from VecGetArray()
1664: Level: beginner
1666: Notes:
1667: For regular PETSc vectors this routine does not involve any copies. For
1668: any special vectors that do not store local vector data in a contiguous
1669: array, this routine will copy the data back into the underlying
1670: vector data structure from the array obtained with VecGetArray().
1672: This routine actually zeros out the a pointer. This is to prevent accidental
1673: us of the array after it has been restored. If you pass null for a it will
1674: not zero the array pointer a.
1676: Fortran Note:
1677: This routine is used differently from Fortran 77
1678: $ Vec x
1679: $ Scalar x_array(1)
1680: $ PetscOffset i_x
1681: $ int ierr
1682: $ call VecGetArray(x,x_array,i_x,ierr)
1683: $
1684: $ Access first local entry in vector with
1685: $ value = x_array(i_x + 1)
1686: $
1687: $ ...... other code
1688: $ call VecRestoreArray(x,x_array,i_x,ierr)
1690: See the Fortran chapter of the users manual and
1691: petsc/src/snes/examples/tutorials/ex5f.F for details.
1692: For Fortran 90 see VecRestoreArrayF90()
1694: .seealso: VecGetArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(), VecRestoreArray2d()
1695: @*/
1696: int VecRestoreArray(Vec x,Scalar *a[])
1697: {
1704: if (x->ops->restorearray) {
1705: (*x->ops->restorearray)(x,a);
1706: }
1707: return(0);
1708: }
1711: /*@C
1712: VecView - Views a vector object.
1714: Collective on Vec
1716: Input Parameters:
1717: + v - the vector
1718: - viewer - an optional visualization context
1720: Notes:
1721: The available visualization contexts include
1722: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
1723: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
1724: output where only the first processor opens
1725: the file. All other processors send their
1726: data to the first processor to print.
1728: You can change the format the vector is printed using the
1729: option PetscViewerSetFormat().
1731: The user can open alternative visualization contexts with
1732: + PetscViewerASCIIOpen() - Outputs vector to a specified file
1733: . PetscViewerBinaryOpen() - Outputs vector in binary to a
1734: specified file; corresponding input uses VecLoad()
1735: . PetscViewerDrawOpen() - Outputs vector to an X window display
1736: - PetscViewerSocketOpen() - Outputs vector to Socket viewer
1738: The user can call PetscViewerSetFormat() to specify the output
1739: format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
1740: PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include
1741: + PETSC_VIEWER_ASCII_DEFAULT - default, prints vector contents
1742: . PETSC_VIEWER_ASCII_MATLAB - prints vector contents in Matlab format
1743: . PETSC_VIEWER_ASCII_INDEX - prints vector contents, including indices of vector elements
1744: - PETSC_VIEWER_ASCII_COMMON - prints vector contents, using a
1745: format common among all vector types
1747: Level: beginner
1749: Concepts: vector^printing
1750: Concepts: vector^saving to disk
1752: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscDrawLGCreate(),
1753: PetscViewerSocketOpen(), PetscViewerBinaryOpen(), VecLoad(), PetscViewerCreate(),
1754: PetscDoubleView(), PetscScalarView(), PetscIntView()
1755: @*/
1756: int VecView(Vec vec,PetscViewer viewer)
1757: {
1758: int ierr;
1759: PetscViewerFormat format;
1764: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(vec->comm);
1767: if (vec->stash.n || vec->bstash.n) SETERRQ(1,"Must call VecAssemblyBegin/End() before viewing this vector");
1769: /*
1770: Check if default viewer has been overridden, but user request it anyways
1771: */
1772: PetscViewerGetFormat(viewer,&format);
1773: if (vec->ops->viewnative && format == PETSC_VIEWER_NATIVE) {
1774: ierr = PetscViewerPopFormat(viewer);
1775: (*vec->ops->viewnative)(vec,viewer);
1776: ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE);
1777: } else {
1778: (*vec->ops->view)(vec,viewer);
1779: }
1780: return(0);
1781: }
1783: /*@
1784: VecGetSize - Returns the global number of elements of the vector.
1786: Not Collective
1788: Input Parameter:
1789: . x - the vector
1791: Output Parameters:
1792: . size - the global length of the vector
1794: Level: beginner
1796: Concepts: vector^local size
1798: .seealso: VecGetLocalSize()
1799: @*/
1800: int VecGetSize(Vec x,int *size)
1801: {
1808: (*x->ops->getsize)(x,size);
1809: return(0);
1810: }
1812: /*@
1813: VecGetLocalSize - Returns the number of elements of the vector stored
1814: in local memory. This routine may be implementation dependent, so use
1815: with care.
1817: Not Collective
1819: Input Parameter:
1820: . x - the vector
1822: Output Parameter:
1823: . size - the length of the local piece of the vector
1825: Level: beginner
1827: Concepts: vector^size
1829: .seealso: VecGetSize()
1830: @*/
1831: int VecGetLocalSize(Vec x,int *size)
1832: {
1839: (*x->ops->getlocalsize)(x,size);
1840: return(0);
1841: }
1843: /*@
1844: VecGetOwnershipRange - Returns the range of indices owned by
1845: this processor, assuming that the vectors are laid out with the
1846: first n1 elements on the first processor, next n2 elements on the
1847: second, etc. For certain parallel layouts this range may not be
1848: well defined.
1850: Not Collective
1852: Input Parameter:
1853: . x - the vector
1855: Output Parameters:
1856: + low - the first local element, pass in PETSC_NULL if not interested
1857: - high - one more than the last local element, pass in PETSC_NULL if not interested
1859: Note:
1860: The high argument is one more than the last element stored locally.
1862: Level: beginner
1864: Concepts: ownership^of vectors
1865: Concepts: vector^ownership of elements
1867: @*/
1868: int VecGetOwnershipRange(Vec x,int *low,int *high)
1869: {
1871: Map map;
1878: (*x->ops->getmap)(x,&map);
1879: MapGetLocalRange(map,low,high);
1880: return(0);
1881: }
1883: /*@C
1884: VecGetMap - Returns the map associated with the vector
1886: Not Collective
1888: Input Parameter:
1889: . x - the vector
1891: Output Parameters:
1892: . map - the map
1894: Level: developer
1896: @*/
1897: int VecGetMap(Vec x,Map *map)
1898: {
1904: (*x->ops->getmap)(x,map);
1905: return(0);
1906: }
1908: /*@
1909: VecSetOption - Sets an option for controling a vector's behavior.
1911: Collective on Vec
1913: Input Parameter:
1914: + x - the vector
1915: - op - the option
1917: Notes:
1918: Currently the only option supported is
1919: VEC_IGNORE_OFF_PROC_ENTRIES, which causes VecSetValues() to ignore
1920: entries destined to be stored on a seperate processor.
1922: Level: intermediate
1924: @*/
1925: int VecSetOption(Vec x,VecOption op)
1926: {
1932: if (x->ops->setoption) {
1933: (*x->ops->setoption)(x,op);
1934: }
1935: return(0);
1936: }
1938: /* Default routines for obtaining and releasing; */
1939: /* may be used by any implementation */
1940: int VecDuplicateVecs_Default(Vec w,int m,Vec *V[])
1941: {
1942: int i,ierr;
1947: if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
1948: PetscMalloc(m*sizeof(Vec*),V);
1949: for (i=0; i<m; i++) {VecDuplicate(w,*V+i);}
1950: return(0);
1951: }
1953: int VecDestroyVecs_Default(const Vec v[], int m)
1954: {
1955: int i,ierr;
1959: if (m <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %d",m);
1960: for (i=0; i<m; i++) {VecDestroy(v[i]);}
1961: PetscFree((Vec*)v);
1962: return(0);
1963: }
1965: /*@
1966: VecPlaceArray - Allows one to replace the array in a vector with an
1967: array provided by the user. This is useful to avoid copying an array
1968: into a vector.
1970: Not Collective
1972: Input Parameters:
1973: + vec - the vector
1974: - array - the array
1976: Notes:
1977: You can return to the original array with a call to VecResetArray()
1979: Level: developer
1981: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecResetArray()
1983: @*/
1984: int VecPlaceArray(Vec vec,const Scalar array[])
1985: {
1991: if (vec->ops->placearray) {
1992: (*vec->ops->placearray)(vec,array);
1993: } else {
1994: SETERRQ(1,"Cannot place array in this type of vector");
1995: }
1996: return(0);
1997: }
1999: /*@
2000: VecResetArray - Resets a vector to use its default memory. Call this
2001: after the use of VecPlaceArray().
2003: Not Collective
2005: Input Parameters:
2006: . vec - the vector
2008: Level: developer
2010: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecPlaceArray()
2012: @*/
2013: int VecResetArray(Vec vec)
2014: {
2020: if (vec->ops->resetarray) {
2021: (*vec->ops->resetarray)(vec);
2022: } else {
2023: SETERRQ(1,"Cannot reset array in this type of vector");
2024: }
2025: return(0);
2026: }
2028: /*@C
2029: VecReplaceArray - Allows one to replace the array in a vector with an
2030: array provided by the user. This is useful to avoid copying an array
2031: into a vector.
2033: Not Collective
2035: Input Parameters:
2036: + vec - the vector
2037: - array - the array
2039: Notes:
2040: This permanently replaces the array and frees the memory associated
2041: with the old array.
2043: Not supported from Fortran
2045: Level: developer
2047: .seealso: VecGetArray(), VecRestoreArray(), VecPlaceArray()
2049: @*/
2050: int VecReplaceArray(Vec vec,const Scalar array[])
2051: {
2057: if (vec->ops->replacearray) {
2058: (*vec->ops->replacearray)(vec,array);
2059: } else {
2060: SETERRQ(1,"Cannot replace array in this type of vector");
2061: }
2062: return(0);
2063: }
2065: /*MC
2066: VecDuplicateVecsF90 - Creates several vectors of the same type as an existing vector
2067: and makes them accessible via a Fortran90 pointer.
2069: Synopsis:
2070: VecGetArrayF90(Vec x,int n,{Scalar, pointer :: y(:)},integer ierr)
2072: Collective on Vec
2074: Input Parameters:
2075: + x - a vector to mimic
2076: - n - the number of vectors to obtain
2078: Output Parameters:
2079: + y - Fortran90 pointer to the array of vectors
2080: - ierr - error code
2082: Example of Usage:
2083: .vb
2084: Vec x
2085: Vec, pointer :: y(:)
2086: ....
2087: call VecDuplicateVecsF90(x,2,y,ierr)
2088: call VecSet(alpha,y(2),ierr)
2089: call VecSet(alpha,y(2),ierr)
2090: ....
2091: call VecDestroyVecsF90(y,2,ierr)
2092: .ve
2094: Notes:
2095: Not yet supported for all F90 compilers
2097: Use VecDestroyVecsF90() to free the space.
2099: Level: beginner
2101: .seealso: VecDestroyVecsF90(), VecDuplicateVecs()
2103: M*/
2105: /*MC
2106: VecRestoreArrayF90 - Restores a vector to a usable state after a call to
2107: VecGetArrayF90().
2109: Synopsis:
2110: VecRestoreArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
2112: Not collective
2114: Input Parameters:
2115: + x - vector
2116: - xx_v - the Fortran90 pointer to the array
2118: Output Parameter:
2119: . ierr - error code
2121: Example of Usage:
2122: .vb
2123: Scalar, pointer :: xx_v(:)
2124: ....
2125: call VecGetArrayF90(x,xx_v,ierr)
2126: a = xx_v(3)
2127: call VecRestoreArrayF90(x,xx_v,ierr)
2128: .ve
2129:
2130: Notes:
2131: Not yet supported for all F90 compilers
2133: Level: beginner
2135: .seealso: VecGetArrayF90(), VecGetArray(), VecRestoreArray()
2137: M*/
2139: /*MC
2140: VecDestroyVecsF90 - Frees a block of vectors obtained with VecDuplicateVecsF90().
2142: Synopsis:
2143: VecDestroyVecsF90({Scalar, pointer :: x(:)},integer n,integer ierr)
2145: Input Parameters:
2146: + x - pointer to array of vector pointers
2147: - n - the number of vectors previously obtained
2149: Output Parameter:
2150: . ierr - error code
2152: Notes:
2153: Not yet supported for all F90 compilers
2155: Level: beginner
2157: .seealso: VecDestroyVecs(), VecDuplicateVecsF90()
2159: M*/
2161: /*MC
2162: VecGetArrayF90 - Accesses a vector array from Fortran90. For default PETSc
2163: vectors, VecGetArrayF90() returns a pointer to the local data array. Otherwise,
2164: this routine is implementation dependent. You MUST call VecRestoreArrayF90()
2165: when you no longer need access to the array.
2167: Synopsis:
2168: VecGetArrayF90(Vec x,{Scalar, pointer :: xx_v(:)},integer ierr)
2170: Not Collective
2172: Input Parameter:
2173: . x - vector
2175: Output Parameters:
2176: + xx_v - the Fortran90 pointer to the array
2177: - ierr - error code
2179: Example of Usage:
2180: .vb
2181: Scalar, pointer :: xx_v(:)
2182: ....
2183: call VecGetArrayF90(x,xx_v,ierr)
2184: a = xx_v(3)
2185: call VecRestoreArrayF90(x,xx_v,ierr)
2186: .ve
2188: Notes:
2189: Not yet supported for all F90 compilers
2191: Level: beginner
2193: .seealso: VecRestoreArrayF90(), VecGetArray(), VecRestoreArray()
2195: M*/
2197: /*@C
2198: VecLoadIntoVector - Loads a vector that has been stored in binary format
2199: with VecView().
2201: Collective on PetscViewer
2203: Input Parameters:
2204: + viewer - binary file viewer, obtained from PetscViewerBinaryOpen()
2205: - vec - vector to contain files values (must be of correct length)
2207: Level: intermediate
2209: Notes:
2210: The input file must contain the full global vector, as
2211: written by the routine VecView().
2213: Use VecLoad() to create the vector as the values are read in
2215: Notes for advanced users:
2216: Most users should not need to know the details of the binary storage
2217: format, since VecLoad() and VecView() completely hide these details.
2218: But for anyone who's interested, the standard binary matrix storage
2219: format is
2220: .vb
2221: int VEC_COOKIE
2222: int number of rows
2223: Scalar *values of all nonzeros
2224: .ve
2226: Note for Cray users, the int's stored in the binary file are 32 bit
2227: integers; not 64 as they are represented in the memory, so if you
2228: write your own routines to read/write these binary files from the Cray
2229: you need to adjust the integer sizes that you read in, see
2230: PetscReadBinary() and PetscWriteBinary() to see how this may be
2231: done.
2233: In addition, PETSc automatically does the byte swapping for
2234: machines that store the bytes reversed, e.g. DEC alpha, freebsd,
2235: linux, nt and the paragon; thus if you write your own binary
2236: read/write routines you have to swap the bytes; see PetscReadBinary()
2237: and PetscWriteBinary() to see how this may be done.
2239: Concepts: vector^loading from file
2241: .seealso: PetscViewerBinaryOpen(), VecView(), MatLoad(), VecLoad()
2242: @*/
2243: int VecLoadIntoVector(PetscViewer viewer,Vec vec)
2244: {
2251: if (!vec->ops->loadintovector) {
2252: SETERRQ(1,"Vector does not support load");
2253: }
2254: (*vec->ops->loadintovector)(viewer,vec);
2255: return(0);
2256: }
2258: /*@
2259: VecReciprocal - Replaces each component of a vector by its reciprocal.
2261: Collective on Vec
2263: Input Parameter:
2264: . v - the vector
2266: Output Parameter:
2267: . v - the vector reciprocal
2269: Level: intermediate
2271: Concepts: vector^reciprocal
2273: @*/
2274: int VecReciprocal(Vec vec)
2275: {
2276: int ierr;
2281: if (!vec->ops->reciprocal) {
2282: SETERRQ(1,"Vector does not support reciprocal operation");
2283: }
2284: (*vec->ops->reciprocal)(vec);
2285: return(0);
2286: }
2288: int VecSetOperation(Vec vec,VecOperation op, void (*f)())
2289: {
2292: /* save the native version of the viewer */
2293: if (op == VECOP_VIEW && !vec->ops->viewnative) {
2294: vec->ops->viewnative = vec->ops->view;
2295: }
2296: (((void(**)())vec->ops)[(int)op]) = f;
2297: return(0);
2298: }
2300: /*@
2301: VecSetStashInitialSize - sets the sizes of the vec-stash, that is
2302: used during the assembly process to store values that belong to
2303: other processors.
2305: Collective on Vec
2307: Input Parameters:
2308: + vec - the vector
2309: . size - the initial size of the stash.
2310: - bsize - the initial size of the block-stash(if used).
2312: Options Database Keys:
2313: + -vecstash_initial_size <size> or <size0,size1,...sizep-1>
2314: - -vecstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1>
2316: Level: intermediate
2318: Notes:
2319: The block-stash is used for values set with VecSetValuesBlocked() while
2320: the stash is used for values set with VecSetValues()
2322: Run with the option -log_info and look for output of the form
2323: VecAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
2324: to determine the appropriate value, MM, to use for size and
2325: VecAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
2326: to determine the value, BMM to use for bsize
2328: Concepts: vector^stash
2329: Concepts: stash^vector
2331: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked(), VecStashView()
2333: @*/
2334: int VecSetStashInitialSize(Vec vec,int size,int bsize)
2335: {
2340: VecStashSetInitialSize_Private(&vec->stash,size);
2341: VecStashSetInitialSize_Private(&vec->bstash,bsize);
2342: return(0);
2343: }
2345: /*@
2346: VecStashView - Prints the entries in the vector stash and block stash.
2348: Collective on Vec
2350: Input Parameters:
2351: + vec - the vector
2352: - viewer - the viewer
2354: Level: advanced
2356: Concepts: vector^stash
2357: Concepts: stash^vector
2359: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked()
2361: @*/
2362: int VecStashView(Vec v,PetscViewer viewer)
2363: {
2364: int ierr,rank,i,j;
2365: PetscTruth match;
2366: VecStash *s;
2367: Scalar val;
2374: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&match);
2375: if (!match) SETERRQ1(1,"Stash viewer only works with ASCII viewer not %sn",((PetscObject)v)->type_name);
2376: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
2377: MPI_Comm_rank(v->comm,&rank);
2378: s = &v->bstash;
2380: /* print block stash */
2381: PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector Block stash size %d block size %dn",rank,s->n,s->bs);
2382: for (i=0; i<s->n; i++) {
2383: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d ",rank,s->idx[i]);
2384: for (j=0; j<s->bs; j++) {
2385: val = s->array[i*s->bs+j];
2386: #if defined(PETSC_USE_COMPLEX)
2387: PetscViewerASCIISynchronizedPrintf(viewer,"(%18.16e %18.16e) ",PetscRealPart(val),PetscImaginaryPart(val));
2388: #else
2389: PetscViewerASCIISynchronizedPrintf(viewer,"%18.16e ",val);
2390: #endif
2391: }
2392: PetscViewerASCIISynchronizedPrintf(viewer,"n");
2393: }
2394: PetscViewerFlush(viewer);
2396: s = &v->stash;
2398: /* print basic stash */
2399: PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector stash size %dn",rank,s->n);
2400: for (i=0; i<s->n; i++) {
2401: val = s->array[i];
2402: #if defined(PETSC_USE_COMPLEX)
2403: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element % (%18.16e %18.16e) ",rank,s->idx[i],PetscRealPart(val),PetscImaginaryPart(val));
2404: #else
2405: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %d %18.16en",rank,s->idx[i],val);
2406: #endif
2407: }
2408: PetscViewerFlush(viewer);
2410: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
2411: return(0);
2412: }
2414: /*@C
2415: VecGetArray2d - Returns a pointer to a 2d contiguous array that contains this
2416: processor's portion of the vector data. You MUST call VecRestoreArray2d()
2417: when you no longer need access to the array.
2419: Not Collective
2421: Input Parameter:
2422: + x - the vector
2423: . m - first dimension of two dimensional array
2424: . n - second dimension of two dimensional array
2425: . mstart - first index you will use in first coordinate direction (often 0)
2426: - nstart - first index in the second coordinate direction (often 0)
2428: Output Parameter:
2429: . a - location to put pointer to the array
2431: Level: beginner
2433: Notes:
2434: For a vector obtained from DACreateLocalVector() mstart and nstart are likely
2435: obtained from the corner indices obtained from DAGetGhostCorners() while for
2436: DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2437: the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray2d().
2438:
2439: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2441: Concepts: vector^accessing local values as 2d array
2443: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2444: VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2445: VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2446: @*/
2447: int VecGetArray2d(Vec x,int m,int n,int mstart,int nstart,Scalar **a[])
2448: {
2449: int i,ierr,N;
2450: Scalar *aa;
2456: VecGetLocalSize(x,&N);
2457: if (m*n != N) SETERRQ3(1,"Local array size %d does not match 2d array dimensions %d by %d",N,m,n);
2458: VecGetArray(x,&aa);
2460: PetscMalloc(m*sizeof(Scalar*),a);
2461: for (i=0; i<m; i++) (*a)[i] = aa + i*n - nstart;
2462: *a -= mstart;
2463: return(0);
2464: }
2466: /*@C
2467: VecRestoreArray2d - Restores a vector after VecGetArray2d() has been called.
2469: Not Collective
2471: Input Parameters:
2472: + x - the vector
2473: . m - first dimension of two dimensional array
2474: . n - second dimension of the two dimensional array
2475: . mstart - first index you will use in first coordinate direction (often 0)
2476: . nstart - first index in the second coordinate direction (often 0)
2477: - a - location of pointer to array obtained from VecGetArray2d()
2479: Level: beginner
2481: Notes:
2482: For regular PETSc vectors this routine does not involve any copies. For
2483: any special vectors that do not store local vector data in a contiguous
2484: array, this routine will copy the data back into the underlying
2485: vector data structure from the array obtained with VecGetArray().
2487: This routine actually zeros out the a pointer.
2489: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2490: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2491: VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2492: @*/
2493: int VecRestoreArray2d(Vec x,int m,int n,int mstart,int nstart,Scalar **a[])
2494: {
2500: PetscFree(*a + mstart);
2501: VecRestoreArray(x,PETSC_NULL);
2502: return(0);
2503: }
2505: /*@C
2506: VecGetArray1d - Returns a pointer to a 1d contiguous array that contains this
2507: processor's portion of the vector data. You MUST call VecRestoreArray1d()
2508: when you no longer need access to the array.
2510: Not Collective
2512: Input Parameter:
2513: + x - the vector
2514: . m - first dimension of two dimensional array
2515: - mstart - first index you will use in first coordinate direction (often 0)
2517: Output Parameter:
2518: . a - location to put pointer to the array
2520: Level: beginner
2522: Notes:
2523: For a vector obtained from DACreateLocalVector() mstart are likely
2524: obtained from the corner indices obtained from DAGetGhostCorners() while for
2525: DACreateGlobalVector() they are the corner indices from DAGetCorners().
2526:
2527: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2529: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2530: VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2531: VecGetarray2d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2532: @*/
2533: int VecGetArray1d(Vec x,int m,int mstart,Scalar *a[])
2534: {
2535: int ierr,N;
2541: VecGetLocalSize(x,&N);
2542: if (m != N) SETERRQ2(1,"Local array size %d does not match 1d array dimensions %d",N,m);
2543: VecGetArray(x,a);
2544: *a -= mstart;
2545: return(0);
2546: }
2548: /*@C
2549: VecRestoreArray1d - Restores a vector after VecGetArray1d() has been called.
2551: Not Collective
2553: Input Parameters:
2554: + x - the vector
2555: . m - first dimension of two dimensional array
2556: . mstart - first index you will use in first coordinate direction (often 0)
2557: - a - location of pointer to array obtained from VecGetArray21()
2559: Level: beginner
2561: Notes:
2562: For regular PETSc vectors this routine does not involve any copies. For
2563: any special vectors that do not store local vector data in a contiguous
2564: array, this routine will copy the data back into the underlying
2565: vector data structure from the array obtained with VecGetArray1d().
2567: This routine actually zeros out the a pointer.
2569: Concepts: vector^accessing local values as 1d array
2571: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2572: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2573: VecGetarray1d(), VecRestoreArray2d(), VecGetArray4d(), VecRestoreArray4d()
2574: @*/
2575: int VecRestoreArray1d(Vec x,int m,int mstart,Scalar *a[])
2576: {
2582: VecRestoreArray(x,PETSC_NULL);
2583: return(0);
2584: }
2586: /*@C
2587: VecConjugate - Conjugates a vector.
2589: Collective on Vec
2591: Input Parameters:
2592: . x - the vector
2594: Level: intermediate
2596: Concepts: vector^conjugate
2598: @*/
2599: int VecConjugate(Vec x)
2600: {
2601: #ifdef PETSC_USE_COMPLEX
2607: (*x->ops->conjugate)(x);
2608: return(0);
2609: #else
2610: return(0);
2611: #endif
2612: }
2614: /*@C
2615: VecGetArray3d - Returns a pointer to a 3d contiguous array that contains this
2616: processor's portion of the vector data. You MUST call VecRestoreArray3d()
2617: when you no longer need access to the array.
2619: Not Collective
2621: Input Parameter:
2622: + x - the vector
2623: . m - first dimension of three dimensional array
2624: . n - second dimension of three dimensional array
2625: . p - third dimension of three dimensional array
2626: . mstart - first index you will use in first coordinate direction (often 0)
2627: . nstart - first index in the second coordinate direction (often 0)
2628: - pstart - first index in the third coordinate direction (often 0)
2630: Output Parameter:
2631: . a - location to put pointer to the array
2633: Level: beginner
2635: Notes:
2636: For a vector obtained from DACreateLocalVector() mstart, nstart, and pstart are likely
2637: obtained from the corner indices obtained from DAGetGhostCorners() while for
2638: DACreateGlobalVector() they are the corner indices from DAGetCorners(). In both cases
2639: the arguments from DAGet[Ghost}Corners() are reversed in the call to VecGetArray3d().
2640:
2641: For standard PETSc vectors this is an inexpensive call; it does not copy the vector values.
2643: Concepts: vector^accessing local values as 3d array
2645: .seealso: VecGetArray(), VecRestoreArray(), VecGetArrays(), VecGetArrayF90(), VecPlaceArray(),
2646: VecRestoreArray2d(), DAVecGetarray(), DAVecRestoreArray(), VecGetArray3d(), VecRestoreArray3d(),
2647: VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d()
2648: @*/
2649: int VecGetArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,Scalar ***a[])
2650: {
2651: int i,ierr,N,j;
2652: Scalar *aa,**b;
2658: VecGetLocalSize(x,&N);
2659: if (m*n*p != N) SETERRQ4(1,"Local array size %d does not match 3d array dimensions %d by %d by %d",N,m,n,p);
2660: VecGetArray(x,&aa);
2662: PetscMalloc(m*sizeof(Scalar**)+m*n*sizeof(Scalar*),a);
2663: b = (Scalar **)((*a) + m);
2664: for (i=0; i<m; i++) (*a)[i] = b + i*n - nstart;
2665: for (i=0; i<m; i++) {
2666: for (j=0; j<n; j++) {
2667: b[i*n+j] = aa + i*n*p + j*p - pstart;
2668: CHKMEMQ;
2669: }
2670: }
2671: *a -= mstart;
2672: CHKMEMQ;
2673: return(0);
2674: }
2676: /*@C
2677: VecRestoreArray3d - Restores a vector after VecGetArray3d() has been called.
2679: Not Collective
2681: Input Parameters:
2682: + x - the vector
2683: . m - first dimension of three dimensional array
2684: . n - second dimension of the three dimensional array
2685: . p - third dimension of the three dimensional array
2686: . mstart - first index you will use in first coordinate direction (often 0)
2687: . nstart - first index in the second coordinate direction (often 0)
2688: . pstart - first index in the third coordinate direction (often 0)
2689: - a - location of pointer to array obtained from VecGetArray3d()
2691: Level: beginner
2693: Notes:
2694: For regular PETSc vectors this routine does not involve any copies. For
2695: any special vectors that do not store local vector data in a contiguous
2696: array, this routine will copy the data back into the underlying
2697: vector data structure from the array obtained with VecGetArray().
2699: This routine actually zeros out the a pointer.
2701: .seealso: VecGetArray(), VecRestoreArray(), VecRestoreArrays(), VecRestoreArrayF90(), VecPlaceArray(),
2702: VecGetArray2d(), VecGetArray3d(), VecRestoreArray3d(), DAVecGetArray(), DAVecRestoreArray()
2703: VecGetarray1d(), VecRestoreArray1d(), VecGetArray4d(), VecRestoreArray4d(), VecGet
2704: @*/
2705: int VecRestoreArray3d(Vec x,int m,int n,int p,int mstart,int nstart,int pstart,Scalar ***a[])
2706: {
2712: PetscFree(*a + mstart);
2713: VecRestoreArray(x,PETSC_NULL);
2714: return(0);
2715: }