Actual source code: vecimpl.h

  1: #pragma once

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

  8: #include <petscvec.h>
  9: #include <petsc/private/petscimpl.h>

 11: PETSC_INTERN PetscBool VecRegisterAllCalled;
 12: PETSC_EXTERN MPI_Op    MPIU_MAXLOC;
 13: PETSC_EXTERN MPI_Op    MPIU_MINLOC;

 15: /* ----------------------------------------------------------------------------*/

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

112: #if defined(offsetof) && (defined(__cplusplus) || (PETSC_C_VERSION >= 23))
113: static_assert(offsetof(struct _VecOps, duplicate) == sizeof(void (*)(void)) * VECOP_DUPLICATE, "");
114: static_assert(offsetof(struct _VecOps, set) == sizeof(void (*)(void)) * VECOP_SET, "");
115: static_assert(offsetof(struct _VecOps, view) == sizeof(void (*)(void)) * VECOP_VIEW, "");
116: static_assert(offsetof(struct _VecOps, load) == sizeof(void (*)(void)) * VECOP_LOAD, "");
117: static_assert(offsetof(struct _VecOps, viewnative) == sizeof(void (*)(void)) * VECOP_VIEWNATIVE, "");
118: static_assert(offsetof(struct _VecOps, loadnative) == sizeof(void (*)(void)) * VECOP_LOADNATIVE, "");
119: #endif

121: /*
122:     The stash is used to temporarily store inserted vec values that
123:   belong to another processor. During the assembly phase the stashed
124:   values are moved to the correct processor and
125: */

127: typedef struct {
128:   PetscInt     nmax;     /* maximum stash size */
129:   PetscInt     umax;     /* max stash size user wants */
130:   PetscInt     oldnmax;  /* the nmax value used previously */
131:   PetscInt     n;        /* stash size */
132:   PetscInt     bs;       /* block size of the stash */
133:   PetscInt     reallocs; /* preserve the no of mallocs invoked */
134:   PetscInt    *idx;      /* global row numbers in stash */
135:   PetscScalar *array;    /* array to hold stashed values */
136:   /* The following variables are used for communication */
137:   MPI_Comm     comm;
138:   PetscMPIInt  size, rank;
139:   PetscMPIInt  tag1, tag2;
140:   MPI_Request *send_waits;        /* array of send requests */
141:   MPI_Request *recv_waits;        /* array of receive requests */
142:   MPI_Status  *send_status;       /* array of send status */
143:   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
144:   PetscScalar *svalues, *rvalues; /* sending and receiving data */
145:   PetscInt    *sindices, *rindices;
146:   PetscInt     rmax;       /* maximum message length */
147:   PetscInt    *nprocs;     /* tmp data used both during scatterbegin and end */
148:   PetscInt     nprocessed; /* number of messages already processed */
149:   PetscBool    donotstash;
150:   PetscBool    ignorenegidx; /* ignore negative indices passed into VecSetValues/VetGetValues */
151:   InsertMode   insertmode;
152:   PetscInt    *bowners;
153: } VecStash;

155: struct _p_Vec {
156:   PETSCHEADER(struct _VecOps);
157:   PetscLayout map;
158:   void       *data; /* implementation-specific data */
159:   PetscBool   array_gotten;
160:   VecStash    stash, bstash; /* used for storing off-proc values during assembly */
161:   PetscBool   petscnative;   /* means the ->data starts with VECHEADER and can use VecGetArrayFast()*/
162: #if PetscDefined(USE_DEBUG)
163:   PetscStack lockstack; /* the file,func,line of where locks are added */
164:   PetscInt   lock;      /* lock state. vector can be free (=0), locked for read (>0) or locked for write(<0) */
165: #endif
166:   PetscOffloadMask offloadmask; /* a mask which indicates where the valid vector data is (GPU, CPU or both) */
167: #if defined(PETSC_HAVE_DEVICE)
168:   void     *spptr; /* this is the special pointer to the array on the GPU */
169:   PetscBool boundtocpu;
170:   PetscBool bindingpropagates;
171:   size_t    minimum_bytes_pinned_memory; /* minimum data size in bytes for which pinned memory will be allocated */
172:   PetscBool pinned_memory;               /* PETSC_TRUE if the current host allocation has been made from pinned memory. */
173: #endif
174:   char *defaultrandtype;
175: };

