Actual source code: petscvec.h
  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
 13: /*S
 14:      Vec - Abstract PETSc vector object
 16:    Level: beginner
 18:   Concepts: field variables, unknowns, arrays
 20: .seealso:  VecCreate(), VecType, VecSetType()
 21: S*/
 22: typedef struct _p_Vec*         Vec;
 24: /*S
 25:      VecScatter - Object used to manage communication of data
 26:        between vectors in parallel. Manages both scatters and gathers
 28:    Level: beginner
 30:   Concepts: scatter
 32: .seealso:  VecScatterCreate(), VecScatterBegin(), VecScatterEnd()
 33: S*/
 34: typedef struct _p_VecScatter*  VecScatter;
 36: /*E
 37:   ScatterMode - Determines the direction of a scatter
 39:   Level: beginner
 41: .seealso: VecScatter, VecScatterBegin(), VecScatterEnd()
 42: E*/
 43: typedef enum {SCATTER_FORWARD=0, SCATTER_REVERSE=1, SCATTER_FORWARD_LOCAL=2, SCATTER_REVERSE_LOCAL=3, SCATTER_LOCAL=2} ScatterMode;
 45: /*MC
 46:     SCATTER_FORWARD - Scatters the values as dictated by the VecScatterCreate() call
 48:     Level: beginner
 50: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD_LOCAL,
 51:           SCATTER_REVERSE_LOCAL
 53: M*/
 55: /*MC
 56:     SCATTER_REVERSE - Moves the values in the opposite direction then the directions indicated in
 57:          in the VecScatterCreate()
 59:     Level: beginner
 61: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
 62:           SCATTER_REVERSE_LOCAL
 64: M*/
 66: /*MC
 67:     SCATTER_FORWARD_LOCAL - Scatters the values as dictated by the VecScatterCreate() call except NO parallel communication
 68:        is done. Any variables that have be moved between processes are ignored
 70:     Level: developer
 72: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_REVERSE, SCATTER_FORWARD,
 73:           SCATTER_REVERSE_LOCAL
 75: M*/
 77: /*MC
 78:     SCATTER_REVERSE_LOCAL - Moves the values in the opposite direction then the directions indicated in
 79:          in the VecScatterCreate()  except NO parallel communication
 80:        is done. Any variables that have be moved between processes are ignored
 82:     Level: developer
 84: .seealso: VecScatter, ScatterMode, VecScatterCreate(), VecScatterBegin(), VecScatterEnd(), SCATTER_FORWARD, SCATTER_FORWARD_LOCAL,
 85:           SCATTER_REVERSE
 87: M*/
 89: /*J
 90:     VecType - String with the name of a PETSc vector or the creation function
 91:        with an optional dynamic library name, for example
 92:        http://www.mcs.anl.gov/petsc/lib.a:myveccreate()
 94:    Level: beginner
 96: .seealso: VecSetType(), Vec
 97: J*/
 98: #define VecType char*
 99: #define VECSEQ         "seq"
