Actual source code: petscvec.h

petsc-3.3-p7 2013-05-11
  1: /* 
  2:     Defines the vector component of PETSc. Vectors generally represent 
  3:   degrees of freedom for finite element/finite difference functions
  4:   on a grid. They have more mathematical structure then simple arrays.
  5: */

  7: #ifndef __PETSCVEC_H 
 9:  #include petscis.h


 12: /*S
 13:      Vec - Abstract PETSc vector object

 15:    Level: beginner

 17:   Concepts: field variables, unknowns, arrays

 19: .seealso:  VecCreate(), VecType, VecSetType()
 20: S*/
 21: typedef struct _p_Vec*         Vec;

 23: /*S
 24:      VecScatter - Object used to manage communication of data
 25:        between vectors in parallel. Manages both scatters and gathers

 27:    Level: beginner

 29:   Concepts: scatter

 31: .seealso:  VecScatterCreate(), VecScatterBegin(), VecScatterEnd()
 32: S*/
 33: typedef struct _p_VecScatter*  VecScatter;

 35: /*E
 36:   ScatterMode - Determines the direction of a scatter

 38:   Level: beginner

 40: .seealso: VecScatter, VecScatterBegin(), VecScatterEnd()
 41: E*/
 42: typedef enum {SCATTER_FORWARD=0, SCATTER_REVERSE=1, SCATTER_FORWARD_LOCAL=2, SCATTER_REVERSE_LOCAL=3, SCATTER_LOCAL=2} ScatterMode;

 44: /*MC
 45:     SCATTER_FORWARD - Scatters the values as dictated by the VecScatterCreate() call

 47:     Level: beginner

 49: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD_LOCAL,
 50:           SCATTER_REVERSE_LOCAL

 52: M*/

 54: /*MC
 55:     SCATTER_REVERSE - Moves the values in the opposite direction then the directions indicated in
 56:          in the VecScatterCreate()

 58:     Level: beginner

 60: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
 61:           SCATTER_REVERSE_LOCAL

 63: M*/

 65: /*MC
 66:     SCATTER_FORWARD_LOCAL - Scatters the values as dictated by the VecScatterCreate() call except NO parallel communication
 67:        is done. Any variables that have be moved between processes are ignored

 69:     Level: developer

 71: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD,
 72:           SCATTER_REVERSE_LOCAL

 74: M*/

 76: /*MC
 77:     SCATTER_REVERSE_LOCAL - Moves the values in the opposite direction then the directions indicated in
 78:          in the VecScatterCreate()  except NO parallel communication
 79:        is done. Any variables that have be moved between processes are ignored

 81:     Level: developer

 83: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
 84:           SCATTER_REVERSE

 86: M*/

 88: /*J
 89:     VecType - String with the name of a PETSc vector or the creation function
 90:        with an optional dynamic library name, for example
 91:        http://www.mcs.anl.gov/petsc/lib.a:myveccreate()

 93:    Level: beginner

 95: .seealso: VecSetType(), Vec
 96: J*/
 97: #define VecType char*
 98: #define VECSEQ         "seq"
 99: #define VECMPI         "mpi"
100: #define VECSTANDARD    "standard"   /* seq on one process and mpi on several */
101: #define VECSHARED      "shared"
102: #define VECSIEVE       "sieve"
103: #define VECSEQCUSP     "seqcusp"
104: #define VECMPICUSP     "mpicusp"
105: #define VECCUSP        "cusp"       /* seqcusp on one process and mpicusp on several */
106: #define VECNEST        "nest"
107: #define VECSEQPTHREAD  "seqpthread"
108: #define VECMPIPTHREAD  "mpipthread"
109: #define VECPTHREAD     "pthread"    /* seqpthread on one process and mpipthread on several */


112: /* Logging support */
113: #define    VEC_FILE_CLASSID 1211214
114: PETSC_EXTERN PetscClassId VEC_CLASSID;
115: PETSC_EXTERN PetscClassId VEC_SCATTER_CLASSID;