177: PETSC_EXTERN PetscLogEvent VEC_SetRandom;
178: PETSC_EXTERN PetscLogEvent VEC_View;
179: PETSC_EXTERN PetscLogEvent VEC_Max;
180: PETSC_EXTERN PetscLogEvent VEC_Min;
181: PETSC_EXTERN PetscLogEvent VEC_Dot;
182: PETSC_EXTERN PetscLogEvent VEC_MDot;
183: PETSC_EXTERN PetscLogEvent VEC_TDot;
184: PETSC_EXTERN PetscLogEvent VEC_MTDot;
185: PETSC_EXTERN PetscLogEvent VEC_Norm;
186: PETSC_EXTERN PetscLogEvent VEC_Normalize;
187: PETSC_EXTERN PetscLogEvent VEC_Scale;
188: PETSC_EXTERN PetscLogEvent VEC_Shift;
189: PETSC_EXTERN PetscLogEvent VEC_Copy;
190: PETSC_EXTERN PetscLogEvent VEC_Set;
191: PETSC_EXTERN PetscLogEvent VEC_AXPY;
192: PETSC_EXTERN PetscLogEvent VEC_AYPX;
193: PETSC_EXTERN PetscLogEvent VEC_WAXPY;
194: PETSC_EXTERN PetscLogEvent VEC_MAXPY;
195: PETSC_EXTERN PetscLogEvent VEC_AssemblyEnd;
196: PETSC_EXTERN PetscLogEvent VEC_PointwiseMult;
197: PETSC_EXTERN PetscLogEvent VEC_SetValues;
198: PETSC_EXTERN PetscLogEvent VEC_SetPreallocateCOO;
199: PETSC_EXTERN PetscLogEvent VEC_SetValuesCOO;
200: PETSC_EXTERN PetscLogEvent VEC_Load;
201: PETSC_EXTERN PetscLogEvent VEC_ScatterBegin;
202: PETSC_EXTERN PetscLogEvent VEC_ScatterEnd;
203: PETSC_EXTERN PetscLogEvent VEC_ReduceArithmetic;
204: PETSC_EXTERN PetscLogEvent VEC_ReduceCommunication;
205: PETSC_EXTERN PetscLogEvent VEC_ReduceBegin;
206: PETSC_EXTERN PetscLogEvent VEC_ReduceEnd;
207: PETSC_EXTERN PetscLogEvent VEC_Swap;
208: PETSC_EXTERN PetscLogEvent VEC_AssemblyBegin;
209: PETSC_EXTERN PetscLogEvent VEC_DotNorm2;
210: PETSC_EXTERN PetscLogEvent VEC_AXPBYPCZ;
211: PETSC_EXTERN PetscLogEvent VEC_Ops;
212: PETSC_EXTERN PetscLogEvent VEC_ViennaCLCopyToGPU;
213: PETSC_EXTERN PetscLogEvent VEC_ViennaCLCopyFromGPU;
214: PETSC_EXTERN PetscLogEvent VEC_CUDACopyToGPU;
215: PETSC_EXTERN PetscLogEvent VEC_CUDACopyFromGPU;
216: PETSC_EXTERN PetscLogEvent VEC_HIPCopyToGPU;
217: PETSC_EXTERN PetscLogEvent VEC_HIPCopyFromGPU;

219: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecView_Seq(Vec, PetscViewer);
220: #if defined(PETSC_HAVE_VIENNACL)
221: PETSC_EXTERN PetscErrorCode VecViennaCLAllocateCheckHost(Vec v);
222: PETSC_EXTERN PetscErrorCode VecViennaCLCopyFromGPU(Vec v);
223: #endif

225: /*
226:      Common header shared by array based vectors,
227:    currently Vec_Seq and Vec_MPI
228: */
229: #define VECHEADER \
230:   PetscScalar *array; \
231:   PetscScalar *array_allocated; /* if the array was allocated by PETSc this is its pointer */ \
232:   PetscScalar *unplacedarray;   /* if one called VecPlaceArray(), this is where it stashed the original */

