Actual source code: vector.c
petsc-3.13.6 2020-09-29
1: /*
2: Provides the interface functions for vector operations that do NOT have PetscScalar/PetscReal in the signature
3: These are the vector functions the user calls.
4: */
5: #include <petsc/private/vecimpl.h>
7: /* Logging support */
8: PetscClassId VEC_CLASSID;
9: PetscLogEvent VEC_View, VEC_Max, VEC_Min, VEC_Dot, VEC_MDot, VEC_TDot;
10: PetscLogEvent VEC_Norm, VEC_Normalize, VEC_Scale, VEC_Copy, VEC_Set, VEC_AXPY, VEC_AYPX, VEC_WAXPY;
11: PetscLogEvent VEC_MTDot, VEC_MAXPY, VEC_Swap, VEC_AssemblyBegin, VEC_ScatterBegin, VEC_ScatterEnd;
12: PetscLogEvent VEC_AssemblyEnd, VEC_PointwiseMult, VEC_SetValues, VEC_Load;
13: PetscLogEvent VEC_SetRandom, VEC_ReduceArithmetic, VEC_ReduceCommunication,VEC_ReduceBegin,VEC_ReduceEnd,VEC_Ops;
14: PetscLogEvent VEC_DotNorm2, VEC_AXPBYPCZ;
15: PetscLogEvent VEC_ViennaCLCopyFromGPU, VEC_ViennaCLCopyToGPU;
16: PetscLogEvent VEC_CUDACopyFromGPU, VEC_CUDACopyToGPU;
17: PetscLogEvent VEC_CUDACopyFromGPUSome, VEC_CUDACopyToGPUSome;
19: /*@
20: VecStashGetInfo - Gets how many values are currently in the vector stash, i.e. need
21: to be communicated to other processors during the VecAssemblyBegin/End() process
23: Not collective
25: Input Parameter:
26: . vec - the vector
28: Output Parameters:
29: + nstash - the size of the stash
30: . reallocs - the number of additional mallocs incurred.
31: . bnstash - the size of the block stash
32: - breallocs - the number of additional mallocs incurred.in the block stash
34: Level: advanced
36: .seealso: VecAssemblyBegin(), VecAssemblyEnd(), Vec, VecStashSetInitialSize(), VecStashView()
38: @*/
39: PetscErrorCode VecStashGetInfo(Vec vec,PetscInt *nstash,PetscInt *reallocs,PetscInt *bnstash,PetscInt *breallocs)
40: {
44: VecStashGetInfo_Private(&vec->stash,nstash,reallocs);
45: VecStashGetInfo_Private(&vec->bstash,bnstash,breallocs);
46: return(0);
47: }
49: /*@
50: VecSetLocalToGlobalMapping - Sets a local numbering to global numbering used
51: by the routine VecSetValuesLocal() to allow users to insert vector entries
52: using a local (per-processor) numbering.
54: Logically Collective on Vec
56: Input Parameters:
57: + x - vector
58: - mapping - mapping created with ISLocalToGlobalMappingCreate() or ISLocalToGlobalMappingCreateIS()
60: Notes:
61: All vectors obtained with VecDuplicate() from this vector inherit the same mapping.
63: Level: intermediate
65: seealso: VecAssemblyBegin(), VecAssemblyEnd(), VecSetValues(), VecSetValuesLocal(),
66: VecSetLocalToGlobalMapping(), VecSetValuesBlockedLocal()
67: @*/
68: PetscErrorCode VecSetLocalToGlobalMapping(Vec x,ISLocalToGlobalMapping mapping)
69: {
76: if (x->ops->setlocaltoglobalmapping) {
77: (*x->ops->setlocaltoglobalmapping)(x,mapping);
78: } else {
79: PetscLayoutSetISLocalToGlobalMapping(x->map,mapping);
80: }
81: return(0);
82: }
84: /*@
85: VecGetLocalToGlobalMapping - Gets the local-to-global numbering set by VecSetLocalToGlobalMapping()
87: Not Collective
89: Input Parameter:
90: . X - the vector
92: Output Parameter:
93: . mapping - the mapping
95: Level: advanced
98: .seealso: VecSetValuesLocal()
99: @*/
100: PetscErrorCode VecGetLocalToGlobalMapping(Vec X,ISLocalToGlobalMapping *mapping)
101: {
106: *mapping = X->map->mapping;
107: return(0);
108: }
110: /*@
111: VecAssemblyBegin - Begins assembling the vector. This routine should
112: be called after completing all calls to VecSetValues().
114: Collective on Vec
116: Input Parameter:
117: . vec - the vector
119: Level: beginner
121: .seealso: VecAssemblyEnd(), VecSetValues()
122: @*/
123: PetscErrorCode VecAssemblyBegin(Vec vec)
124: {
130: VecStashViewFromOptions(vec,NULL,"-vec_view_stash");
131: PetscLogEventBegin(VEC_AssemblyBegin,vec,0,0,0);
132: if (vec->ops->assemblybegin) {
133: (*vec->ops->assemblybegin)(vec);
134: }
135: PetscLogEventEnd(VEC_AssemblyBegin,vec,0,0,0);
136: PetscObjectStateIncrease((PetscObject)vec);
137: return(0);
138: }
140: /*@
141: VecAssemblyEnd - Completes assembling the vector. This routine should
142: be called after VecAssemblyBegin().
144: Collective on Vec
146: Input Parameter:
147: . vec - the vector
149: Options Database Keys:
150: + -vec_view - Prints vector in ASCII format
151: . -vec_view ::ascii_matlab - Prints vector in ASCII MATLAB format to stdout
152: . -vec_view matlab:filename - Prints vector in MATLAB format to matlaboutput.mat
153: . -vec_view draw - Activates vector viewing using drawing tools
154: . -display <name> - Sets display name (default is host)
155: . -draw_pause <sec> - Sets number of seconds to pause after display
156: - -vec_view socket - Activates vector viewing using a socket
158: Level: beginner
160: .seealso: VecAssemblyBegin(), VecSetValues()
161: @*/
162: PetscErrorCode VecAssemblyEnd(Vec vec)
163: {
168: PetscLogEventBegin(VEC_AssemblyEnd,vec,0,0,0);
170: if (vec->ops->assemblyend) {
171: (*vec->ops->assemblyend)(vec);
172: }
173: PetscLogEventEnd(VEC_AssemblyEnd,vec,0,0,0);
174: VecViewFromOptions(vec,NULL,"-vec_view");
175: return(0);
176: }
178: /*@
179: VecPointwiseMax - Computes the componentwise maximum w_i = max(x_i, y_i).
181: Logically Collective on Vec
183: Input Parameters:
184: . x, y - the vectors
186: Output Parameter:
187: . w - the result
189: Level: advanced
191: Notes:
192: any subset of the x, y, and w may be the same vector.
193: For complex numbers compares only the real part
195: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMin(), VecPointwiseMaxAbs(), VecMaxPointwiseDivide()
196: @*/
197: PetscErrorCode VecPointwiseMax(Vec w,Vec x,Vec y)
198: {
210: VecCheckSameSize(w,1,x,2);
211: VecCheckSameSize(w,1,y,3);
212: VecSetErrorIfLocked(w,1);
213: (*w->ops->pointwisemax)(w,x,y);
214: PetscObjectStateIncrease((PetscObject)w);
215: return(0);
216: }
218: /*@
219: VecPointwiseMin - Computes the componentwise minimum w_i = min(x_i, y_i).
221: Logically Collective on Vec
223: Input Parameters:
224: . x, y - the vectors
226: Output Parameter:
227: . w - the result
229: Level: advanced
231: Notes:
232: any subset of the x, y, and w may be the same vector.
233: For complex numbers compares only the real part
235: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMin(), VecPointwiseMaxAbs(), VecMaxPointwiseDivide()
236: @*/
237: PetscErrorCode VecPointwiseMin(Vec w,Vec x,Vec y)
238: {
250: VecCheckSameSize(w,1,x,2);
251: VecCheckSameSize(w,1,y,3);
252: VecSetErrorIfLocked(w,1);
253: (*w->ops->pointwisemin)(w,x,y);
254: PetscObjectStateIncrease((PetscObject)w);
255: return(0);
256: }
258: /*@
259: VecPointwiseMaxAbs - Computes the componentwise maximum of the absolute values w_i = max(abs(x_i), abs(y_i)).
261: Logically Collective on Vec
263: Input Parameters:
264: . x, y - the vectors
266: Output Parameter:
267: . w - the result
269: Level: advanced
271: Notes:
272: any subset of the x, y, and w may be the same vector.
274: .seealso: VecPointwiseDivide(), VecPointwiseMult(), VecPointwiseMin(), VecPointwiseMax(), VecMaxPointwiseDivide()
275: @*/
276: PetscErrorCode VecPointwiseMaxAbs(Vec w,Vec x,Vec y)
277: {
289: VecCheckSameSize(w,1,x,2);
290: VecCheckSameSize(w,1,y,3);
291: VecSetErrorIfLocked(w,1);
292: (*w->ops->pointwisemaxabs)(w,x,y);
293: PetscObjectStateIncrease((PetscObject)w);
294: return(0);
295: }
297: /*@
298: VecPointwiseDivide - Computes the componentwise division w = x/y.
300: Logically Collective on Vec
302: Input Parameters:
303: . x, y - the vectors
305: Output Parameter:
306: . w - the result
308: Level: advanced
310: Notes:
311: any subset of the x, y, and w may be the same vector.
313: .seealso: VecPointwiseMult(), VecPointwiseMax(), VecPointwiseMin(), VecPointwiseMaxAbs(), VecMaxPointwiseDivide()
314: @*/
315: PetscErrorCode VecPointwiseDivide(Vec w,Vec x,Vec y)
316: {
328: VecCheckSameSize(w,1,x,2);
329: VecCheckSameSize(w,1,y,3);
330: VecSetErrorIfLocked(w,1);
331: (*w->ops->pointwisedivide)(w,x,y);
332: PetscObjectStateIncrease((PetscObject)w);
333: return(0);
334: }
337: /*@
338: VecDuplicate - Creates a new vector of the same type as an existing vector.
340: Collective on Vec
342: Input Parameters:
343: . v - a vector to mimic
345: Output Parameter:
346: . newv - location to put new vector
348: Notes:
349: VecDuplicate() DOES NOT COPY the vector entries, but rather allocates storage
350: for the new vector. Use VecCopy() to copy a vector.
352: Use VecDestroy() to free the space. Use VecDuplicateVecs() to get several
353: vectors.
355: Level: beginner
357: .seealso: VecDestroy(), VecDuplicateVecs(), VecCreate(), VecCopy()
358: @*/
359: PetscErrorCode VecDuplicate(Vec v,Vec *newv)
360: {
367: (*v->ops->duplicate)(v,newv);
368: PetscObjectStateIncrease((PetscObject)*newv);
369: return(0);
370: }
372: /*@
373: VecDestroy - Destroys a vector.
375: Collective on Vec
377: Input Parameters:
378: . v - the vector
380: Level: beginner
382: .seealso: VecDuplicate(), VecDestroyVecs()
383: @*/
384: PetscErrorCode VecDestroy(Vec *v)
385: {
389: if (!*v) return(0);
391: if (--((PetscObject)(*v))->refct > 0) {*v = 0; return(0);}
393: PetscObjectSAWsViewOff((PetscObject)*v);
394: /* destroy the internal part */
395: if ((*v)->ops->destroy) {
396: (*(*v)->ops->destroy)(*v);
397: }
398: /* destroy the external/common part */
399: PetscLayoutDestroy(&(*v)->map);
400: PetscHeaderDestroy(v);
401: return(0);
402: }
404: /*@C
405: VecDuplicateVecs - Creates several vectors of the same type as an existing vector.
407: Collective on Vec
409: Input Parameters:
410: + m - the number of vectors to obtain
411: - v - a vector to mimic
413: Output Parameter:
414: . V - location to put pointer to array of vectors
416: Notes:
417: Use VecDestroyVecs() to free the space. Use VecDuplicate() to form a single
418: vector.
420: Fortran Note:
421: The Fortran interface is slightly different from that given below, it
422: requires one to pass in V a Vec (integer) array of size at least m.
423: See the Fortran chapter of the users manual and petsc/src/vec/vec/examples for details.
425: Level: intermediate
427: .seealso: VecDestroyVecs(), VecDuplicate(), VecCreate(), VecDuplicateVecsF90()
428: @*/
429: PetscErrorCode VecDuplicateVecs(Vec v,PetscInt m,Vec *V[])
430: {
437: (*v->ops->duplicatevecs)(v,m,V);
438: return(0);
439: }
441: /*@C
442: VecDestroyVecs - Frees a block of vectors obtained with VecDuplicateVecs().
444: Collective on Vec
446: Input Parameters:
447: + vv - pointer to pointer to array of vector pointers, if NULL no vectors are destroyed
448: - m - the number of vectors previously obtained, if zero no vectors are destroyed
450: Fortran Note:
451: The Fortran interface is slightly different from that given below.
452: See the Fortran chapter of the users manual
454: Level: intermediate
456: .seealso: VecDuplicateVecs(), VecDestroyVecsf90()
457: @*/
458: PetscErrorCode VecDestroyVecs(PetscInt m,Vec *vv[])
459: {
464: if (m < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to destroy negative number of vectors %D",m);
465: if (!m || !*vv) {*vv = NULL; return(0);}
468: (*(**vv)->ops->destroyvecs)(m,*vv);
469: *vv = NULL;
470: return(0);
471: }
473: /*@C
474: VecViewFromOptions - View from Options
476: Collective on Vec
478: Input Parameters:
479: + A - the vector
480: . obj - Optional object
481: - name - command line option
483: Level: intermediate
484: .seealso: Vec, VecView, PetscObjectViewFromOptions(), VecCreate()
485: @*/
486: PetscErrorCode VecViewFromOptions(Vec A,PetscObject obj,const char name[])
487: {
492: PetscObjectViewFromOptions((PetscObject)A,obj,name);
493: return(0);
494: }
496: /*@C
497: VecView - Views a vector object.
499: Collective on Vec
501: Input Parameters:
502: + vec - the vector
503: - viewer - an optional visualization context
505: Notes:
506: The available visualization contexts include
507: + PETSC_VIEWER_STDOUT_SELF - for sequential vectors
508: . PETSC_VIEWER_STDOUT_WORLD - for parallel vectors created on PETSC_COMM_WORLD
509: - PETSC_VIEWER_STDOUT_(comm) - for parallel vectors created on MPI communicator comm
511: You can change the format the vector is printed using the
512: option PetscViewerPushFormat().
514: The user can open alternative viewers with
515: + PetscViewerASCIIOpen() - Outputs vector to a specified file
516: . PetscViewerBinaryOpen() - Outputs vector in binary to a
517: specified file; corresponding input uses VecLoad()
518: . PetscViewerDrawOpen() - Outputs vector to an X window display
519: . PetscViewerSocketOpen() - Outputs vector to Socket viewer
520: - PetscViewerHDF5Open() - Outputs vector to HDF5 file viewer
522: The user can call PetscViewerPushFormat() to specify the output
523: format of ASCII printed objects (when using PETSC_VIEWER_STDOUT_SELF,
524: PETSC_VIEWER_STDOUT_WORLD and PetscViewerASCIIOpen). Available formats include
525: + PETSC_VIEWER_DEFAULT - default, prints vector contents
526: . PETSC_VIEWER_ASCII_MATLAB - prints vector contents in MATLAB format
527: . PETSC_VIEWER_ASCII_INDEX - prints vector contents, including indices of vector elements
528: - PETSC_VIEWER_ASCII_COMMON - prints vector contents, using a
529: format common among all vector types
531: Notes:
532: You can pass any number of vector objects, or other PETSc objects to the same viewer.
534: Notes for binary viewer:
535: If you pass multiple vectors to a binary viewer you can read them back in in the same order
536: with VecLoad().
538: If the blocksize of the vector is greater than one then you must provide a unique prefix to
539: the vector with PetscObjectSetOptionsPrefix((PetscObject)vec,"uniqueprefix"); BEFORE calling VecView() on the
540: vector to be stored and then set that same unique prefix on the vector that you pass to VecLoad(). The blocksize
541: information is stored in an ASCII file with the same name as the binary file plus a ".info" appended to the
542: filename. If you copy the binary file, make sure you copy the associated .info file with it.
544: See the manual page for VecLoad() on the exact format the binary viewer stores
545: the values in the file.
548: Notes for HDF5 Viewer:
549: The name of the Vec (given with PetscObjectSetName() is the name that is used
550: for the object in the HDF5 file. If you wish to store the same Vec into multiple
551: datasets in the same file (typically with different values), you must change its
552: name each time before calling the VecView(). To load the same vector,
553: the name of the Vec object passed to VecLoad() must be the same.
555: If the block size of the vector is greater than 1 then it is used as the first dimension in the HDF5 array.
556: If the function PetscViewerHDF5SetBaseDimension2()is called then even if the block size is one it will
557: be used as the first dimension in the HDF5 array (that is the HDF5 array will always be two dimensional)
558: See also PetscViewerHDF5SetTimestep() which adds an additional complication to reading and writing Vecs
559: with the HDF5 viewer.
561: Level: beginner
564: .seealso: PetscViewerASCIIOpen(), PetscViewerDrawOpen(), PetscDrawLGCreate(),
565: PetscViewerSocketOpen(), PetscViewerBinaryOpen(), VecLoad(), PetscViewerCreate(),
566: PetscRealView(), PetscScalarView(), PetscIntView(), PetscViewerHDF5SetTimestep()
567: @*/
568: PetscErrorCode VecView(Vec vec,PetscViewer viewer)
569: {
570: PetscErrorCode ierr;
571: PetscBool iascii;
572: PetscViewerFormat format;
573: PetscMPIInt size;
578: if (!viewer) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)vec),&viewer);}
580: PetscViewerGetFormat(viewer,&format);
581: MPI_Comm_size(PetscObjectComm((PetscObject)vec),&size);
582: if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) return(0);
584: if (vec->stash.n || vec->bstash.n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call VecAssemblyBegin/End() before viewing this vector");
586: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
587: if (iascii) {
588: PetscInt rows,bs;
590: PetscObjectPrintClassNamePrefixType((PetscObject)vec,viewer);
591: if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
592: PetscViewerASCIIPushTab(viewer);
593: VecGetSize(vec,&rows);
594: VecGetBlockSize(vec,&bs);
595: if (bs != 1) {
596: PetscViewerASCIIPrintf(viewer,"length=%D, bs=%D\n",rows,bs);
597: } else {
598: PetscViewerASCIIPrintf(viewer,"length=%D\n",rows);
599: }
600: PetscViewerASCIIPopTab(viewer);
601: }
602: }
603: VecLockReadPush(vec);
604: PetscLogEventBegin(VEC_View,vec,viewer,0,0);
605: if ((format == PETSC_VIEWER_NATIVE || format == PETSC_VIEWER_LOAD_BALANCE) && vec->ops->viewnative) {
606: (*vec->ops->viewnative)(vec,viewer);
607: } else {
608: (*vec->ops->view)(vec,viewer);
609: }
610: VecLockReadPop(vec);
611: PetscLogEventEnd(VEC_View,vec,viewer,0,0);
612: return(0);
613: }
615: #if defined(PETSC_USE_DEBUG)
616: #include <../src/sys/totalview/tv_data_display.h>
617: PETSC_UNUSED static int TV_display_type(const struct _p_Vec *v)
618: {
619: const PetscScalar *values;
620: char type[32];
621: PetscErrorCode ierr;
624: TV_add_row("Local rows", "int", &v->map->n);
625: TV_add_row("Global rows", "int", &v->map->N);
626: TV_add_row("Typename", TV_ascii_string_type, ((PetscObject)v)->type_name);
627: VecGetArrayRead((Vec)v,&values);
628: PetscSNPrintf(type,32,"double[%d]",v->map->n);
629: TV_add_row("values",type, values);
630: VecRestoreArrayRead((Vec)v,&values);
631: return TV_format_OK;
632: }
633: #endif
635: /*@
636: VecGetSize - Returns the global number of elements of the vector.
638: Not Collective
640: Input Parameter:
641: . x - the vector
643: Output Parameters:
644: . size - the global length of the vector
646: Level: beginner
648: .seealso: VecGetLocalSize()
649: @*/
650: PetscErrorCode VecGetSize(Vec x,PetscInt *size)
651: {
658: (*x->ops->getsize)(x,size);
659: return(0);
660: }
662: /*@
663: VecGetLocalSize - Returns the number of elements of the vector stored
664: in local memory.
666: Not Collective
668: Input Parameter:
669: . x - the vector
671: Output Parameter:
672: . size - the length of the local piece of the vector
674: Level: beginner
676: .seealso: VecGetSize()
677: @*/
678: PetscErrorCode VecGetLocalSize(Vec x,PetscInt *size)
679: {
686: (*x->ops->getlocalsize)(x,size);
687: return(0);
688: }
690: /*@C
691: VecGetOwnershipRange - Returns the range of indices owned by
692: this processor, assuming that the vectors are laid out with the
693: first n1 elements on the first processor, next n2 elements on the
694: second, etc. For certain parallel layouts this range may not be
695: well defined.
697: Not Collective
699: Input Parameter:
700: . x - the vector
702: Output Parameters:
703: + low - the first local element, pass in NULL if not interested
704: - high - one more than the last local element, pass in NULL if not interested
706: Note:
707: The high argument is one more than the last element stored locally.
709: Fortran: PETSC_NULL_INTEGER should be used instead of NULL
711: Level: beginner
714: .seealso: MatGetOwnershipRange(), MatGetOwnershipRanges(), VecGetOwnershipRanges()
715: @*/
716: PetscErrorCode VecGetOwnershipRange(Vec x,PetscInt *low,PetscInt *high)
717: {
723: if (low) *low = x->map->rstart;
724: if (high) *high = x->map->rend;
725: return(0);
726: }
728: /*@C
729: VecGetOwnershipRanges - Returns the range of indices owned by EACH processor,
730: assuming that the vectors are laid out with the
731: first n1 elements on the first processor, next n2 elements on the
732: second, etc. For certain parallel layouts this range may not be
733: well defined.
735: Not Collective
737: Input Parameter:
738: . x - the vector
740: Output Parameters:
741: . range - array of length size+1 with the start and end+1 for each process
743: Note:
744: The high argument is one more than the last element stored locally.
746: Fortran: You must PASS in an array of length size+1
748: Level: beginner
751: .seealso: MatGetOwnershipRange(), MatGetOwnershipRanges(), VecGetOwnershipRange()
752: @*/
753: PetscErrorCode VecGetOwnershipRanges(Vec x,const PetscInt *ranges[])
754: {
760: PetscLayoutGetRanges(x->map,ranges);
761: return(0);
762: }
764: /*@
765: VecSetOption - Sets an option for controling a vector's behavior.
767: Collective on Vec
769: Input Parameter:
770: + x - the vector
771: . op - the option
772: - flag - turn the option on or off
774: Supported Options:
775: + VEC_IGNORE_OFF_PROC_ENTRIES, which causes VecSetValues() to ignore
776: entries destined to be stored on a separate processor. This can be used
777: to eliminate the global reduction in the VecAssemblyXXXX() if you know
778: that you have only used VecSetValues() to set local elements
779: . VEC_IGNORE_NEGATIVE_INDICES, which means you can pass negative indices
780: in ix in calls to VecSetValues() or VecGetValues(). These rows are simply
781: ignored.
782: - VEC_SUBSET_OFF_PROC_ENTRIES, which causes VecAssemblyBegin() to assume that the off-process
783: entries will always be a subset (possibly equal) of the off-process entries set on the
784: first assembly which had a true VEC_SUBSET_OFF_PROC_ENTRIES and the vector has not
785: changed this flag afterwards. If this assembly is not such first assembly, then this
786: assembly can reuse the communication pattern setup in that first assembly, thus avoiding
787: a global reduction. Subsequent assemblies setting off-process values should use the same
788: InsertMode as the first assembly.
790: Developer Note:
791: The InsertMode restriction could be removed by packing the stash messages out of place.
793: Level: intermediate
795: @*/
796: PetscErrorCode VecSetOption(Vec x,VecOption op,PetscBool flag)
797: {
803: if (x->ops->setoption) {
804: (*x->ops->setoption)(x,op,flag);
805: }
806: return(0);
807: }
809: /* Default routines for obtaining and releasing; */
810: /* may be used by any implementation */
811: PetscErrorCode VecDuplicateVecs_Default(Vec w,PetscInt m,Vec *V[])
812: {
814: PetscInt i;
819: if (m <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"m must be > 0: m = %D",m);
820: PetscMalloc1(m,V);
821: for (i=0; i<m; i++) {VecDuplicate(w,*V+i);}
822: return(0);
823: }
825: PetscErrorCode VecDestroyVecs_Default(PetscInt m,Vec v[])
826: {
828: PetscInt i;
832: for (i=0; i<m; i++) {VecDestroy(&v[i]);}
833: PetscFree(v);
834: return(0);
835: }
837: /*@
838: VecResetArray - Resets a vector to use its default memory. Call this
839: after the use of VecPlaceArray().
841: Not Collective
843: Input Parameters:
844: . vec - the vector
846: Level: developer
848: .seealso: VecGetArray(), VecRestoreArray(), VecReplaceArray(), VecPlaceArray()
850: @*/
851: PetscErrorCode VecResetArray(Vec vec)
852: {
858: if (vec->ops->resetarray) {
859: (*vec->ops->resetarray)(vec);
860: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot reset array in this type of vector");
861: PetscObjectStateIncrease((PetscObject)vec);
862: return(0);
863: }
865: /*@C
866: VecLoad - Loads a vector that has been stored in binary or HDF5 format
867: with VecView().
869: Collective on PetscViewer
871: Input Parameters:
872: + vec - the newly loaded vector, this needs to have been created with VecCreate() or
873: some related function before a call to VecLoad().
874: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
875: HDF5 file viewer, obtained from PetscViewerHDF5Open()
877: Level: intermediate
879: Notes:
880: Defaults to the standard Seq or MPI Vec, if you want some other type of Vec call VecSetFromOptions()
881: before calling this.
883: The input file must contain the full global vector, as
884: written by the routine VecView().
886: If the type or size of vec is not set before a call to VecLoad, PETSc
887: sets the type and the local and global sizes. If type and/or
888: sizes are already set, then the same are used.
890: If using the binary viewer and the blocksize of the vector is greater than one then you must provide a unique prefix to
891: the vector with PetscObjectSetOptionsPrefix((PetscObject)vec,"uniqueprefix"); BEFORE calling VecView() on the
892: vector to be stored and then set that same unique prefix on the vector that you pass to VecLoad(). The blocksize
893: information is stored in an ASCII file with the same name as the binary file plus a ".info" appended to the
894: filename. If you copy the binary file, make sure you copy the associated .info file with it.
896: If using HDF5, you must assign the Vec the same name as was used in the Vec
897: that was stored in the file using PetscObjectSetName(). Otherwise you will
898: get the error message: "Cannot H5DOpen2() with Vec name NAMEOFOBJECT".
900: If the HDF5 file contains a two dimensional array the first dimension is treated as the block size
901: in loading the vector. Hence, for example, using Matlab notation h5create('vector.dat','/Test_Vec',[27 1]);
902: will load a vector of size 27 and block size 27 thus resulting in all 27 entries being on the first process of
903: vectors communicator and the rest of the processes having zero entries
905: Notes for advanced users when using the binary viewer:
906: Most users should not need to know the details of the binary storage
907: format, since VecLoad() and VecView() completely hide these details.
908: But for anyone who's interested, the standard binary vector storage
909: format is
910: .vb
911: PetscInt VEC_FILE_CLASSID
912: PetscInt number of rows
913: PetscScalar *values of all entries
914: .ve
916: In addition, PETSc automatically uses byte swapping to work on all machines; the files
917: are written ALWAYS using big-endian ordering. On small-endian machines the numbers
918: are converted to the small-endian format when they are read in from the file.
919: See PetscBinaryRead() and PetscBinaryWrite() to see how this may be done.
921: .seealso: PetscViewerBinaryOpen(), VecView(), MatLoad(), VecLoad()
922: @*/
923: PetscErrorCode VecLoad(Vec vec, PetscViewer viewer)
924: {
925: PetscErrorCode ierr;
926: PetscBool isbinary,ishdf5,isadios,isadios2;
927: PetscViewerFormat format;
933: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
934: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERHDF5,&ishdf5);
935: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERADIOS,&isadios);
936: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERADIOS,&isadios2);
937: if (!isbinary && !ishdf5 && !isadios && !isadios2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
939: VecSetErrorIfLocked(vec,1);
940: if (!((PetscObject)vec)->type_name && !vec->ops->create) {
941: VecSetType(vec, VECSTANDARD);
942: }
943: PetscLogEventBegin(VEC_Load,viewer,0,0,0);
944: PetscViewerGetFormat(viewer,&format);
945: if (format == PETSC_VIEWER_NATIVE && vec->ops->loadnative) {
946: (*vec->ops->loadnative)(vec,viewer);
947: } else {
948: (*vec->ops->load)(vec,viewer);
949: }
950: PetscLogEventEnd(VEC_Load,viewer,0,0,0);
951: return(0);
952: }
955: /*@
956: VecReciprocal - Replaces each component of a vector by its reciprocal.
958: Logically Collective on Vec
960: Input Parameter:
961: . vec - the vector
963: Output Parameter:
964: . vec - the vector reciprocal
966: Level: intermediate
968: .seealso: VecLog(), VecExp(), VecSqrtAbs()
970: @*/
971: PetscErrorCode VecReciprocal(Vec vec)
972: {
978: if (vec->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
979: if (!vec->ops->reciprocal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Vector does not support reciprocal operation");
980: VecSetErrorIfLocked(vec,1);
981: (*vec->ops->reciprocal)(vec);
982: PetscObjectStateIncrease((PetscObject)vec);
983: return(0);
984: }
986: /*@C
987: VecSetOperation - Allows user to set a vector operation.
989: Logically Collective on Vec
991: Input Parameters:
992: + vec - the vector
993: . op - the name of the operation
994: - f - the function that provides the operation.
996: Level: advanced
998: Usage:
999: $ PetscErrorCode userview(Vec,PetscViewer);
1000: $ VecCreateMPI(comm,m,M,&x);
1001: $ VecSetOperation(x,VECOP_VIEW,(void(*)(void))userview);
1003: Notes:
1004: See the file include/petscvec.h for a complete list of matrix
1005: operations, which all have the form VECOP_<OPERATION>, where
1006: <OPERATION> is the name (in all capital letters) of the
1007: user interface routine (e.g., VecView() -> VECOP_VIEW).
1009: This function is not currently available from Fortran.
1011: .seealso: VecCreate(), MatShellSetOperation()
1012: @*/
1013: PetscErrorCode VecSetOperation(Vec vec,VecOperation op, void (*f)(void))
1014: {
1017: if (op == VECOP_VIEW && !vec->ops->viewnative) {
1018: vec->ops->viewnative = vec->ops->view;
1019: } else if (op == VECOP_LOAD && !vec->ops->loadnative) {
1020: vec->ops->loadnative = vec->ops->load;
1021: }
1022: (((void(**)(void))vec->ops)[(int)op]) = f;
1023: return(0);
1024: }
1027: /*@
1028: VecStashSetInitialSize - sets the sizes of the vec-stash, that is
1029: used during the assembly process to store values that belong to
1030: other processors.
1032: Not Collective, different processes can have different size stashes
1034: Input Parameters:
1035: + vec - the vector
1036: . size - the initial size of the stash.
1037: - bsize - the initial size of the block-stash(if used).
1039: Options Database Keys:
1040: + -vecstash_initial_size <size> or <size0,size1,...sizep-1>
1041: - -vecstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1>
1043: Level: intermediate
1045: Notes:
1046: The block-stash is used for values set with VecSetValuesBlocked() while
1047: the stash is used for values set with VecSetValues()
1049: Run with the option -info and look for output of the form
1050: VecAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs.
1051: to determine the appropriate value, MM, to use for size and
1052: VecAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs.
1053: to determine the value, BMM to use for bsize
1056: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked(), VecStashView()
1058: @*/
1059: PetscErrorCode VecStashSetInitialSize(Vec vec,PetscInt size,PetscInt bsize)
1060: {
1065: VecStashSetInitialSize_Private(&vec->stash,size);
1066: VecStashSetInitialSize_Private(&vec->bstash,bsize);
1067: return(0);
1068: }
1070: /*@
1071: VecConjugate - Conjugates a vector.
1073: Logically Collective on Vec
1075: Input Parameters:
1076: . x - the vector
1078: Level: intermediate
1080: @*/
1081: PetscErrorCode VecConjugate(Vec x)
1082: {
1083: #if defined(PETSC_USE_COMPLEX)
1089: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
1090: VecSetErrorIfLocked(x,1);
1091: (*x->ops->conjugate)(x);
1092: /* we need to copy norms here */
1093: PetscObjectStateIncrease((PetscObject)x);
1094: return(0);
1095: #else
1096: return(0);
1097: #endif
1098: }
1100: /*@
1101: VecPointwiseMult - Computes the componentwise multiplication w = x*y.
1103: Logically Collective on Vec
1105: Input Parameters:
1106: . x, y - the vectors
1108: Output Parameter:
1109: . w - the result
1111: Level: advanced
1113: Notes:
1114: any subset of the x, y, and w may be the same vector.
1116: .seealso: VecPointwiseDivide(), VecPointwiseMax(), VecPointwiseMin(), VecPointwiseMaxAbs(), VecMaxPointwiseDivide()
1117: @*/
1118: PetscErrorCode VecPointwiseMult(Vec w, Vec x,Vec y)
1119: {
1131: VecCheckSameSize(w,1,x,2);
1132: VecCheckSameSize(w,2,y,3);
1133: VecSetErrorIfLocked(w,1);
1134: PetscLogEventBegin(VEC_PointwiseMult,x,y,w,0);
1135: (*w->ops->pointwisemult)(w,x,y);
1136: PetscLogEventEnd(VEC_PointwiseMult,x,y,w,0);
1137: PetscObjectStateIncrease((PetscObject)w);
1138: return(0);
1139: }
1141: /*@
1142: VecSetRandom - Sets all components of a vector to random numbers.
1144: Logically Collective on Vec
1146: Input Parameters:
1147: + x - the vector
1148: - rctx - the random number context, formed by PetscRandomCreate(), or NULL and
1149: it will create one internally.
1151: Output Parameter:
1152: . x - the vector
1154: Example of Usage:
1155: .vb
1156: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
1157: VecSetRandom(x,rctx);
1158: PetscRandomDestroy(rctx);
1159: .ve
1161: Level: intermediate
1164: .seealso: VecSet(), VecSetValues(), PetscRandomCreate(), PetscRandomDestroy()
1165: @*/
1166: PetscErrorCode VecSetRandom(Vec x,PetscRandom rctx)
1167: {
1169: PetscRandom randObj = NULL;
1175: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
1176: VecSetErrorIfLocked(x,1);
1178: if (!rctx) {
1179: MPI_Comm comm;
1180: PetscObjectGetComm((PetscObject)x,&comm);
1181: PetscRandomCreate(comm,&randObj);
1182: PetscRandomSetFromOptions(randObj);
1183: rctx = randObj;
1184: }
1186: PetscLogEventBegin(VEC_SetRandom,x,rctx,0,0);
1187: (*x->ops->setrandom)(x,rctx);
1188: PetscLogEventEnd(VEC_SetRandom,x,rctx,0,0);
1190: PetscRandomDestroy(&randObj);
1191: PetscObjectStateIncrease((PetscObject)x);
1192: return(0);
1193: }
1195: /*@
1196: VecZeroEntries - puts a 0.0 in each element of a vector
1198: Logically Collective on Vec
1200: Input Parameter:
1201: . vec - The vector
1203: Level: beginner
1205: .seealso: VecCreate(), VecSetOptionsPrefix(), VecSet(), VecSetValues()
1206: @*/
1207: PetscErrorCode VecZeroEntries(Vec vec)
1208: {
1212: VecSet(vec,0);
1213: return(0);
1214: }
1216: /*
1217: VecSetTypeFromOptions_Private - Sets the type of vector from user options. Defaults to a PETSc sequential vector on one
1218: processor and a PETSc MPI vector on more than one processor.
1220: Collective on Vec
1222: Input Parameter:
1223: . vec - The vector
1225: Level: intermediate
1227: .seealso: VecSetFromOptions(), VecSetType()
1228: */
1229: static PetscErrorCode VecSetTypeFromOptions_Private(PetscOptionItems *PetscOptionsObject,Vec vec)
1230: {
1231: PetscBool opt;
1232: VecType defaultType;
1233: char typeName[256];
1234: PetscMPIInt size;
1238: if (((PetscObject)vec)->type_name) defaultType = ((PetscObject)vec)->type_name;
1239: else {
1240: MPI_Comm_size(PetscObjectComm((PetscObject)vec), &size);
1241: if (size > 1) defaultType = VECMPI;
1242: else defaultType = VECSEQ;
1243: }
1245: VecRegisterAll();
1246: PetscOptionsFList("-vec_type","Vector type","VecSetType",VecList,defaultType,typeName,256,&opt);
1247: if (opt) {
1248: VecSetType(vec, typeName);
1249: } else {
1250: VecSetType(vec, defaultType);
1251: }
1252: return(0);
1253: }
1255: /*@
1256: VecSetFromOptions - Configures the vector from the options database.
1258: Collective on Vec
1260: Input Parameter:
1261: . vec - The vector
1263: Notes:
1264: To see all options, run your program with the -help option, or consult the users manual.
1265: Must be called after VecCreate() but before the vector is used.
1267: Level: beginner
1270: .seealso: VecCreate(), VecSetOptionsPrefix()
1271: @*/
1272: PetscErrorCode VecSetFromOptions(Vec vec)
1273: {
1279: PetscObjectOptionsBegin((PetscObject)vec);
1280: /* Handle vector type options */
1281: VecSetTypeFromOptions_Private(PetscOptionsObject,vec);
1283: /* Handle specific vector options */
1284: if (vec->ops->setfromoptions) {
1285: (*vec->ops->setfromoptions)(PetscOptionsObject,vec);
1286: }
1288: /* process any options handlers added with PetscObjectAddOptionsHandler() */
1289: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)vec);
1290: PetscOptionsEnd();
1291: return(0);
1292: }
1294: /*@
1295: VecSetSizes - Sets the local and global sizes, and checks to determine compatibility
1297: Collective on Vec
1299: Input Parameters:
1300: + v - the vector
1301: . n - the local size (or PETSC_DECIDE to have it set)
1302: - N - the global size (or PETSC_DECIDE)
1304: Notes:
1305: n and N cannot be both PETSC_DECIDE
1306: If one processor calls this with N of PETSC_DECIDE then all processors must, otherwise the program will hang.
1308: Level: intermediate
1310: .seealso: VecGetSize(), PetscSplitOwnership()
1311: @*/
1312: PetscErrorCode VecSetSizes(Vec v, PetscInt n, PetscInt N)
1313: {
1319: if (N >= 0 && n > N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local size %D cannot be larger than global size %D",n,N);
1320: if ((v->map->n >= 0 || v->map->N >= 0) && (v->map->n != n || v->map->N != N)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vector sizes to %D local %D global after previously setting them to %D local %D global",n,N,v->map->n,v->map->N);
1321: v->map->n = n;
1322: v->map->N = N;
1323: if (v->ops->create) {
1324: (*v->ops->create)(v);
1325: v->ops->create = 0;
1326: }
1327: return(0);
1328: }
1330: /*@
1331: VecSetBlockSize - Sets the blocksize for future calls to VecSetValuesBlocked()
1332: and VecSetValuesBlockedLocal().
1334: Logically Collective on Vec
1336: Input Parameter:
1337: + v - the vector
1338: - bs - the blocksize
1340: Notes:
1341: All vectors obtained by VecDuplicate() inherit the same blocksize.
1343: Level: advanced
1345: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMapping(), VecGetBlockSize()
1347: @*/
1348: PetscErrorCode VecSetBlockSize(Vec v,PetscInt bs)
1349: {
1355: PetscLayoutSetBlockSize(v->map,bs);
1356: v->bstash.bs = bs; /* use the same blocksize for the vec's block-stash */
1357: return(0);
1358: }
1360: /*@
1361: VecGetBlockSize - Gets the blocksize for the vector, i.e. what is used for VecSetValuesBlocked()
1362: and VecSetValuesBlockedLocal().
1364: Not Collective
1366: Input Parameter:
1367: . v - the vector
1369: Output Parameter:
1370: . bs - the blocksize
1372: Notes:
1373: All vectors obtained by VecDuplicate() inherit the same blocksize.
1375: Level: advanced
1377: .seealso: VecSetValuesBlocked(), VecSetLocalToGlobalMapping(), VecSetBlockSize()
1380: @*/
1381: PetscErrorCode VecGetBlockSize(Vec v,PetscInt *bs)
1382: {
1388: PetscLayoutGetBlockSize(v->map,bs);
1389: return(0);
1390: }
1392: /*@C
1393: VecSetOptionsPrefix - Sets the prefix used for searching for all
1394: Vec options in the database.
1396: Logically Collective on Vec
1398: Input Parameter:
1399: + v - the Vec context
1400: - prefix - the prefix to prepend to all option names
1402: Notes:
1403: A hyphen (-) must NOT be given at the beginning of the prefix name.
1404: The first character of all runtime options is AUTOMATICALLY the hyphen.
1406: Level: advanced
1408: .seealso: VecSetFromOptions()
1409: @*/
1410: PetscErrorCode VecSetOptionsPrefix(Vec v,const char prefix[])
1411: {
1416: PetscObjectSetOptionsPrefix((PetscObject)v,prefix);
1417: return(0);
1418: }
1420: /*@C
1421: VecAppendOptionsPrefix - Appends to the prefix used for searching for all
1422: Vec options in the database.
1424: Logically Collective on Vec
1426: Input Parameters:
1427: + v - the Vec context
1428: - prefix - the prefix to prepend to all option names
1430: Notes:
1431: A hyphen (-) must NOT be given at the beginning of the prefix name.
1432: The first character of all runtime options is AUTOMATICALLY the hyphen.
1434: Level: advanced
1436: .seealso: VecGetOptionsPrefix()
1437: @*/
1438: PetscErrorCode VecAppendOptionsPrefix(Vec v,const char prefix[])
1439: {
1444: PetscObjectAppendOptionsPrefix((PetscObject)v,prefix);
1445: return(0);
1446: }
1448: /*@C
1449: VecGetOptionsPrefix - Sets the prefix used for searching for all
1450: Vec options in the database.
1452: Not Collective
1454: Input Parameter:
1455: . v - the Vec context
1457: Output Parameter:
1458: . prefix - pointer to the prefix string used
1460: Notes:
1461: On the fortran side, the user should pass in a string 'prefix' of
1462: sufficient length to hold the prefix.
1464: Level: advanced
1466: .seealso: VecAppendOptionsPrefix()
1467: @*/
1468: PetscErrorCode VecGetOptionsPrefix(Vec v,const char *prefix[])
1469: {
1474: PetscObjectGetOptionsPrefix((PetscObject)v,prefix);
1475: return(0);
1476: }
1478: /*@
1479: VecSetUp - Sets up the internal vector data structures for the later use.
1481: Collective on Vec
1483: Input Parameters:
1484: . v - the Vec context
1486: Notes:
1487: For basic use of the Vec classes the user need not explicitly call
1488: VecSetUp(), since these actions will happen automatically.
1490: Level: advanced
1492: .seealso: VecCreate(), VecDestroy()
1493: @*/
1494: PetscErrorCode VecSetUp(Vec v)
1495: {
1496: PetscMPIInt size;
1501: if (v->map->n < 0 && v->map->N < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Sizes not set");
1502: if (!((PetscObject)v)->type_name) {
1503: MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);
1504: if (size == 1) {
1505: VecSetType(v, VECSEQ);
1506: } else {
1507: VecSetType(v, VECMPI);
1508: }
1509: }
1510: return(0);
1511: }
1513: /*
1514: These currently expose the PetscScalar/PetscReal in updating the
1515: cached norm. If we push those down into the implementation these
1516: will become independent of PetscScalar/PetscReal
1517: */
1519: /*@
1520: VecCopy - Copies a vector. y <- x
1522: Logically Collective on Vec
1524: Input Parameter:
1525: . x - the vector
1527: Output Parameter:
1528: . y - the copy
1530: Notes:
1531: For default parallel PETSc vectors, both x and y must be distributed in
1532: the same manner; local copies are done.
1534: Developer Notes:
1536: of the vectors to be sequential and one to be parallel so long as both have the same
1537: local sizes. This is used in some internal functions in PETSc.
1539: Level: beginner
1541: .seealso: VecDuplicate()
1542: @*/
1543: PetscErrorCode VecCopy(Vec x,Vec y)
1544: {
1545: PetscBool flgs[4];
1546: PetscReal norms[4] = {0.0,0.0,0.0,0.0};
1548: PetscInt i;
1555: if (x == y) return(0);
1556: VecCheckSameLocalSize(x,1,y,2);
1557: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
1558: VecSetErrorIfLocked(y,2);
1560: #if !defined(PETSC_USE_MIXED_PRECISION)
1561: for (i=0; i<4; i++) {
1562: PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],norms[i],flgs[i]);
1563: }
1564: #endif
1566: PetscLogEventBegin(VEC_Copy,x,y,0,0);
1567: #if defined(PETSC_USE_MIXED_PRECISION)
1568: extern PetscErrorCode VecGetArray(Vec,double**);
1569: extern PetscErrorCode VecRestoreArray(Vec,double**);
1570: extern PetscErrorCode VecGetArray(Vec,float**);
1571: extern PetscErrorCode VecRestoreArray(Vec,float**);
1572: extern PetscErrorCode VecGetArrayRead(Vec,const double**);
1573: extern PetscErrorCode VecRestoreArrayRead(Vec,const double**);
1574: extern PetscErrorCode VecGetArrayRead(Vec,const float**);
1575: extern PetscErrorCode VecRestoreArrayRead(Vec,const float**);
1576: if ((((PetscObject)x)->precision == PETSC_PRECISION_SINGLE) && (((PetscObject)y)->precision == PETSC_PRECISION_DOUBLE)) {
1577: PetscInt i,n;
1578: const float *xx;
1579: double *yy;
1580: VecGetArrayRead(x,&xx);
1581: VecGetArray(y,&yy);
1582: VecGetLocalSize(x,&n);
1583: for (i=0; i<n; i++) yy[i] = xx[i];
1584: VecRestoreArrayRead(x,&xx);
1585: VecRestoreArray(y,&yy);
1586: } else if ((((PetscObject)x)->precision == PETSC_PRECISION_DOUBLE) && (((PetscObject)y)->precision == PETSC_PRECISION_SINGLE)) {
1587: PetscInt i,n;
1588: float *yy;
1589: const double *xx;
1590: VecGetArrayRead(x,&xx);
1591: VecGetArray(y,&yy);
1592: VecGetLocalSize(x,&n);
1593: for (i=0; i<n; i++) yy[i] = (float) xx[i];
1594: VecRestoreArrayRead(x,&xx);
1595: VecRestoreArray(y,&yy);
1596: } else {
1597: (*x->ops->copy)(x,y);
1598: }
1599: #else
1600: (*x->ops->copy)(x,y);
1601: #endif
1603: PetscObjectStateIncrease((PetscObject)y);
1604: #if !defined(PETSC_USE_MIXED_PRECISION)
1605: for (i=0; i<4; i++) {
1606: if (flgs[i]) {
1607: PetscObjectComposedDataSetReal((PetscObject)y,NormIds[i],norms[i]);
1608: }
1609: }
1610: #endif
1612: PetscLogEventEnd(VEC_Copy,x,y,0,0);
1613: return(0);
1614: }
1616: /*@
1617: VecSwap - Swaps the vectors x and y.
1619: Logically Collective on Vec
1621: Input Parameters:
1622: . x, y - the vectors
1624: Level: advanced
1626: @*/
1627: PetscErrorCode VecSwap(Vec x,Vec y)
1628: {
1629: PetscReal normxs[4]={0.0,0.0,0.0,0.0},normys[4]={0.0,0.0,0.0,0.0};
1630: PetscBool flgxs[4],flgys[4];
1632: PetscInt i;
1640: VecCheckSameSize(x,1,y,2);
1641: if (x->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
1642: if (y->stash.insertmode != NOT_SET_VALUES) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled vector");
1643: VecSetErrorIfLocked(x,1);
1644: VecSetErrorIfLocked(y,2);
1646: PetscLogEventBegin(VEC_Swap,x,y,0,0);
1647: for (i=0; i<4; i++) {
1648: PetscObjectComposedDataGetReal((PetscObject)x,NormIds[i],normxs[i],flgxs[i]);
1649: PetscObjectComposedDataGetReal((PetscObject)y,NormIds[i],normys[i],flgys[i]);
1650: }
1651: (*x->ops->swap)(x,y);
1652: PetscObjectStateIncrease((PetscObject)x);
1653: PetscObjectStateIncrease((PetscObject)y);
1654: for (i=0; i<4; i++) {
1655: if (flgxs[i]) {
1656: PetscObjectComposedDataSetReal((PetscObject)y,NormIds[i],normxs[i]);
1657: }
1658: if (flgys[i]) {
1659: PetscObjectComposedDataSetReal((PetscObject)x,NormIds[i],normys[i]);
1660: }
1661: }
1662: PetscLogEventEnd(VEC_Swap,x,y,0,0);
1663: return(0);
1664: }
1666: /*
1667: VecStashViewFromOptions - Processes command line options to determine if/how an VecStash object is to be viewed.
1669: Collective on VecStash
1671: Input Parameters:
1672: + obj - the VecStash object
1673: . bobj - optional other object that provides the prefix
1674: - optionname - option to activate viewing
1676: Level: intermediate
1678: Developer Note: This cannot use PetscObjectViewFromOptions() because it takes a Vec as an argument but does not use VecView
1680: */
1681: PetscErrorCode VecStashViewFromOptions(Vec obj,PetscObject bobj,const char optionname[])
1682: {
1683: PetscErrorCode ierr;
1684: PetscViewer viewer;
1685: PetscBool flg;
1686: PetscViewerFormat format;
1687: char *prefix;
1690: prefix = bobj ? bobj->prefix : ((PetscObject)obj)->prefix;
1691: PetscOptionsGetViewer(PetscObjectComm((PetscObject)obj),((PetscObject)obj)->options,prefix,optionname,&viewer,&format,&flg);
1692: if (flg) {
1693: PetscViewerPushFormat(viewer,format);
1694: VecStashView(obj,viewer);
1695: PetscViewerPopFormat(viewer);
1696: PetscViewerDestroy(&viewer);
1697: }
1698: return(0);
1699: }
1701: /*@
1702: VecStashView - Prints the entries in the vector stash and block stash.
1704: Collective on Vec
1706: Input Parameters:
1707: + v - the vector
1708: - viewer - the viewer
1710: Level: advanced
1713: .seealso: VecSetBlockSize(), VecSetValues(), VecSetValuesBlocked()
1715: @*/
1716: PetscErrorCode VecStashView(Vec v,PetscViewer viewer)
1717: {
1719: PetscMPIInt rank;
1720: PetscInt i,j;
1721: PetscBool match;
1722: VecStash *s;
1723: PetscScalar val;
1730: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&match);
1731: if (!match) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Stash viewer only works with ASCII viewer not %s\n",((PetscObject)v)->type_name);
1732: PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);
1733: MPI_Comm_rank(PetscObjectComm((PetscObject)v),&rank);
1734: s = &v->bstash;
1736: /* print block stash */
1737: PetscViewerASCIIPushSynchronized(viewer);
1738: PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector Block stash size %D block size %D\n",rank,s->n,s->bs);
1739: for (i=0; i<s->n; i++) {
1740: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %D ",rank,s->idx[i]);
1741: for (j=0; j<s->bs; j++) {
1742: val = s->array[i*s->bs+j];
1743: #if defined(PETSC_USE_COMPLEX)
1744: PetscViewerASCIISynchronizedPrintf(viewer,"(%18.16e %18.16e) ",PetscRealPart(val),PetscImaginaryPart(val));
1745: #else
1746: PetscViewerASCIISynchronizedPrintf(viewer,"%18.16e ",val);
1747: #endif
1748: }
1749: PetscViewerASCIISynchronizedPrintf(viewer,"\n");
1750: }
1751: PetscViewerFlush(viewer);
1753: s = &v->stash;
1755: /* print basic stash */
1756: PetscViewerASCIISynchronizedPrintf(viewer,"[%d]Vector stash size %D\n",rank,s->n);
1757: for (i=0; i<s->n; i++) {
1758: val = s->array[i];
1759: #if defined(PETSC_USE_COMPLEX)
1760: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %D (%18.16e %18.16e) ",rank,s->idx[i],PetscRealPart(val),PetscImaginaryPart(val));
1761: #else
1762: PetscViewerASCIISynchronizedPrintf(viewer,"[%d] Element %D %18.16e\n",rank,s->idx[i],val);
1763: #endif
1764: }
1765: PetscViewerFlush(viewer);
1766: PetscViewerASCIIPopSynchronized(viewer);
1767: PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);
1768: return(0);
1769: }
1771: PetscErrorCode PetscOptionsGetVec(PetscOptions options,const char prefix[],const char key[],Vec v,PetscBool *set)
1772: {
1773: PetscInt i,N,rstart,rend;
1775: PetscScalar *xx;
1776: PetscReal *xreal;
1777: PetscBool iset;
1780: VecGetOwnershipRange(v,&rstart,&rend);
1781: VecGetSize(v,&N);
1782: PetscCalloc1(N,&xreal);
1783: PetscOptionsGetRealArray(options,prefix,key,xreal,&N,&iset);
1784: if (iset) {
1785: VecGetArray(v,&xx);
1786: for (i=rstart; i<rend; i++) xx[i-rstart] = xreal[i];
1787: VecRestoreArray(v,&xx);
1788: }
1789: PetscFree(xreal);
1790: if (set) *set = iset;
1791: return(0);
1792: }
1794: /*@
1795: VecGetLayout - get PetscLayout describing vector layout
1797: Not Collective
1799: Input Arguments:
1800: . x - the vector
1802: Output Arguments:
1803: . map - the layout
1805: Level: developer
1807: .seealso: VecGetSizes(), VecGetOwnershipRange(), VecGetOwnershipRanges()
1808: @*/
1809: PetscErrorCode VecGetLayout(Vec x,PetscLayout *map)
1810: {
1814: *map = x->map;
1815: return(0);
1816: }
1818: /*@
1819: VecSetLayout - set PetscLayout describing vector layout
1821: Not Collective
1823: Input Arguments:
1824: + x - the vector
1825: - map - the layout
1827: Notes:
1828: It is normally only valid to replace the layout with a layout known to be equivalent.
1830: Level: developer
1832: .seealso: VecGetLayout(), VecGetSizes(), VecGetOwnershipRange(), VecGetOwnershipRanges()
1833: @*/
1834: PetscErrorCode VecSetLayout(Vec x,PetscLayout map)
1835: {
1840: PetscLayoutReference(map,&x->map);
1841: return(0);
1842: }
1844: PetscErrorCode VecSetInf(Vec xin)
1845: {
1846: PetscInt i,n = xin->map->n;
1847: PetscScalar *xx;
1848: PetscScalar zero=0.0,one=1.0,inf=one/zero;
1852: VecGetArray(xin,&xx);
1853: for (i=0; i<n; i++) xx[i] = inf;
1854: VecRestoreArray(xin,&xx);
1855: return(0);
1856: }
1858: /*@
1859: VecBindToCPU - marks a vector to temporarily stay on the CPU and perform computations on the CPU
1861: Input Parameters:
1862: + v - the vector
1863: - flg - bind to the CPU if value of PETSC_TRUE
1865: Level: intermediate
1866: @*/
1867: PetscErrorCode VecBindToCPU(Vec v,PetscBool flg)
1868: {
1869: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1873: if (v->boundtocpu == flg) return(0);
1874: v->boundtocpu = flg;
1875: if (v->ops->bindtocpu) {
1876: (*v->ops->bindtocpu)(v,flg);
1877: }
1878: return(0);
1879: #else
1880: return 0;
1881: #endif
1882: }
1884: /*@C
1885: VecSetPinnedMemoryMin - Set the minimum data size for which pinned memory will be used for host (CPU) allocations.
1887: Logically Collective on Vec
1889: Input Parameters:
1890: + v - the vector
1891: - mbytes - minimum data size in bytes
1893: Options Database Keys:
1895: . -vec_pinned_memory_min <size> - minimum size (in bytes) for an allocation to use pinned memory on host.
1896: Note that this takes a PetscScalar, to accommodate large values;
1897: specifying -1 ensures that pinned memory will always be used.
1899: Level: developer
1901: .seealso: VecGetPinnedMemoryMin()
1902: @*/
1903: PetscErrorCode VecSetPinnedMemoryMin(Vec v,size_t mbytes)
1904: {
1905: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1907: v->minimum_bytes_pinned_memory = mbytes;
1908: return(0);
1909: #else
1910: return 0;
1911: #endif
1912: }
1914: /*@C
1915: VecGetPinnedMemoryMin - Get the minimum data size for which pinned memory will be used for host (CPU) allocations.
1917: Logically Collective on Vec
1919: Input Parameters:
1920: . v - the vector
1922: Output Parameters:
1923: . mbytes - minimum data size in bytes
1925: Level: developer
1927: .seealso: VecSetPinnedMemoryMin()
1928: @*/
1929: PetscErrorCode VecGetPinnedMemoryMin(Vec v,size_t *mbytes)
1930: {
1931: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
1933: *mbytes = v->minimum_bytes_pinned_memory;
1934: return(0);
1935: #else
1936: return 0;
1937: #endif
1938: }