118: PETSC_EXTERN PetscErrorCode VecInitializePackage(const char[]);
119: PETSC_EXTERN PetscErrorCode VecFinalizePackage(void);

121: PETSC_EXTERN PetscErrorCode VecCreate(MPI_Comm,Vec*);
122: PETSC_EXTERN PetscErrorCode VecCreateSeq(MPI_Comm,PetscInt,Vec*);
123: PETSC_EXTERN PetscErrorCode VecCreateMPI(MPI_Comm,PetscInt,PetscInt,Vec*);
124: PETSC_EXTERN PetscErrorCode VecCreateSeqWithArray(MPI_Comm,PetscInt,PetscInt,const PetscScalar[],Vec*);
125: PETSC_EXTERN PetscErrorCode VecCreateMPIWithArray(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscScalar[],Vec*);
126: PETSC_EXTERN PetscErrorCode VecCreateShared(MPI_Comm,PetscInt,PetscInt,Vec*);
127: PETSC_EXTERN PetscErrorCode VecSetFromOptions(Vec);
128: PETSC_EXTERN PetscErrorCode VecSetUp(Vec);
129: PETSC_EXTERN PetscErrorCode VecDestroy(Vec*);
130: PETSC_EXTERN PetscErrorCode VecZeroEntries(Vec);
131: PETSC_EXTERN PetscErrorCode VecSetOptionsPrefix(Vec,const char[]);
132: PETSC_EXTERN PetscErrorCode VecAppendOptionsPrefix(Vec,const char[]);
133: PETSC_EXTERN PetscErrorCode VecGetOptionsPrefix(Vec,const char*[]);

135: PETSC_EXTERN PetscErrorCode VecSetSizes(Vec,PetscInt,PetscInt);

137: PETSC_EXTERN PetscErrorCode VecDotNorm2(Vec,Vec,PetscScalar*,PetscScalar*);
138: PETSC_EXTERN PetscErrorCode VecDot(Vec,Vec,PetscScalar*);
139: PETSC_EXTERN PetscErrorCode VecTDot(Vec,Vec,PetscScalar*);
140: PETSC_EXTERN PetscErrorCode VecMDot(Vec,PetscInt,const Vec[],PetscScalar[]);
141: PETSC_EXTERN PetscErrorCode VecMTDot(Vec,PetscInt,const Vec[],PetscScalar[]);
142: PETSC_EXTERN PetscErrorCode VecGetSubVector(Vec,IS,Vec*);
143: PETSC_EXTERN PetscErrorCode VecRestoreSubVector(Vec,IS,Vec*);

145: /*E
146:     NormType - determines what type of norm to compute

148:     Level: beginner

150: .seealso: VecNorm(), VecNormBegin(), VecNormEnd(), MatNorm()
151: E*/
152: typedef enum {NORM_1=0,NORM_2=1,NORM_FROBENIUS=2,NORM_INFINITY=3,NORM_1_AND_2=4} NormType;
153: PETSC_EXTERN const char *NormTypes[];
154: #define NORM_MAX NORM_INFINITY

156: /*MC
157:      NORM_1 - the one norm, ||v|| = sum_i | v_i |. ||A|| = max_j || v_*j ||, maximum column sum

159:    Level: beginner

161: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_2, NORM_FROBENIUS, 
162:            NORM_INFINITY, NORM_1_AND_2

164: M*/

166: /*MC
167:      NORM_2 - the two norm, ||v|| = sqrt(sum_i (v_i)^2) (vectors only)

169:    Level: beginner

171: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_FROBENIUS, 
172:            NORM_INFINITY, NORM_1_AND_2

174: M*/

176: /*MC
177:      NORM_FROBENIUS - ||A|| = sqrt(sum_ij (A_ij)^2), same as NORM_2 for vectors

179:    Level: beginner

181: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
182:            NORM_INFINITY, NORM_1_AND_2

184: M*/

186: /*MC
187:      NORM_INFINITY - ||v|| = max_i |v_i|. ||A|| = max_i || v_i* ||, maximum row sum

189:    Level: beginner

191: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
192:            NORM_FROBINIUS, NORM_1_AND_2

194: M*/