100: #define VECMPI         "mpi"
101: #define VECSTANDARD    "standard"   /* seq on one process and mpi on several */
102: #define VECSHARED      "shared"
103: #define VECSIEVE       "sieve"
104: #define VECSEQCUSP     "seqcusp"
105: #define VECMPICUSP     "mpicusp"
106: #define VECCUSP        "cusp"       /* seqcusp on one process and mpicusp on several */
107: #define VECNEST        "nest"
108: #define VECSEQPTHREAD  "seqpthread"
109: #define VECPTHREAD     "pthread"    /* seqpthread on one process and mpipthread on several */
111: /* Logging support */
112: #define    VEC_FILE_CLASSID 1211214
121: PetscPolymorphicSubroutine(VecCreate,(Vec *x),(PETSC_COMM_SELF,x))
123: PetscPolymorphicSubroutine(VecCreateSeq,(PetscInt n,Vec *x),(PETSC_COMM_SELF,n,x))
125: PetscPolymorphicSubroutine(VecCreateMPI,(PetscInt n,PetscInt N,Vec *x),(PETSC_COMM_WORLD,n,N,x))
127: PetscPolymorphicSubroutine(VecCreateSeqWithArray,(PetscInt n,PetscScalar s[],Vec *x),(PETSC_COMM_SELF,n,s,x))
129: PetscPolymorphicSubroutine(VecCreateMPIWithArray,(PetscInt n,PetscInt N,PetscScalar s[],Vec *x),(PETSC_COMM_WORLD,n,N,s,x))
143: PetscPolymorphicFunction(VecDot,(Vec x,Vec y),(x,y,&s),PetscScalar,s)
145: PetscPolymorphicFunction(VecTDot,(Vec x,Vec y),(x,y,&s),PetscScalar,s)
151: /*E
152:     NormType - determines what type of norm to compute
154:     Level: beginner
156: .seealso: VecNorm(), VecNormBegin(), VecNormEnd(), MatNorm()
157: E*/
158: typedef enum {NORM_1=0,NORM_2=1,NORM_FROBENIUS=2,NORM_INFINITY=3,NORM_1_AND_2=4} NormType;
160: #define NORM_MAX NORM_INFINITY
162: /*MC
163:      NORM_1 - the one norm, ||v|| = sum_i | v_i |. ||A|| = max_j || v_*j ||, maximum column sum
165:    Level: beginner
167: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_2, NORM_FROBENIUS, 
168:            NORM_INFINITY, NORM_1_AND_2
170: M*/
172: /*MC
173:      NORM_2 - the two norm, ||v|| = sqrt(sum_i (v_i)^2) (vectors only)
175:    Level: beginner
177: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_FROBENIUS, 
178:            NORM_INFINITY, NORM_1_AND_2
180: M*/
182: /*MC
183:      NORM_FROBENIUS - ||A|| = sqrt(sum_ij (A_ij)^2), same as NORM_2 for vectors
185:    Level: beginner
187: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
188:            NORM_INFINITY, NORM_1_AND_2
190: M*/
192: /*MC
193:      NORM_INFINITY - ||v|| = max_i |v_i|. ||A|| = max_i || v_i* ||, maximum row sum
195:    Level: beginner
197: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
198:            NORM_FROBINIUS, NORM_1_AND_2
200: M*/
202: /*MC
203:      NORM_1_AND_2 - computes both the 1 and 2 norm of a vector
205:    Level: beginner
207: .seealso:  NormType, MatNorm(), VecNorm(), VecNormBegin(), VecNormEnd(), NORM_1, NORM_2, 
208:            NORM_FROBINIUS, NORM_INFINITY
210: M*/
212: /*MC
213:      NORM_MAX - see NORM_INFINITY
215:    Level: beginner
217: M*/
221: PetscPolymorphicSubroutine(VecNorm,(Vec x,PetscReal *r),(x,NORM_2,r))
222: PetscPolymorphicFunction(VecNorm,(Vec x,NormType t),(x,t,&r),PetscReal,r)
223: PetscPolymorphicFunction(VecNorm,(Vec x),(x,NORM_2,&r),PetscReal,r)
227: PetscPolymorphicSubroutine(VecMax,(Vec x,PetscReal *r),(x,PETSC_NULL,r))
229: PetscPolymorphicSubroutine(VecMin,(Vec x,PetscReal *r),(x,PETSC_NULL,r))
242: PetscPolymorphicSubroutine(VecPointwiseMax,(Vec x,Vec y),(x,y,y))
244: PetscPolymorphicSubroutine(VecPointwiseMaxAbs,(Vec x,Vec y),(x,y,y))
246: PetscPolymorphicSubroutine(VecPointwiseMin,(Vec x,Vec y),(x,y,y))
248: PetscPolymorphicSubroutine(VecPointwiseMult,(Vec x,Vec y),(x,x,y))
250: PetscPolymorphicSubroutine(VecPointwiseDivide,(Vec x,Vec y),(x,x,y))
268: PetscPolymorphicFunction(VecStrideNorm,(Vec x,PetscInt i),(x,i,NORM_2,&r),PetscReal,r)
269: PetscPolymorphicFunction(VecStrideNorm,(Vec x,PetscInt i,NormType t),(x,i,t,&r),PetscReal,r)
271: PetscPolymorphicFunction(VecStrideMax,(Vec x,PetscInt i),(x,i,PETSC_NULL,&r),PetscReal,r)
273: PetscPolymorphicFunction(VecStrideMin,(Vec x,PetscInt i),(x,i,PETSC_NULL,&r),PetscReal,r)
291: /*MC
292:    VecSetValue - Set a single entry into a vector.
294:    Synopsis:
295:    PetscErrorCode VecSetValue(Vec v,int row,PetscScalar value, InsertMode mode);
297:    Not Collective
299:    Input Parameters:
300: +  v - the vector
301: .  row - the row location of the entry
302: .  value - the value to insert
303: -  mode - either INSERT_VALUES or ADD_VALUES
305:    Notes:
306:    For efficiency one should use VecSetValues() and set several or 
307:    many values simultaneously if possible.
309:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
310:    MUST be called after all calls to VecSetValues() have been completed.
312:    VecSetValues() uses 0-based indices in Fortran as well as in C.
314:    Level: beginner
316: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(), VecSetValueLocal()
317: M*/
318: PETSC_STATIC_INLINE PetscErrorCode VecSetValue(Vec v,PetscInt i,PetscScalar va,InsertMode mode) {return VecSetValues(v,1,&i,&va,mode);}
323: PetscPolymorphicFunction(VecGetBlockSize,(Vec x),(x,&i),PetscInt,i)
326: /* Dynamic creation and loading functions */
335: /*MC
336:   VecRegisterDynamic - Adds a new vector component implementation
338:   Synopsis:
339:   PetscErrorCode VecRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(Vec))
341:   Not Collective
343:   Input Parameters:
344: + name        - The name of a new user-defined creation routine
345: . path        - The path (either absolute or relative) of the library containing this routine
346: . func_name   - The name of routine to create method context
347: - create_func - The creation routine itself
349:   Notes:
350:   VecRegisterDynamic() may be called multiple times to add several user-defined vectors
352:   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.
354:   Sample usage:
355: .vb
356:     VecRegisterDynamic("my_vec","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyVectorCreate", MyVectorCreate);
357: .ve
359:   Then, your vector type can be chosen with the procedural interface via
360: .vb
361:     VecCreate(MPI_Comm, Vec *);
362:     VecSetType(Vec,"my_vector_name");
363: .ve
364:    or at runtime via the option
365: .vb
366:     -vec_type my_vector_name
367: .ve
369:   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
370:          If your function is not being put into a shared library then use VecRegister() instead
371:         
372:   Level: advanced
374: .keywords: Vec, register
375: .seealso: VecRegisterAll(), VecRegisterDestroy(), VecRegister()
376: M*/
377: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
378: #define VecRegisterDynamic(a,b,c,d) VecRegister(a,b,c,0)
379: #else
380: #define VecRegisterDynamic(a,b,c,d) VecRegister(a,b,c,d)
381: #endif
385: PetscPolymorphicFunction(VecScatterCreate,(Vec x,IS is1,Vec y,IS is2),(x,is1,y,is2,&s),VecScatter,s)
386: PetscPolymorphicSubroutine(VecScatterCreate,(Vec x,IS is,Vec y,VecScatter *s),(x,is,y,PETSC_NULL,s))
387: PetscPolymorphicFunction(VecScatterCreate,(Vec x,IS is,Vec y),(x,is,y,PETSC_NULL,&s),VecScatter,s)
388: PetscPolymorphicSubroutine(VecScatterCreate,(Vec x,Vec y,IS is,VecScatter *s),(x,PETSC_NULL,y,is,s))
389: PetscPolymorphicFunction(VecScatterCreate,(Vec x,Vec y,IS is),(x,PETSC_NULL,y,is,&s),VecScatter,s)
418: PetscPolymorphicFunction(VecEqual,(Vec x,Vec y),(x,y,&s),PetscBool ,s)
422: PetscPolymorphicFunction(VecGetSize,(Vec x),(x,&s),PetscInt,s)
424: PetscPolymorphicFunction(VecGetLocalSize,(Vec x),(x,&s),PetscInt,s)
431: /*MC
432:    VecSetValueLocal - Set a single entry into a vector using the local numbering
434:    Synopsis:
435:    PetscErrorCode VecSetValueLocal(Vec v,int row,PetscScalar value, InsertMode mode);
437:    Not Collective
439:    Input Parameters:
440: +  v - the vector
441: .  row - the row location of the entry
442: .  value - the value to insert
443: -  mode - either INSERT_VALUES or ADD_VALUES
445:    Notes:
446:    For efficiency one should use VecSetValues() and set several or 
447:    many values simultaneously if possible.
449:    These values may be cached, so VecAssemblyBegin() and VecAssemblyEnd() 
450:    MUST be called after all calls to VecSetValues() have been completed.
452:    VecSetValues() uses 0-based indices in Fortran as well as in C.
454:    Level: beginner
456: .seealso: VecSetValues(), VecAssemblyBegin(), VecAssemblyEnd(), VecSetValuesBlockedLocal(), VecSetValue()
457: M*/
458: PETSC_STATIC_INLINE PetscErrorCode VecSetValueLocal(Vec v,PetscInt i,PetscScalar va,InsertMode mode) {return VecSetValuesLocal(v,1,&i,&va,mode);}
466: PetscPolymorphicSubroutine(VecDotBegin,(Vec x,Vec y),(x,y,PETSC_NULL))
468: PetscPolymorphicFunction(VecDotEnd,(Vec x,Vec y),(x,y,&s),PetscScalar,s)
470: PetscPolymorphicSubroutine(VecTDotBegin,(Vec x,Vec y),(x,y,PETSC_NULL))
472: PetscPolymorphicFunction(VecTDotEnd,(Vec x,Vec y),(x,y,&s),PetscScalar,s)
474: PetscPolymorphicSubroutine(VecNormBegin,(Vec x,NormType t),(x,t,PETSC_NULL))
475: PetscPolymorphicSubroutine(VecNormBegin,(Vec x),(x,NORM_2,PETSC_NULL))
477: PetscPolymorphicFunction(VecNormEnd,(Vec x,NormType t),(x,t,&s),PetscReal,s)
478: PetscPolymorphicFunction(VecNormEnd,(Vec x),(x,NORM_2,&s),PetscReal,s)
486: typedef enum {VEC_IGNORE_OFF_PROC_ENTRIES,VEC_IGNORE_NEGATIVE_INDICES} VecOption;
489: /*
490:    Expose VecGetArray()/VecRestoreArray() to users. Allows this to work without any function
491:    call overhead on any 'native' Vecs.
492: */
494:  #include private/vecimpl.h
498: /*
499:     These numbers need to match the entries in 
500:   the function table in vecimpl.h
501: */
502: typedef enum { VECOP_VIEW = 33, VECOP_LOAD = 41, VECOP_DUPLICATE = 0} VecOperation;
505: /*
506:      Routines for dealing with ghosted vectors:
507:   vectors with ghost elements at the end of the array.
508: */
514: PetscPolymorphicFunction(VecGhostGetLocalForm,(Vec x),(x,&s),Vec,s)
527: /*S
528:      Vecs - Collection of vectors where the data for the vectors is stored in 
529:             one contiguous memory
531:    Level: advanced
533:    Notes:
534:     Temporary construct for handling multiply right hand side solves
536:     This is faked by storing a single vector that has enough array space for 
537:     n vectors
539:   Concepts: parallel decomposition
541: S*/
542:         struct _n_Vecs  {PetscInt n; Vec v;};
543: typedef struct _n_Vecs* Vecs;
544: #define VecsDestroy(x)            (VecDestroy(&(x)->v)         || PetscFree(x))
545: #define VecsCreateSeq(comm,p,m,x) (PetscNew(struct _n_Vecs,x) || VecCreateSeq(comm,p*m,&(*(x))->v) || (-1 == ((*(x))->n = (m))))
546: #define VecsCreateSeqWithArray(comm,p,m,a,x) (PetscNew(struct _n_Vecs,x) || VecCreateSeqWithArray(comm,p*m,a,&(*(x))->v) || (-1 == ((*(x))->n = (m))))
547: #define VecsDuplicate(x,y)        (PetscNew(struct _n_Vecs,y) || VecDuplicate(x->v,&(*(y))->v) || (-1 == ((*(y))->n = (x)->n)))
549: #if defined(PETSC_HAVE_CUSP)
550: typedef struct _p_PetscCUSPIndices* PetscCUSPIndices;
558: #endif
566: #endif