Actual source code: vecimpl.h

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:    This private file should not be included in users' code.
  4:    Defines the fields shared by all vector implementations.

  6: */

  8: #ifndef __VECIMPL_H

 11:  #include <petscvec.h>
 12:  #include <petsc/private/petscimpl.h>
 13:  #include <petscviewer.h>

 15: PETSC_EXTERN PetscBool VecRegisterAllCalled;
 16: PETSC_EXTERN PetscErrorCode VecRegisterAll(void);

 18: /* ----------------------------------------------------------------------------*/

 20: typedef struct _VecOps *VecOps;
 21: struct _VecOps {
 22:   PetscErrorCode (*duplicate)(Vec,Vec*);         /* get single vector */
 23:   PetscErrorCode (*duplicatevecs)(Vec,PetscInt,Vec**);     /* get array of vectors */
 24:   PetscErrorCode (*destroyvecs)(PetscInt,Vec[]);           /* free array of vectors */
 25:   PetscErrorCode (*dot)(Vec,Vec,PetscScalar*);             /* z = x^H * y */
 26:   PetscErrorCode (*mdot)(Vec,PetscInt,const Vec[],PetscScalar*); /* z[j] = x dot y[j] */
 27:   PetscErrorCode (*norm)(Vec,NormType,PetscReal*);        /* z = sqrt(x^H * x) */
 28:   PetscErrorCode (*tdot)(Vec,Vec,PetscScalar*);             /* x'*y */
 29:   PetscErrorCode (*mtdot)(Vec,PetscInt,const Vec[],PetscScalar*);/* z[j] = x dot y[j] */
 30:   PetscErrorCode (*scale)(Vec,PetscScalar);                 /* x = alpha * x   */
 31:   PetscErrorCode (*copy)(Vec,Vec);                     /* y = x */
 32:   PetscErrorCode (*set)(Vec,PetscScalar);                        /* y = alpha  */
 33:   PetscErrorCode (*swap)(Vec,Vec);                               /* exchange x and y */
 34:   PetscErrorCode (*axpy)(Vec,PetscScalar,Vec);                   /* y = y + alpha * x */
 35:   PetscErrorCode (*axpby)(Vec,PetscScalar,PetscScalar,Vec);      /* y = alpha * x + beta * y*/
 36:   PetscErrorCode (*maxpy)(Vec,PetscInt,const PetscScalar*,Vec*); /* y = y + alpha[j] x[j] */
 37:   PetscErrorCode (*aypx)(Vec,PetscScalar,Vec);                   /* y = x + alpha * y */
 38:   PetscErrorCode (*waxpy)(Vec,PetscScalar,Vec,Vec);         /* w = y + alpha * x */
 39:   PetscErrorCode (*axpbypcz)(Vec,PetscScalar,PetscScalar,PetscScalar,Vec,Vec);   /* z = alpha * x + beta *y + gamma *z*/
 40:   PetscErrorCode (*pointwisemult)(Vec,Vec,Vec);        /* w = x .* y */
 41:   PetscErrorCode (*pointwisedivide)(Vec,Vec,Vec);      /* w = x ./ y */
 42:   PetscErrorCode (*setvalues)(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 43:   PetscErrorCode (*assemblybegin)(Vec);                /* start global assembly */
 44:   PetscErrorCode (*assemblyend)(Vec);                  /* end global assembly */
 45:   PetscErrorCode (*getarray)(Vec,PetscScalar**);            /* get data array */
 46:   PetscErrorCode (*getsize)(Vec,PetscInt*);
 47:   PetscErrorCode (*getlocalsize)(Vec,PetscInt*);
 48:   PetscErrorCode (*restorearray)(Vec,PetscScalar**);        /* restore data array */
 49:   PetscErrorCode (*max)(Vec,PetscInt*,PetscReal*);      /* z = max(x); idx=index of max(x) */
 50:   PetscErrorCode (*min)(Vec,PetscInt*,PetscReal*);      /* z = min(x); idx=index of min(x) */
 51:   PetscErrorCode (*setrandom)(Vec,PetscRandom);         /* set y[j] = random numbers */
 52:   PetscErrorCode (*setoption)(Vec,VecOption,PetscBool );
 53:   PetscErrorCode (*setvaluesblocked)(Vec,PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 54:   PetscErrorCode (*destroy)(Vec);
 55:   PetscErrorCode (*view)(Vec,PetscViewer);
 56:   PetscErrorCode (*placearray)(Vec,const PetscScalar*);     /* place data array */
 57:   PetscErrorCode (*replacearray)(Vec,const PetscScalar*);     /* replace data array */
 58:   PetscErrorCode (*dot_local)(Vec,Vec,PetscScalar*);
 59:   PetscErrorCode (*tdot_local)(Vec,Vec,PetscScalar*);
 60:   PetscErrorCode (*norm_local)(Vec,NormType,PetscReal*);
 61:   PetscErrorCode (*mdot_local)(Vec,PetscInt,const Vec[],PetscScalar*);
 62:   PetscErrorCode (*mtdot_local)(Vec,PetscInt,const Vec[],PetscScalar*);
 63:   PetscErrorCode (*load)(Vec,PetscViewer);
 64:   PetscErrorCode (*reciprocal)(Vec);
 65:   PetscErrorCode (*conjugate)(Vec);
 66:   PetscErrorCode (*setlocaltoglobalmapping)(Vec,ISLocalToGlobalMapping);
 67:   PetscErrorCode (*setvalueslocal)(Vec,PetscInt,const PetscInt *,const PetscScalar *,InsertMode);
 68:   PetscErrorCode (*resetarray)(Vec);      /* vector points to its original array, i.e. undoes any VecPlaceArray() */
 69:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,Vec);
 70:   PetscErrorCode (*maxpointwisedivide)(Vec,Vec,PetscReal*);      /* m = max abs(x ./ y) */
 71:   PetscErrorCode (*pointwisemax)(Vec,Vec,Vec);
 72:   PetscErrorCode (*pointwisemaxabs)(Vec,Vec,Vec);
 73:   PetscErrorCode (*pointwisemin)(Vec,Vec,Vec);
 74:   PetscErrorCode (*getvalues)(Vec,PetscInt,const PetscInt[],PetscScalar[]);
 75:   PetscErrorCode (*sqrt)(Vec);
 76:   PetscErrorCode (*abs)(Vec);
 77:   PetscErrorCode (*exp)(Vec);
 78:   PetscErrorCode (*log)(Vec);
 79:   PetscErrorCode (*shift)(Vec,PetscScalar);
 80:   PetscErrorCode (*create)(Vec);
 81:   PetscErrorCode (*stridegather)(Vec,PetscInt,Vec,InsertMode);
 82:   PetscErrorCode (*stridescatter)(Vec,PetscInt,Vec,InsertMode);
 83:   PetscErrorCode (*dotnorm2)(Vec,Vec,PetscScalar*,PetscScalar*);
 84:   PetscErrorCode (*getsubvector)(Vec,IS,Vec*);
 85:   PetscErrorCode (*restoresubvector)(Vec,IS,Vec*);
 86:   PetscErrorCode (*getarrayread)(Vec,const PetscScalar**);
 87:   PetscErrorCode (*restorearrayread)(Vec,const PetscScalar**);
 88:   PetscErrorCode (*stridesubsetgather)(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);
 89:   PetscErrorCode (*stridesubsetscatter)(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);
 90:   PetscErrorCode (*viewnative)(Vec,PetscViewer);
 91:   PetscErrorCode (*loadnative)(Vec,PetscViewer);
 92:   PetscErrorCode (*getlocalvector)(Vec,Vec);
 93:   PetscErrorCode (*restorelocalvector)(Vec,Vec);
 94:   PetscErrorCode (*getlocalvectorread)(Vec,Vec);
 95:   PetscErrorCode (*restorelocalvectorread)(Vec,Vec);
 96:   PetscErrorCode (*bindtocpu)(Vec,PetscBool);
 97:   PetscErrorCode (*getarraywrite)(Vec,PetscScalar**);
 98:   PetscErrorCode (*restorearraywrite)(Vec,PetscScalar**);
 99: };

101: /*
102:     The stash is used to temporarily store inserted vec values that
103:   belong to another processor. During the assembly phase the stashed
104:   values are moved to the correct processor and
105: */

107: typedef struct {
108:   PetscInt      nmax;                   /* maximum stash size */
109:   PetscInt      umax;                   /* max stash size user wants */
110:   PetscInt      oldnmax;                /* the nmax value used previously */
111:   PetscInt      n;                      /* stash size */
112:   PetscInt      bs;                     /* block size of the stash */
113:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
114:   PetscInt      *idx;                   /* global row numbers in stash */
115:   PetscScalar   *array;                 /* array to hold stashed values */
116:   /* The following variables are used for communication */
117:   MPI_Comm      comm;
118:   PetscMPIInt   size,rank;
119:   PetscMPIInt   tag1,tag2;
120:   MPI_Request   *send_waits;            /* array of send requests */
121:   MPI_Request   *recv_waits;            /* array of receive requests */
122:   MPI_Status    *send_status;           /* array of send status */
123:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
124:   PetscScalar   *svalues,*rvalues;      /* sending and receiving data */
125:   PetscInt      *sindices,*rindices;
126:   PetscInt      rmax;                   /* maximum message length */
127:   PetscInt      *nprocs;                /* tmp data used both during scatterbegin and end */
128:   PetscInt      nprocessed;             /* number of messages already processed */
129:   PetscBool     donotstash;
130:   PetscBool     ignorenegidx;           /* ignore negative indices passed into VecSetValues/VetGetValues */
131:   InsertMode    insertmode;
132:   PetscInt      *bowners;
133: } VecStash;

135: struct _p_Vec {
136:   PETSCHEADER(struct _VecOps);
137:   PetscLayout            map;
138:   void                   *data;     /* implementation-specific data */
139:   PetscBool              array_gotten;
140:   VecStash               stash,bstash; /* used for storing off-proc values during assembly */
141:   PetscBool              petscnative;  /* means the ->data starts with VECHEADER and can use VecGetArrayFast()*/
142:   PetscInt               lock;         /* lock state. vector can be free (=0), locked for read (>0) or locked for write(<0) */
143: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
144:   void                   *spptr; /* this is the special pointer to the array on the GPU */
145:   PetscOffloadMask       offloadmask;  /* a mask which indicates where the valid vector data is (GPU, CPU or both) */
146:   PetscBool              boundtocpu;
147:   size_t                 minimum_bytes_pinned_memory; /* minimum data size in bytes for which pinned memory will be allocated */
148:   PetscBool              pinned_memory; /* PETSC_TRUE if the current host allocation has been made from pinned memory. */
149: #endif
150: };

152: PETSC_EXTERN PetscLogEvent VEC_SetRandom;
153: PETSC_EXTERN PetscLogEvent VEC_View;
154: PETSC_EXTERN PetscLogEvent VEC_Max;
155: PETSC_EXTERN PetscLogEvent VEC_Min;
156: PETSC_EXTERN PetscLogEvent VEC_Dot;
157: PETSC_EXTERN PetscLogEvent VEC_MDot;
158: PETSC_EXTERN PetscLogEvent VEC_TDot;
159: PETSC_EXTERN PetscLogEvent VEC_MTDot;
160: PETSC_EXTERN PetscLogEvent VEC_Norm;
161: PETSC_EXTERN PetscLogEvent VEC_Normalize;
162: PETSC_EXTERN PetscLogEvent VEC_Scale;
163: PETSC_EXTERN PetscLogEvent VEC_Copy;
164: PETSC_EXTERN PetscLogEvent VEC_Set;
165: PETSC_EXTERN PetscLogEvent VEC_AXPY;
166: PETSC_EXTERN PetscLogEvent VEC_AYPX;
167: PETSC_EXTERN PetscLogEvent VEC_WAXPY;
168: PETSC_EXTERN PetscLogEvent VEC_MAXPY;
169: PETSC_EXTERN PetscLogEvent VEC_AssemblyEnd;
170: PETSC_EXTERN PetscLogEvent VEC_PointwiseMult;
171: PETSC_EXTERN PetscLogEvent VEC_SetValues;
172: PETSC_EXTERN PetscLogEvent VEC_Load;
173: PETSC_EXTERN PetscLogEvent VEC_ScatterBegin;
174: PETSC_EXTERN PetscLogEvent VEC_ScatterEnd;
175: PETSC_EXTERN PetscLogEvent VEC_ReduceArithmetic;
176: PETSC_EXTERN PetscLogEvent VEC_ReduceCommunication;
177: PETSC_EXTERN PetscLogEvent VEC_ReduceBegin;
178: PETSC_EXTERN PetscLogEvent VEC_ReduceEnd;
179: PETSC_EXTERN PetscLogEvent VEC_Swap;
180: PETSC_EXTERN PetscLogEvent VEC_AssemblyBegin;
181: PETSC_EXTERN PetscLogEvent VEC_DotNorm2;
182: PETSC_EXTERN PetscLogEvent VEC_AXPBYPCZ;
183: PETSC_EXTERN PetscLogEvent VEC_Ops;
184: PETSC_EXTERN PetscLogEvent VEC_ViennaCLCopyToGPU;
185: PETSC_EXTERN PetscLogEvent VEC_ViennaCLCopyFromGPU;
186: PETSC_EXTERN PetscLogEvent VEC_CUDACopyToGPU;
187: PETSC_EXTERN PetscLogEvent VEC_CUDACopyFromGPU;
188: PETSC_EXTERN PetscLogEvent VEC_CUDACopyToGPUSome;
189: PETSC_EXTERN PetscLogEvent VEC_CUDACopyFromGPUSome;

191: PETSC_EXTERN PetscErrorCode VecView_Seq(Vec,PetscViewer);
192: #if defined(PETSC_HAVE_VIENNACL)
193: PETSC_EXTERN PetscErrorCode VecViennaCLAllocateCheckHost(Vec v);
194: PETSC_EXTERN PetscErrorCode VecViennaCLCopyFromGPU(Vec v);
195: #endif
196: #if defined(PETSC_HAVE_CUDA)
197: PETSC_EXTERN PetscErrorCode VecCUDAAllocateCheckHost(Vec v);
198: PETSC_EXTERN PetscErrorCode VecCUDACopyFromGPU(Vec v);
199: #endif


202: /*
203:      Common header shared by array based vectors,
204:    currently Vec_Seq and Vec_MPI
205: */
206: #define VECHEADER                          \
207:   PetscScalar *array;                      \
208:   PetscScalar *array_allocated;                        /* if the array was allocated by PETSc this is its pointer */  \
209:   PetscScalar *unplacedarray;                           /* if one called VecPlaceArray(), this is where it stashed the original */

211: /* Lock a vector for exclusive read&write access */
212: #if defined(PETSC_USE_DEBUG)
213: PETSC_INTERN PetscErrorCode VecLockWriteSet_Private(Vec,PetscBool);
214: #else
215: #define VecLockWriteSet_Private(x,flg) 0
216: #endif

218: /* Default obtain and release vectors; can be used by any implementation */
219: PETSC_EXTERN PetscErrorCode VecDuplicateVecs_Default(Vec,PetscInt,Vec *[]);
220: PETSC_EXTERN PetscErrorCode VecDestroyVecs_Default(PetscInt,Vec []);
221: PETSC_EXTERN PetscErrorCode VecView_Binary(Vec, PetscViewer);
222: PETSC_EXTERN PetscErrorCode VecLoad_Binary(Vec, PetscViewer);
223: PETSC_EXTERN PetscErrorCode VecLoad_Default(Vec, PetscViewer);

225: PETSC_EXTERN PetscInt  NormIds[7];  /* map from NormType to IDs used to cache/retreive values of norms */

227: PETSC_INTERN PetscErrorCode VecStashCreate_Private(MPI_Comm,PetscInt,VecStash*);
228: PETSC_INTERN PetscErrorCode VecStashDestroy_Private(VecStash*);
229: PETSC_EXTERN PetscErrorCode VecStashExpand_Private(VecStash*,PetscInt);
230: PETSC_INTERN PetscErrorCode VecStashScatterEnd_Private(VecStash*);
231: PETSC_INTERN PetscErrorCode VecStashSetInitialSize_Private(VecStash*,PetscInt);
232: PETSC_INTERN PetscErrorCode VecStashGetInfo_Private(VecStash*,PetscInt*,PetscInt*);
233: PETSC_INTERN PetscErrorCode VecStashScatterBegin_Private(VecStash*,PetscInt*);
234: PETSC_INTERN PetscErrorCode VecStashScatterGetMesg_Private(VecStash*,PetscMPIInt*,PetscInt**,PetscScalar**,PetscInt*);
235: PETSC_INTERN PetscErrorCode VecStashSortCompress_Private(VecStash*);
236: PETSC_INTERN PetscErrorCode VecStashGetOwnerList_Private(VecStash*,PetscLayout,PetscMPIInt*,PetscMPIInt**);

238: /*
239:   VecStashValue_Private - inserts a single value into the stash.

241:   Input Parameters:
242:   stash  - the stash
243:   idx    - the global of the inserted value
244:   values - the value inserted
245: */
246: PETSC_STATIC_INLINE PetscErrorCode VecStashValue_Private(VecStash *stash,PetscInt row,PetscScalar value)
247: {
249:   /* Check and see if we have sufficient memory */
250:   if (((stash)->n + 1) > (stash)->nmax) {
251:     VecStashExpand_Private(stash,1);
252:   }
253:   (stash)->idx[(stash)->n]   = row;
254:   (stash)->array[(stash)->n] = value;
255:   (stash)->n++;
256:   return 0;
257: }

259: /*
260:   VecStashValuesBlocked_Private - inserts 1 block of values into the stash.

262:   Input Parameters:
263:   stash  - the stash
264:   idx    - the global block index
265:   values - the values inserted
266: */
267: PETSC_STATIC_INLINE PetscErrorCode VecStashValuesBlocked_Private(VecStash *stash,PetscInt row,PetscScalar *values)
268: {
269:   PetscInt       jj,stash_bs=(stash)->bs;
270:   PetscScalar    *array;
272:   if (((stash)->n+1) > (stash)->nmax) {
273:     VecStashExpand_Private(stash,1);
274:   }
275:   array = (stash)->array + stash_bs*(stash)->n;
276:   (stash)->idx[(stash)->n]   = row;
277:   for (jj=0; jj<stash_bs; jj++) { array[jj] = values[jj];}
278:   (stash)->n++;
279:   return 0;
280: }

282: PETSC_INTERN PetscErrorCode VecStrideGather_Default(Vec,PetscInt,Vec,InsertMode);
283: PETSC_INTERN PetscErrorCode VecStrideScatter_Default(Vec,PetscInt,Vec,InsertMode);
284: PETSC_INTERN PetscErrorCode VecReciprocal_Default(Vec);
285: PETSC_INTERN PetscErrorCode VecStrideSubSetGather_Default(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);
286: PETSC_INTERN PetscErrorCode VecStrideSubSetScatter_Default(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);

288: #if defined(PETSC_HAVE_MATLAB_ENGINE)
289: PETSC_EXTERN PetscErrorCode VecMatlabEnginePut_Default(PetscObject,void*);
290: PETSC_EXTERN PetscErrorCode VecMatlabEngineGet_Default(PetscObject,void*);
291: #endif

293: PETSC_EXTERN PetscErrorCode PetscSectionGetField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);
294: PETSC_EXTERN PetscErrorCode PetscSectionRestoreField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);