196: /*MC
197:      NORM_1_AND_2 - computes both the 1 and 2 norm of a vector

199:    Level: beginner

201: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
202:            NORM_FROBINIUS, NORM_INFINITY

204: M*/

206: /*MC
207:      NORM_MAX - see NORM_INFINITY

209:    Level: beginner

211: M*/

213: PETSC_EXTERN PetscErrorCode VecNorm(Vec,NormType,PetscReal *);
214: PETSC_EXTERN PetscErrorCode VecNormAvailable(Vec,NormType,PetscBool *,PetscReal *);
215: PETSC_EXTERN PetscErrorCode VecNormalize(Vec,PetscReal *);
216: PETSC_EXTERN PetscErrorCode VecSum(Vec,PetscScalar*);
217: PETSC_EXTERN PetscErrorCode VecMax(Vec,PetscInt*,PetscReal *);
218: PETSC_EXTERN PetscErrorCode VecMin(Vec,PetscInt*,PetscReal *);
219: PETSC_EXTERN PetscErrorCode VecScale(Vec,PetscScalar);
220: PETSC_EXTERN PetscErrorCode VecCopy(Vec,Vec);
221: PETSC_EXTERN PetscErrorCode VecSetRandom(Vec,PetscRandom);
222: PETSC_EXTERN PetscErrorCode VecSet(Vec,PetscScalar);
223: PETSC_EXTERN PetscErrorCode VecSwap(Vec,Vec);
224: PETSC_EXTERN PetscErrorCode VecAXPY(Vec,PetscScalar,Vec);
225: PETSC_EXTERN PetscErrorCode VecAXPBY(Vec,PetscScalar,PetscScalar,Vec);
226: PETSC_EXTERN PetscErrorCode VecMAXPY(Vec,PetscInt,const PetscScalar[],Vec[]);
227: PETSC_EXTERN PetscErrorCode VecAYPX(Vec,PetscScalar,Vec);
228: PETSC_EXTERN PetscErrorCode VecWAXPY(Vec,PetscScalar,Vec,Vec);
229: PETSC_EXTERN PetscErrorCode VecAXPBYPCZ(Vec,PetscScalar,PetscScalar,PetscScalar,Vec,Vec);
230: PETSC_EXTERN PetscErrorCode VecPointwiseMax(Vec,Vec,Vec);
231: PETSC_EXTERN PetscErrorCode VecPointwiseMaxAbs(Vec,Vec,Vec);
232: PETSC_EXTERN PetscErrorCode VecPointwiseMin(Vec,Vec,Vec);
233: PETSC_EXTERN PetscErrorCode VecPointwiseMult(Vec,Vec,Vec);
234: PETSC_EXTERN PetscErrorCode VecPointwiseDivide(Vec,Vec,Vec);
235: PETSC_EXTERN PetscErrorCode VecMaxPointwiseDivide(Vec,Vec,PetscReal*);
236: PETSC_EXTERN PetscErrorCode VecShift(Vec,PetscScalar);
237: PETSC_EXTERN PetscErrorCode VecReciprocal(Vec);
238: PETSC_EXTERN PetscErrorCode VecPermute(Vec, IS, PetscBool );
239: PETSC_EXTERN PetscErrorCode VecSqrtAbs(Vec);
240: PETSC_EXTERN PetscErrorCode VecLog(Vec);
241: PETSC_EXTERN PetscErrorCode VecExp(Vec);
242: PETSC_EXTERN PetscErrorCode VecAbs(Vec);
243: PETSC_EXTERN PetscErrorCode VecDuplicate(Vec,Vec*);
244: PETSC_EXTERN PetscErrorCode VecDuplicateVecs(Vec,PetscInt,Vec*[]);
245: PETSC_EXTERN PetscErrorCode VecDestroyVecs(PetscInt, Vec*[]);
246: PETSC_EXTERN PetscErrorCode VecStrideNormAll(Vec,NormType,PetscReal[]);
247: PETSC_EXTERN PetscErrorCode VecStrideMaxAll(Vec,PetscInt [],PetscReal []);
248: PETSC_EXTERN PetscErrorCode VecStrideMinAll(Vec,PetscInt [],PetscReal []);
249: PETSC_EXTERN PetscErrorCode VecStrideScaleAll(Vec,const PetscScalar[]);