234: /* Get Root type of vector. e.g. VECSEQ -> VECSTANDARD, VECMPICUDA -> VECCUDA */
235: PETSC_EXTERN PetscErrorCode VecGetRootType_Private(Vec, VecType *);

237: /* Default obtain and release vectors; can be used by any implementation */
238: PETSC_INTERN PetscErrorCode VecDuplicateVecs_Default(Vec, PetscInt, Vec *[]);
239: PETSC_INTERN PetscErrorCode VecDestroyVecs_Default(PetscInt, Vec[]);
240: PETSC_INTERN PetscErrorCode VecView_Binary(Vec, PetscViewer);

242: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecLoad_Default(Vec, PetscViewer);

244: PETSC_INTERN PetscInt NormIds[4]; /* map from NormType to IDs used to cache/retrieve values of norms, 1_AND_2 is excluded */

246: PETSC_INTERN PetscErrorCode VecStashCreate_Private(MPI_Comm, PetscInt, VecStash *);
247: PETSC_INTERN PetscErrorCode VecStashDestroy_Private(VecStash *);
248: PETSC_INTERN PetscErrorCode VecStashScatterEnd_Private(VecStash *);
249: PETSC_INTERN PetscErrorCode VecStashSetInitialSize_Private(VecStash *, PetscInt);
250: PETSC_INTERN PetscErrorCode VecStashGetInfo_Private(VecStash *, PetscInt *, PetscInt *);
251: PETSC_INTERN PetscErrorCode VecStashScatterBegin_Private(VecStash *, const PetscInt *);
252: PETSC_INTERN PetscErrorCode VecStashScatterGetMesg_Private(VecStash *, PetscMPIInt *, PetscInt **, PetscScalar **, PetscInt *);
253: PETSC_INTERN PetscErrorCode VecStashSortCompress_Private(VecStash *);
254: PETSC_INTERN PetscErrorCode VecStashGetOwnerList_Private(VecStash *, PetscLayout, PetscMPIInt *, PetscMPIInt **);

256: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecStashExpand_Private(VecStash *, PetscInt);

258: /*
259:   VecStashValue_Private - inserts a single value into the stash.

261:   Input Parameters:
262:   stash  - the stash
263:   idx    - the global of the inserted value
264:   values - the value inserted
265: */
266: static inline PetscErrorCode VecStashValue_Private(VecStash *stash, PetscInt row, PetscScalar value)
267: {
268:   /* Check and see if we have sufficient memory */
269:   PetscFunctionBegin;
270:   if (((stash)->n + 1) > (stash)->nmax) PetscCall(VecStashExpand_Private(stash, 1));
271:   (stash)->idx[(stash)->n]   = row;
272:   (stash)->array[(stash)->n] = value;
273:   (stash)->n++;
274:   PetscFunctionReturn(PETSC_SUCCESS);
275: }

277: /*
278:   VecStashValuesBlocked_Private - inserts 1 block of values into the stash.

280:   Input Parameters:
281:   stash  - the stash
282:   idx    - the global block index
283:   values - the values inserted
284: */
285: static inline PetscErrorCode VecStashValuesBlocked_Private(VecStash *stash, PetscInt row, PetscScalar *values)
286: {
287:   PetscInt     stash_bs = (stash)->bs;
288:   PetscScalar *array;

290:   PetscFunctionBegin;
291:   if (((stash)->n + 1) > (stash)->nmax) PetscCall(VecStashExpand_Private(stash, 1));
292:   array                    = (stash)->array + stash_bs * (stash)->n;
293:   (stash)->idx[(stash)->n] = row;
294:   PetscCall(PetscArraycpy(array, values, stash_bs));
295:   (stash)->n++;
296:   PetscFunctionReturn(PETSC_SUCCESS);
297: }

299: PETSC_INTERN PetscErrorCode VecStrideGather_Default(Vec, PetscInt, Vec, InsertMode);
300: PETSC_INTERN PetscErrorCode VecStrideScatter_Default(Vec, PetscInt, Vec, InsertMode);
301: PETSC_INTERN PetscErrorCode VecStrideSubSetGather_Default(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);
302: PETSC_INTERN PetscErrorCode VecStrideSubSetScatter_Default(Vec, PetscInt, const PetscInt[], const PetscInt[], Vec, InsertMode);