296: #define VecCheckSameLocalSize(x,ar1,y,ar2) do { \
297:     if ((x)->map->n != (y)->map->n) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths parameter # %d local size %D != parameter # %d local size %D", ar1,(x)->map->n, ar2,(y)->map->n); \
298:   } while (0)

300: #define VecCheckSameSize(x,ar1,y,ar2) do { \
301:     if ((x)->map->N != (y)->map->N) SETERRQ4(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths parameter # %d global size %D != parameter # %d global size %D", ar1,(x)->map->N, ar2,(y)->map->N); \
302:     VecCheckSameLocalSize(x,ar1,y,ar2); \
303:   } while(0)

305: #define VecCheckLocalSize(x,ar1,n) do { \
306:     if ((x)->map->n != n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect vector local size: parameter # %d local size %D != %D",ar1,(x)->map->n,n); \
307:   } while (0)

309: #define VecCheckSize(x,ar1,n,N) do { \
310:     if ((x)->map->N != N) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_INCOMP,"Incorrect vector global size: parameter # %d global size %D != %D",ar1,(x)->map->N, N); \
311:     VecCheckLocalSize(x,ar1,n); \
312:   } while (0)

314: typedef struct _VecTaggerOps *VecTaggerOps;
315: struct _VecTaggerOps {
316:   PetscErrorCode (*create) (VecTagger);
317:   PetscErrorCode (*destroy) (VecTagger);
318:   PetscErrorCode (*setfromoptions) (PetscOptionItems*,VecTagger);
319:   PetscErrorCode (*setup) (VecTagger);
320:   PetscErrorCode (*view) (VecTagger,PetscViewer);
321:   PetscErrorCode (*computeboxes) (VecTagger,Vec,PetscInt *,VecTaggerBox **);
322:   PetscErrorCode (*computeis) (VecTagger,Vec,IS *);
323: };
324: struct _p_VecTagger {
325:   PETSCHEADER(struct _VecTaggerOps);
326:   void      *data;
327:   PetscInt  blocksize;
328:   PetscBool invert;
329:   PetscBool setupcalled;
330: };

332: PETSC_EXTERN PetscBool      VecTaggerRegisterAllCalled;
333: PETSC_EXTERN PetscErrorCode VecTaggerRegisterAll(void);
334: PETSC_EXTERN PetscErrorCode VecTaggerComputeIS_FromBoxes(VecTagger,Vec,IS*);
335: PETSC_EXTERN PetscMPIInt Petsc_Reduction_keyval;

337: #endif