251: PETSC_EXTERN PetscErrorCode VecStrideNorm(Vec,PetscInt,NormType,PetscReal*);
252: PETSC_EXTERN PetscErrorCode VecStrideMax(Vec,PetscInt,PetscInt *,PetscReal *);
253: PETSC_EXTERN PetscErrorCode VecStrideMin(Vec,PetscInt,PetscInt *,PetscReal *);
254: PETSC_EXTERN PetscErrorCode VecStrideScale(Vec,PetscInt,PetscScalar);
255: PETSC_EXTERN PetscErrorCode VecStrideSet(Vec,PetscInt,PetscScalar);


258: PETSC_EXTERN PetscErrorCode VecStrideGather(Vec,PetscInt,Vec,InsertMode);
259: PETSC_EXTERN PetscErrorCode VecStrideScatter(Vec,PetscInt,Vec,InsertMode);
260: PETSC_EXTERN PetscErrorCode VecStrideGatherAll(Vec,Vec[],InsertMode);
261: PETSC_EXTERN PetscErrorCode VecStrideScatterAll(Vec[],Vec,InsertMode);

263: PETSC_EXTERN PetscErrorCode VecSetValues(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
264: PETSC_EXTERN PetscErrorCode VecGetValues(Vec,PetscInt,const PetscInt[],PetscScalar[]);
265: PETSC_EXTERN PetscErrorCode VecAssemblyBegin(Vec);
266: PETSC_EXTERN PetscErrorCode VecAssemblyEnd(Vec);
267: PETSC_EXTERN PetscErrorCode VecStashSetInitialSize(Vec,PetscInt,PetscInt);
268: PETSC_EXTERN PetscErrorCode VecStashView(Vec,PetscViewer);
269: PETSC_EXTERN PetscErrorCode VecStashGetInfo(Vec,PetscInt*,PetscInt*,PetscInt*,PetscInt*);

271: /*MC
272:    VecSetValue - Set a single entry into a vector.

274:    Synopsis:
275:    PetscErrorCode VecSetValue(Vec v,int row,PetscScalar value, InsertMode mode);

277:    Not Collective

279:    Input Parameters:
280: +  v - the vector
281: .  row - the row location of the entry
282: .  value - the value to insert
283: -  mode - either INSERT_VALUES or ADD_VALUES

285:    Notes:
286:    For efficiency one should use VecSetValues() and set several or 
287:    many values simultaneously if possible.

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

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

294:    Level: beginner

296: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(), VecSetValueLocal()
297: M*/
298: PETSC_STATIC_INLINE PetscErrorCode VecSetValue(Vec v,PetscInt i,PetscScalar va,InsertMode mode) {return VecSetValues(v,1,&i,&va,mode);}


301: PETSC_EXTERN PetscErrorCode VecSetBlockSize(Vec,PetscInt);
302: PETSC_EXTERN PetscErrorCode VecGetBlockSize(Vec,PetscInt*);
303: PETSC_EXTERN PetscErrorCode VecSetValuesBlocked(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

305: /* Dynamic creation and loading functions */
306: PETSC_EXTERN PetscFList VecList;
307: PETSC_EXTERN PetscBool VecRegisterAllCalled;
308: PETSC_EXTERN PetscErrorCode VecSetType(Vec, const VecType);
309: PETSC_EXTERN PetscErrorCode VecGetType(Vec, const VecType *);
310: PETSC_EXTERN PetscErrorCode VecRegister(const char[],const char[],const char[],PetscErrorCode (*)(Vec));
311: PETSC_EXTERN PetscErrorCode VecRegisterAll(const char []);
312: PETSC_EXTERN PetscErrorCode VecRegisterDestroy(void);

314: /*MC
315:   VecRegisterDynamic - Adds a new vector component implementation

317:   Synopsis:
318:   PetscErrorCode VecRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(Vec))

320:   Not Collective

322:   Input Parameters:
323: + name        - The name of a new user-defined creation routine
324: . path        - The path (either absolute or relative) of the library containing this routine
325: . func_name   - The name of routine to create method context
326: - create_func - The creation routine itself

328:   Notes:
329:   VecRegisterDynamic() may be called multiple times to add several user-defined vectors

331:   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.

333:   Sample usage:
334: .vb
335:     VecRegisterDynamic("my_vec","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyVectorCreate", MyVectorCreate);
336: .ve

338:   Then, your vector type can be chosen with the procedural interface via
339: .vb
340:     VecCreate(MPI_Comm, Vec *);
341:     VecSetType(Vec,"my_vector_name");
342: .ve
343:    or at runtime via the option
344: .vb
345:     -vec_type my_vector_name
346: .ve

348:   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
349:          If your function is not being put into a shared library then use VecRegister() instead
350:         
351:   Level: advanced

353: .keywords: Vec, register
354: .seealso: VecRegisterAll(), VecRegisterDestroy(), VecRegister()
355: M*/
356: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
357: #define VecRegisterDynamic(a,b,c,d) VecRegister(a,b,c,0)
358: #else
359: #define VecRegisterDynamic(a,b,c,d) VecRegister(a,b,c,d)
360: #endif


363: PETSC_EXTERN PetscErrorCode VecScatterCreate(Vec,IS,Vec,IS,VecScatter *);
364: PETSC_EXTERN PetscErrorCode VecScatterCreateEmpty(MPI_Comm,VecScatter *);
365: PETSC_EXTERN PetscErrorCode VecScatterCreateLocal(VecScatter,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],PetscInt);
366: PETSC_EXTERN PetscErrorCode VecScatterBegin(VecScatter,Vec,Vec,InsertMode,ScatterMode);
367: PETSC_EXTERN PetscErrorCode VecScatterEnd(VecScatter,Vec,Vec,InsertMode,ScatterMode);
368: PETSC_EXTERN PetscErrorCode VecScatterDestroy(VecScatter*);
369: PETSC_EXTERN PetscErrorCode VecScatterCopy(VecScatter,VecScatter *);
370: PETSC_EXTERN PetscErrorCode VecScatterView(VecScatter,PetscViewer);
371: PETSC_EXTERN PetscErrorCode VecScatterRemap(VecScatter,PetscInt *,PetscInt*);
372: PETSC_EXTERN PetscErrorCode VecScatterGetMerged(VecScatter,PetscBool *);

374: PETSC_EXTERN PetscErrorCode VecGetArray4d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar****[]);
375: PETSC_EXTERN PetscErrorCode VecRestoreArray4d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar****[]);
376: PETSC_EXTERN PetscErrorCode VecGetArray3d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar***[]);
377: PETSC_EXTERN PetscErrorCode VecRestoreArray3d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar***[]);
378: PETSC_EXTERN PetscErrorCode VecGetArray2d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar**[]);
379: PETSC_EXTERN PetscErrorCode VecRestoreArray2d(Vec,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar**[]);
380: PETSC_EXTERN PetscErrorCode VecGetArray1d(Vec,PetscInt,PetscInt,PetscScalar *[]);
381: PETSC_EXTERN PetscErrorCode VecRestoreArray1d(Vec,PetscInt,PetscInt,PetscScalar *[]);