304: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecReciprocal_Default(Vec);
305: #if defined(PETSC_HAVE_MATLAB)
306: PETSC_EXTERN PetscErrorCode VecMatlabEnginePut_Default(PetscObject, void *);
307: PETSC_EXTERN PetscErrorCode VecMatlabEngineGet_Default(PetscObject, void *);
308: #endif

310: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscSectionGetField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);
311: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscSectionRestoreField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);

313: #define VecCheckSameLocalSize(x, ar1, y, ar2) \
314:   do { \
315:     PetscCheck((x)->map->n == (y)->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible vector local lengths parameter # %d local size %" PetscInt_FMT " != parameter # %d local size %" PetscInt_FMT, ar1, (x)->map->n, ar2, (y)->map->n); \
316:   } while (0)

318: #define VecCheckSameSize(x, ar1, y, ar2) \
319:   do { \
320:     PetscCheck((x)->map->N == (y)->map->N, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_INCOMP, "Incompatible vector global lengths parameter # %d global size %" PetscInt_FMT " != parameter # %d global size %" PetscInt_FMT, ar1, (x)->map->N, ar2, \
321:                (y)->map->N); \
322:     VecCheckSameLocalSize(x, ar1, y, ar2); \
323:   } while (0)

325: #define VecCheckLocalSize(x, ar1, n) \
326:   do { \
327:     PetscCheck((x)->map->n == (n), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incorrect vector local size: parameter # %d local size %" PetscInt_FMT " != %" PetscInt_FMT, ar1, (x)->map->n, n); \
328:   } while (0)

330: #define VecCheckSize(x, ar1, n, N) \
331:   do { \
332:     PetscCheck((x)->map->N == (N), PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_INCOMP, "Incorrect vector global size: parameter # %d global size %" PetscInt_FMT " != %" PetscInt_FMT, ar1, (x)->map->N, N); \
333:     VecCheckLocalSize(x, ar1, n); \
334:   } while (0)

336: typedef struct _VecTaggerOps *VecTaggerOps;
337: struct _VecTaggerOps {
338:   PetscErrorCode (*create)(VecTagger);
339:   PetscErrorCode (*destroy)(VecTagger);
340:   PetscErrorCode (*setfromoptions)(VecTagger, PetscOptionItems *);
341:   PetscErrorCode (*setup)(VecTagger);
342:   PetscErrorCode (*view)(VecTagger, PetscViewer);
343:   PetscErrorCode (*computeboxes)(VecTagger, Vec, PetscInt *, VecTaggerBox **, PetscBool *);
344:   PetscErrorCode (*computeis)(VecTagger, Vec, IS *, PetscBool *);
345: };
346: struct _p_VecTagger {
347:   PETSCHEADER(struct _VecTaggerOps);
348:   void     *data;
349:   PetscInt  blocksize;
350:   PetscBool invert;
351:   PetscBool setupcalled;
352: };

354: PETSC_INTERN PetscBool      VecTaggerRegisterAllCalled;
355: PETSC_INTERN PetscErrorCode VecTaggerComputeIS_FromBoxes(VecTagger, Vec, IS *, PetscBool *);
356: PETSC_INTERN PetscMPIInt    Petsc_Reduction_keyval;

358: PETSC_INTERN PetscInt       VecGetSubVectorSavedStateId;
359: PETSC_INTERN PetscErrorCode VecGetSubVectorContiguityAndBS_Private(Vec, IS, PetscBool *, PetscInt *, PetscInt *);
360: PETSC_INTERN PetscErrorCode VecGetSubVectorThroughVecScatter_Private(Vec, IS, PetscInt, Vec *);

362: #if PetscDefined(HAVE_CUDA)
363: PETSC_INTERN PetscErrorCode VecCreate_CUDA(Vec);
364: PETSC_INTERN PetscErrorCode VecCreate_SeqCUDA(Vec);
365: PETSC_INTERN PetscErrorCode VecCreate_MPICUDA(Vec);
366: PETSC_INTERN PetscErrorCode VecCUDAGetArrays_Private(Vec, const PetscScalar **, const PetscScalar **, PetscOffloadMask *);
367: PETSC_INTERN PetscErrorCode VecConvert_Seq_SeqCUDA_inplace(Vec);
368: PETSC_INTERN PetscErrorCode VecConvert_MPI_MPICUDA_inplace(Vec);
369: #endif

371: #if PetscDefined(HAVE_HIP)
372: PETSC_INTERN PetscErrorCode VecCreate_HIP(Vec);
373: PETSC_INTERN PetscErrorCode VecCreate_SeqHIP(Vec);
374: PETSC_INTERN PetscErrorCode VecCreate_MPIHIP(Vec);
375: PETSC_INTERN PetscErrorCode VecHIPGetArrays_Private(Vec, const PetscScalar **, const PetscScalar **, PetscOffloadMask *);
376: PETSC_INTERN PetscErrorCode VecConvert_Seq_SeqHIP_inplace(Vec);
377: PETSC_INTERN PetscErrorCode VecConvert_MPI_MPIHIP_inplace(Vec);
378: #endif

380: #if defined(PETSC_HAVE_KOKKOS)
381: PETSC_INTERN PetscErrorCode VecCreateMPIKokkosWithArrays_Private(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscScalar *, const PetscScalar *, Vec *);
382: PETSC_INTERN PetscErrorCode VecConvert_Seq_SeqKokkos_inplace(Vec);
383: PETSC_INTERN PetscErrorCode VecConvert_MPI_MPIKokkos_inplace(Vec);
384: #endif

386: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecCreateWithLayout_Private(PetscLayout, Vec *);
387: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecCreateSeqWithLayoutAndArray_Private(PetscLayout, const PetscScalar *, Vec *);
388: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecCreateMPIWithLayoutAndArray_Private(PetscLayout, const PetscScalar *, Vec *);

390: /* std::upper_bound(): Given a sorted array, return index of the first element in range [first,last) whose value
391:    is greater than value, or last if there is no such element.
392: */
393: static inline PetscErrorCode PetscSortedIntUpperBound(const PetscInt *array, PetscCount first, PetscCount last, PetscInt value, PetscCount *upper)
394: {
395:   PetscCount count = last - first;

397:   PetscFunctionBegin;
398:   while (count > 0) {
399:     const PetscCount step = count / 2;
400:     PetscCount       it   = first + step;

402:     if (value < array[it]) {
403:       count = step;
404:     } else {
405:       first = it + 1;
406:       count -= step + 1;
407:     }
408:   }
409:   *upper = first;
410:   PetscFunctionReturn(PETSC_SUCCESS);
411: }

413: #define VEC_ASYNC_FN_NAME(Base) "Vec" Base "Async_Private_C"
414: #define VecAsyncFnName(Base)    VEC_##Base##_ASYNC_FN_NAME

416: #define VEC_Abs_ASYNC_FN_NAME             VEC_ASYNC_FN_NAME("Abs")
417: #define VEC_AXPBY_ASYNC_FN_NAME           VEC_ASYNC_FN_NAME("AXPBY")
418: #define VEC_AXPBYPCZ_ASYNC_FN_NAME        VEC_ASYNC_FN_NAME("AXPBYPCZ")
419: #define VEC_AXPY_ASYNC_FN_NAME            VEC_ASYNC_FN_NAME("AXPY")
420: #define VEC_AYPX_ASYNC_FN_NAME            VEC_ASYNC_FN_NAME("AYPX")
421: #define VEC_Conjugate_ASYNC_FN_NAME       VEC_ASYNC_FN_NAME("Conjugate")
422: #define VEC_Copy_ASYNC_FN_NAME            VEC_ASYNC_FN_NAME("Copy")
423: #define VEC_Exp_ASYNC_FN_NAME             VEC_ASYNC_FN_NAME("Exp")
424: #define VEC_Log_ASYNC_FN_NAME             VEC_ASYNC_FN_NAME("Log")
425: #define VEC_MAXPY_ASYNC_FN_NAME           VEC_ASYNC_FN_NAME("MAXPY")
426: #define VEC_PointwiseDivide_ASYNC_FN_NAME VEC_ASYNC_FN_NAME("PointwiseDivide")
427: #define VEC_PointwiseMax_ASYNC_FN_NAME    VEC_ASYNC_FN_NAME("PointwiseMax")
428: #define VEC_PointwiseMaxAbs_ASYNC_FN_NAME VEC_ASYNC_FN_NAME("PointwiseMaxAbs")
429: #define VEC_PointwiseMin_ASYNC_FN_NAME    VEC_ASYNC_FN_NAME("PointwiseMin")
430: #define VEC_PointwiseMult_ASYNC_FN_NAME   VEC_ASYNC_FN_NAME("PointwiseMult")
431: #define VEC_Reciprocal_ASYNC_FN_NAME      VEC_ASYNC_FN_NAME("Reciprocal")
432: #define VEC_Scale_ASYNC_FN_NAME           VEC_ASYNC_FN_NAME("Scale")
433: #define VEC_Set_ASYNC_FN_NAME             VEC_ASYNC_FN_NAME("Set")
434: #define VEC_Shift_ASYNC_FN_NAME           VEC_ASYNC_FN_NAME("Shift")
435: #define VEC_SqrtAbs_ASYNC_FN_NAME         VEC_ASYNC_FN_NAME("SqrtAbs")
436: #define VEC_Swap_ASYNC_FN_NAME            VEC_ASYNC_FN_NAME("Swap")
437: #define VEC_WAXPY_ASYNC_FN_NAME           VEC_ASYNC_FN_NAME("WAXPY")

439: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecAbsAsync_Private(Vec, PetscDeviceContext);
440: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecAXPBYAsync_Private(Vec, PetscScalar, PetscScalar, Vec, PetscDeviceContext);
441: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecAXPBYPCZAsync_Private(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec, PetscDeviceContext);
442: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecAXPYAsync_Private(Vec, PetscScalar, Vec, PetscDeviceContext);
443: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecAYPXAsync_Private(Vec, PetscScalar, Vec, PetscDeviceContext);
444: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecConjugateAsync_Private(Vec, PetscDeviceContext);
445: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecCopyAsync_Private(Vec, Vec, PetscDeviceContext);
446: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecExpAsync_Private(Vec, PetscDeviceContext);
447: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecLogAsync_Private(Vec, PetscDeviceContext);
448: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecMAXPYAsync_Private(Vec, PetscInt, const PetscScalar *, Vec[], PetscDeviceContext);
449: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecPointwiseDivideAsync_Private(Vec, Vec, Vec, PetscDeviceContext);
450: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecPointwiseMaxAsync_Private(Vec, Vec, Vec, PetscDeviceContext);
451: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecPointwiseMaxAbsAsync_Private(Vec, Vec, Vec, PetscDeviceContext);
452: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecPointwiseMinAsync_Private(Vec, Vec, Vec, PetscDeviceContext);
453: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecPointwiseMultAsync_Private(Vec, Vec, Vec, PetscDeviceContext);
454: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecReciprocalAsync_Private(Vec, PetscDeviceContext);
455: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecScaleAsync_Private(Vec, PetscScalar, PetscDeviceContext);
456: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecSetAsync_Private(Vec, PetscScalar, PetscDeviceContext);
457: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecShiftAsync_Private(Vec, PetscScalar, PetscDeviceContext);
458: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecSqrtAbsAsync_Private(Vec, PetscDeviceContext);
459: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecSwapAsync_Private(Vec, Vec, PetscDeviceContext);
460: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode VecWAXPYAsync_Private(Vec, PetscScalar, Vec, Vec, PetscDeviceContext);

462: #define VecMethodDispatch(v, dctx, async_name, name, async_arg_types, ...) \
463:   do { \
464:     PetscErrorCode(*_8_f) async_arg_types = NULL; \
465:     if (dctx) PetscCall(PetscObjectQueryFunction((PetscObject)(v), async_name, &_8_f)); \
466:     if (_8_f) { \
467:       PetscCall((*_8_f)(v, __VA_ARGS__, dctx)); \
468:     } else { \
469:       PetscUseTypeMethod(v, name, __VA_ARGS__); \
470:     } \
471:   } while (0)