383: PETSC_EXTERN PetscErrorCode VecPlaceArray(Vec,const PetscScalar[]);
384: PETSC_EXTERN PetscErrorCode VecResetArray(Vec);
385: PETSC_EXTERN PetscErrorCode VecReplaceArray(Vec,const PetscScalar[]);
386: PETSC_EXTERN PetscErrorCode VecGetArrays(const Vec[],PetscInt,PetscScalar**[]);
387: PETSC_EXTERN PetscErrorCode VecRestoreArrays(const Vec[],PetscInt,PetscScalar**[]);

389: PETSC_EXTERN PetscErrorCode VecView(Vec,PetscViewer);
390: PETSC_EXTERN PetscErrorCode VecViewFromOptions(Vec, const char *);
391: PETSC_EXTERN PetscErrorCode VecEqual(Vec,Vec,PetscBool *);
392: PETSC_EXTERN PetscErrorCode VecLoad(Vec, PetscViewer);

394: PETSC_EXTERN PetscErrorCode VecGetSize(Vec,PetscInt*);
395: PETSC_EXTERN PetscErrorCode VecGetLocalSize(Vec,PetscInt*);
396: PETSC_EXTERN PetscErrorCode VecGetOwnershipRange(Vec,PetscInt*,PetscInt*);
397: PETSC_EXTERN PetscErrorCode VecGetOwnershipRanges(Vec,const PetscInt *[]);

399: PETSC_EXTERN PetscErrorCode VecSetLocalToGlobalMapping(Vec,ISLocalToGlobalMapping);
400: PETSC_EXTERN PetscErrorCode VecSetValuesLocal(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);

402: /*MC
403:    VecSetValueLocal - Set a single entry into a vector using the local numbering

405:    Synopsis:
406:    PetscErrorCode VecSetValueLocal(Vec v,int row,PetscScalar value, InsertMode mode);

408:    Not Collective

410:    Input Parameters:
411: +  v - the vector
412: .  row - the row location of the entry
413: .  value - the value to insert
414: -  mode - either INSERT_VALUES or ADD_VALUES

416:    Notes:
417:    For efficiency one should use VecSetValues() and set several or 
418:    many values simultaneously if possible.

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

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

425:    Level: beginner

427: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(), VecSetValue()
428: M*/
429: PETSC_STATIC_INLINE PetscErrorCode VecSetValueLocal(Vec v,PetscInt i,PetscScalar va,InsertMode mode) {return VecSetValuesLocal(v,1,&i,&va,mode);}

431: PETSC_EXTERN PetscErrorCode VecSetLocalToGlobalMappingBlock(Vec,ISLocalToGlobalMapping);
432: PETSC_EXTERN PetscErrorCode VecSetValuesBlockedLocal(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
433: PETSC_EXTERN PetscErrorCode VecGetLocalToGlobalMappingBlock(Vec,ISLocalToGlobalMapping*);
434: PETSC_EXTERN PetscErrorCode VecGetLocalToGlobalMapping(Vec,ISLocalToGlobalMapping*);

436: PETSC_EXTERN PetscErrorCode VecDotBegin(Vec,Vec,PetscScalar *);
437: PETSC_EXTERN PetscErrorCode VecDotEnd(Vec,Vec,PetscScalar *);
438: PETSC_EXTERN PetscErrorCode VecTDotBegin(Vec,Vec,PetscScalar *);
439: PETSC_EXTERN PetscErrorCode VecTDotEnd(Vec,Vec,PetscScalar *);
440: PETSC_EXTERN PetscErrorCode VecNormBegin(Vec,NormType,PetscReal *);
441: PETSC_EXTERN PetscErrorCode VecNormEnd(Vec,NormType,PetscReal *);

443: PETSC_EXTERN PetscErrorCode VecMDotBegin(Vec,PetscInt,const Vec[],PetscScalar[]);
444: PETSC_EXTERN PetscErrorCode VecMDotEnd(Vec,PetscInt,const Vec[],PetscScalar[]);
445: PETSC_EXTERN PetscErrorCode VecMTDotBegin(Vec,PetscInt,const Vec[],PetscScalar[]);
446: PETSC_EXTERN PetscErrorCode VecMTDotEnd(Vec,PetscInt,const Vec[],PetscScalar[]);
447: PETSC_EXTERN PetscErrorCode PetscCommSplitReductionBegin(MPI_Comm);


450: #if defined(PETSC_USE_DEBUG)
451: #define VecValidValues(vec,argnum,input) do {                           \
452:     PetscErrorCode     _ierr;                                           \
453:     PetscInt          _n,_i;                                            \
454:     const PetscScalar *_x;                                              \
455:                                                                         \
456:     if (vec->petscnative || vec->ops->getarray) {                       \
457:       _VecGetLocalSize(vec,&_n);CHKERRQ(_ierr);                  \
458:       _VecGetArrayRead(vec,&_x);CHKERRQ(_ierr);                  \
459:       for (_i=0; _i<_n; _i++) {                                         \
460:         if (input) {                                                    \
461:           if (PetscIsInfOrNanScalar(_x[_i])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FP,"Vec entry at local location %D is not-a-number or infinite at beginning of function: Parameter number %d",_i,argnum); \
462:         } else {                                                        \
463:           if (PetscIsInfOrNanScalar(_x[_i])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FP,"Vec entry at local location %D is not-a-number or infinite at end of function: Parameter number %d",_i,argnum); \
464:         }                                                               \
465:       }                                                                 \
466:       _VecRestoreArrayRead(vec,&_x);CHKERRQ(_ierr);              \
467:     }                                                                   \
468:   } while (0)
469: #else
470: #define VecValidValues(vec,argnum,input)
471: #endif


474: typedef enum {VEC_IGNORE_OFF_PROC_ENTRIES,VEC_IGNORE_NEGATIVE_INDICES} VecOption;
475: PETSC_EXTERN PetscErrorCode VecSetOption(Vec,VecOption,PetscBool );

477: /*
478:    Expose VecGetArray()/VecRestoreArray() to users. Allows this to work without any function
479:    call overhead on any 'native' Vecs.
480: */

482: #include "petsc-private/vecimpl.h"

484: PETSC_EXTERN PetscErrorCode VecContourScale(Vec,PetscReal,PetscReal);

486: /*
487:     These numbers need to match the entries in 
488:   the function table in vecimpl.h
489: */
490: typedef enum { VECOP_VIEW = 33, VECOP_LOAD = 41, VECOP_DUPLICATE = 0} VecOperation;
491: PETSC_EXTERN PetscErrorCode VecSetOperation(Vec,VecOperation,void(*)(void));

493: /*
494:      Routines for dealing with ghosted vectors:
495:   vectors with ghost elements at the end of the array.
496: */
497: PETSC_EXTERN PetscErrorCode VecMPISetGhost(Vec,PetscInt,const PetscInt[]);
498: PETSC_EXTERN PetscErrorCode VecCreateGhost(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Vec*);
499: PETSC_EXTERN PetscErrorCode VecCreateGhostWithArray(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],Vec*);
500: PETSC_EXTERN PetscErrorCode VecCreateGhostBlock(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Vec*);
501: PETSC_EXTERN PetscErrorCode VecCreateGhostBlockWithArray(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],Vec*);
502: PETSC_EXTERN PetscErrorCode VecGhostGetLocalForm(Vec,Vec*);
503: PETSC_EXTERN PetscErrorCode VecGhostRestoreLocalForm(Vec,Vec*);
504: PETSC_EXTERN PetscErrorCode VecGhostUpdateBegin(Vec,InsertMode,ScatterMode);
505: PETSC_EXTERN PetscErrorCode VecGhostUpdateEnd(Vec,InsertMode,ScatterMode);

507: PETSC_EXTERN PetscErrorCode VecConjugate(Vec);

509: PETSC_EXTERN PetscErrorCode VecScatterCreateToAll(Vec,VecScatter*,Vec*);
510: PETSC_EXTERN PetscErrorCode VecScatterCreateToZero(Vec,VecScatter*,Vec*);

512: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaGetVector(PetscViewer, Vec);
513: PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutVector(PetscViewer, Vec);

515: /*S
516:      Vecs - Collection of vectors where the data for the vectors is stored in 
517:             one contiguous memory

519:    Level: advanced

521:    Notes:
522:     Temporary construct for handling multiply right hand side solves

524:     This is faked by storing a single vector that has enough array space for 
525:     n vectors

527:   Concepts: parallel decomposition

529: S*/
530:         struct _n_Vecs  {PetscInt n; Vec v;};
531: typedef struct _n_Vecs* Vecs;
532: #define VecsDestroy(x)            (VecDestroy(&(x)->v)         || PetscFree(x))
533: #define VecsCreateSeq(comm,p,m,x) (PetscNew(struct _n_Vecs,x) || VecCreateSeq(comm,p*m,&(*(x))->v) || (-1 == ((*(x))->n = (m))))
534: #define VecsCreateSeqWithArray(comm,p,m,a,x) (PetscNew(struct _n_Vecs,x) || VecCreateSeqWithArray(comm,1,p*m,a,&(*(x))->v) || (-1 == ((*(x))->n = (m))))
535: #define VecsDuplicate(x,y)        (PetscNew(struct _n_Vecs,y) || VecDuplicate(x->v,&(*(y))->v) || (-1 == ((*(y))->n = (x)->n)))

537: #if defined(PETSC_HAVE_CUSP)
538: typedef struct _p_PetscCUSPIndices* PetscCUSPIndices;
539: PETSC_EXTERN PetscErrorCode PetscCUSPIndicesCreate(PetscInt, PetscInt*,PetscInt, PetscInt*,PetscCUSPIndices*);
540: PETSC_EXTERN PetscErrorCode PetscCUSPIndicesDestroy(PetscCUSPIndices*);
541: PETSC_EXTERN PetscErrorCode VecCUSPCopyToGPUSome_Public(Vec,PetscCUSPIndices);
542: PETSC_EXTERN PetscErrorCode VecCUSPCopyFromGPUSome_Public(Vec,PetscCUSPIndices);

544: #if defined(PETSC_HAVE_TXPETSCGPU)
545: PETSC_EXTERN PetscErrorCode VecCUSPResetIndexBuffersFlagsGPU_Public(PetscCUSPIndices);
546: PETSC_EXTERN PetscErrorCode VecCUSPCopySomeToContiguousBufferGPU_Public(Vec,PetscCUSPIndices);
547: PETSC_EXTERN PetscErrorCode VecCUSPCopySomeFromContiguousBufferGPU_Public(Vec,PetscCUSPIndices);
548: PETSC_EXTERN PetscErrorCode VecScatterInitializeForGPU(VecScatter,Vec,ScatterMode);
549: PETSC_EXTERN PetscErrorCode VecScatterFinalizeForGPU(VecScatter);
550: #endif

552: PETSC_EXTERN PetscErrorCode VecCreateSeqCUSP(MPI_Comm,PetscInt,Vec*);
553: PETSC_EXTERN PetscErrorCode VecCreateMPICUSP(MPI_Comm,PetscInt,PetscInt,Vec*);
554: #endif

556: #if defined(PETSC_HAVE_PTHREADCLASSES)
557: PETSC_EXTERN PetscErrorCode VecSetNThreads(Vec,PetscInt);
558: PETSC_EXTERN PetscErrorCode VecGetNThreads(Vec,PetscInt*);
559: PETSC_EXTERN PetscErrorCode VecGetThreadOwnershipRange(Vec,PetscInt,PetscInt*,PetscInt*);
560: PETSC_EXTERN PetscErrorCode VecSetThreadAffinities(Vec,const PetscInt[]);
561: PETSC_EXTERN PetscErrorCode VecCreateSeqPThread(MPI_Comm,PetscInt,PetscInt,PetscInt[],Vec*);
562: PETSC_EXTERN PetscErrorCode VecCreateMPIPThread(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],Vec*);
563: #endif

565: PETSC_EXTERN PetscErrorCode VecNestGetSubVecs(Vec,PetscInt*,Vec**);
566: PETSC_EXTERN PetscErrorCode VecNestGetSubVec(Vec,PetscInt,Vec*);
567: PETSC_EXTERN PetscErrorCode VecNestSetSubVecs(Vec,PetscInt,PetscInt*,Vec*);
568: PETSC_EXTERN PetscErrorCode VecNestSetSubVec(Vec,PetscInt,Vec);
569: PETSC_EXTERN PetscErrorCode VecCreateNest(MPI_Comm,PetscInt,IS*,Vec*,Vec*);
570: PETSC_EXTERN PetscErrorCode VecNestGetSize(Vec,PetscInt*);

572: #endif