Actual source code: vecimpl.h
petsc-3.8.4 2018-03-24
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: };
98: /*
99: The stash is used to temporarily store inserted vec values that
100: belong to another processor. During the assembly phase the stashed
101: values are moved to the correct processor and
102: */
104: typedef struct {
105: PetscInt nmax; /* maximum stash size */
106: PetscInt umax; /* max stash size user wants */
107: PetscInt oldnmax; /* the nmax value used previously */
108: PetscInt n; /* stash size */
109: PetscInt bs; /* block size of the stash */
110: PetscInt reallocs; /* preserve the no of mallocs invoked */
111: PetscInt *idx; /* global row numbers in stash */
112: PetscScalar *array; /* array to hold stashed values */
113: /* The following variables are used for communication */
114: MPI_Comm comm;
115: PetscMPIInt size,rank;
116: PetscMPIInt tag1,tag2;
117: MPI_Request *send_waits; /* array of send requests */
118: MPI_Request *recv_waits; /* array of receive requests */
119: MPI_Status *send_status; /* array of send status */
120: PetscInt nsends,nrecvs; /* numbers of sends and receives */
121: PetscScalar *svalues,*rvalues; /* sending and receiving data */
122: PetscInt *sindices,*rindices;
123: PetscInt rmax; /* maximum message length */
124: PetscInt *nprocs; /* tmp data used both during scatterbegin and end */
125: PetscInt nprocessed; /* number of messages already processed */
126: PetscBool donotstash;
127: PetscBool ignorenegidx; /* ignore negative indices passed into VecSetValues/VetGetValues */
128: InsertMode insertmode;
129: PetscInt *bowners;
130: } VecStash;
132: struct _p_Vec {
133: PETSCHEADER(struct _VecOps);
134: PetscLayout map;
135: void *data; /* implementation-specific data */
136: PetscBool array_gotten;
137: VecStash stash,bstash; /* used for storing off-proc values during assembly */
138: PetscBool petscnative; /* means the ->data starts with VECHEADER and can use VecGetArrayFast()*/
139: PetscInt lock; /* vector is locked to read only */
140: #if defined(PETSC_HAVE_CUSP)
141: PetscCUSPFlag valid_GPU_array; /* indicates where the most recently modified vector data is (GPU or CPU) */
142: void *spptr; /* if we're using CUSP, then this is the special pointer to the array on the GPU */
143: #elif defined(PETSC_HAVE_VIENNACL)
144: PetscViennaCLFlag valid_GPU_array; /* indicates where the most recently modified vector data is (GPU or CPU) */
145: void *spptr; /* if we're using ViennaCL, then this is the special pointer to the array on the GPU */
146: #elif defined(PETSC_HAVE_VECCUDA)
147: PetscCUDAFlag valid_GPU_array; /* indicates where the most recently modified vector data is (GPU or CPU) */
148: void *spptr; /* if we're using CUDA, then this is the special pointer to the array on the GPU */
149: #endif
150: };
152: PETSC_EXTERN PetscLogEvent VEC_SetRandom;
153: PETSC_EXTERN PetscLogEvent VEC_View, VEC_Max, VEC_Min, VEC_DotBarrier, VEC_Dot, VEC_MDotBarrier, VEC_MDot, VEC_TDot, VEC_MTDot;
154: PETSC_EXTERN PetscLogEvent VEC_Norm, VEC_Normalize, VEC_Scale, VEC_Copy, VEC_Set, VEC_AXPY, VEC_AYPX, VEC_WAXPY, VEC_MAXPY;
155: PETSC_EXTERN PetscLogEvent VEC_AssemblyEnd, VEC_PointwiseMult, VEC_SetValues, VEC_Load, VEC_ScatterBarrier, VEC_ScatterBegin, VEC_ScatterEnd;
156: PETSC_EXTERN PetscLogEvent VEC_ReduceArithmetic, VEC_ReduceBarrier, VEC_ReduceCommunication;
157: PETSC_EXTERN PetscLogEvent VEC_ReduceBegin,VEC_ReduceEnd;
158: PETSC_EXTERN PetscLogEvent VEC_Swap, VEC_AssemblyBegin, VEC_NormBarrier, VEC_DotNormBarrier, VEC_DotNorm, VEC_AXPBYPCZ, VEC_Ops;
159: PETSC_EXTERN PetscLogEvent VEC_CUSPCopyToGPU, VEC_CUSPCopyFromGPU;
160: PETSC_EXTERN PetscLogEvent VEC_CUSPCopyToGPUSome, VEC_CUSPCopyFromGPUSome;
161: PETSC_EXTERN PetscLogEvent VEC_ViennaCLCopyToGPU, VEC_ViennaCLCopyFromGPU;
162: PETSC_EXTERN PetscLogEvent VEC_CUDACopyToGPU, VEC_CUDACopyFromGPU;
163: PETSC_EXTERN PetscLogEvent VEC_CUDACopyToGPUSome, VEC_CUDACopyFromGPUSome;
165: PETSC_EXTERN PetscErrorCode VecView_Seq(Vec,PetscViewer);
166: #if defined(PETSC_HAVE_CUSP)
167: PETSC_EXTERN PetscErrorCode VecCUSPAllocateCheckHost(Vec v);
168: PETSC_EXTERN PetscErrorCode VecCUSPCopyFromGPU(Vec v);
169: #elif defined(PETSC_HAVE_VIENNACL)
170: PETSC_EXTERN PetscErrorCode VecViennaCLAllocateCheckHost(Vec v);
171: PETSC_EXTERN PetscErrorCode VecViennaCLCopyFromGPU(Vec v);
172: #elif defined(PETSC_HAVE_VECCUDA)
173: PETSC_EXTERN PetscErrorCode VecCUDAAllocateCheckHost(Vec v);
174: PETSC_EXTERN PetscErrorCode VecCUDACopyFromGPU(Vec v);
175: #endif
178: /*
179: Common header shared by array based vectors,
180: currently Vec_Seq and Vec_MPI
181: */
182: #define VECHEADER \
183: PetscScalar *array; \
184: PetscScalar *array_allocated; /* if the array was allocated by PETSc this is its pointer */ \
185: PetscScalar *unplacedarray; /* if one called VecPlaceArray(), this is where it stashed the original */
187: /* Default obtain and release vectors; can be used by any implementation */
188: PETSC_EXTERN PetscErrorCode VecDuplicateVecs_Default(Vec,PetscInt,Vec *[]);
189: PETSC_EXTERN PetscErrorCode VecDestroyVecs_Default(PetscInt,Vec []);
190: PETSC_INTERN PetscErrorCode VecLoad_Binary(Vec, PetscViewer);
191: PETSC_EXTERN PetscErrorCode VecLoad_Default(Vec, PetscViewer);
193: PETSC_EXTERN PetscInt NormIds[7]; /* map from NormType to IDs used to cache/retreive values of norms */
195: /* --------------------------------------------------------------------*/
196: /* */
197: /* Defines the data structures used in the Vec Scatter operations */
199: typedef enum { VEC_SCATTER_SEQ_GENERAL,VEC_SCATTER_SEQ_STRIDE,
200: VEC_SCATTER_MPI_GENERAL,VEC_SCATTER_MPI_TOALL,
201: VEC_SCATTER_MPI_TOONE} VecScatterType;
203: #define VECSCATTER_IMPL_HEADER \
204: VecScatterType type;
206: typedef struct {
207: VECSCATTER_IMPL_HEADER
208: } VecScatter_Common;
210: /*
211: These scatters are for the purely local case.
212: */
213: typedef struct {
214: VECSCATTER_IMPL_HEADER
215: PetscInt n; /* number of components to scatter */
216: PetscInt *vslots; /* locations of components */
217: /*
218: The next three fields are used in parallel scatters, they contain
219: optimization in the special case that the "to" vector and the "from"
220: vector are the same, so one only needs copy components that truly
221: copies instead of just y[idx[i]] = y[jdx[i]] where idx[i] == jdx[i].
222: */
223: PetscBool nonmatching_computed;
224: PetscInt n_nonmatching; /* number of "from"s != "to"s */
225: PetscInt *slots_nonmatching; /* locations of "from"s != "to"s */
226: PetscBool is_copy;
227: PetscInt copy_start; /* local scatter is a copy starting at copy_start */
228: PetscInt copy_length;
229: } VecScatter_Seq_General;
231: typedef struct {
232: VECSCATTER_IMPL_HEADER
233: PetscInt n;
234: PetscInt first;
235: PetscInt step;
236: } VecScatter_Seq_Stride;
238: /*
239: This scatter is for a global vector copied (completely) to each processor (or all to one)
240: */
241: typedef struct {
242: VECSCATTER_IMPL_HEADER
243: PetscMPIInt *count; /* elements of vector on each processor */
244: PetscMPIInt *displx;
245: PetscScalar *work1;
246: PetscScalar *work2;
247: } VecScatter_MPI_ToAll;
249: /*
250: This is the general parallel scatter
251: */
252: typedef struct {
253: VECSCATTER_IMPL_HEADER
254: PetscInt n; /* number of processors to send/receive */
255: PetscInt *starts; /* starting point in indices and values for each proc*/
256: PetscInt *indices; /* list of all components sent or received */
257: PetscMPIInt *procs; /* processors we are communicating with in scatter */
258: MPI_Request *requests,*rev_requests;
259: PetscScalar *values; /* buffer for all sends or receives */
260: VecScatter_Seq_General local; /* any part that happens to be local */
261: MPI_Status *sstatus,*rstatus;
262: PetscBool use_readyreceiver;
263: PetscInt bs;
264: PetscBool sendfirst;
265: PetscBool contiq;
266: /* for MPI_Alltoallv() approach */
267: PetscBool use_alltoallv;
268: PetscMPIInt *counts,*displs;
269: /* for MPI_Alltoallw() approach */
270: PetscBool use_alltoallw;
271: #if defined(PETSC_HAVE_MPI_ALLTOALLW)
272: PetscMPIInt *wcounts,*wdispls;
273: MPI_Datatype *types;
274: #endif
275: PetscBool use_window;
276: #if defined(PETSC_HAVE_MPI_WIN_CREATE)
277: MPI_Win window;
278: PetscInt *winstarts; /* displacements in the processes I am putting to */
279: #endif
280: } VecScatter_MPI_General;
283: PETSC_INTERN PetscErrorCode VecScatterGetTypes_Private(VecScatter,VecScatterType*,VecScatterType*);
284: PETSC_INTERN PetscErrorCode VecScatterIsSequential_Private(VecScatter_Common*,PetscBool*);
286: typedef struct _VecScatterOps *VecScatterOps;
287: struct _VecScatterOps {
288: PetscErrorCode (*begin)(VecScatter,Vec,Vec,InsertMode,ScatterMode);
289: PetscErrorCode (*end)(VecScatter,Vec,Vec,InsertMode,ScatterMode);
290: PetscErrorCode (*copy)(VecScatter,VecScatter);
291: PetscErrorCode (*destroy)(VecScatter);
292: PetscErrorCode (*view)(VecScatter,PetscViewer);
293: PetscErrorCode (*viewfromoptions)(VecScatter,const char prefix[],const char name[]);
294: PetscErrorCode (*remap)(VecScatter,PetscInt *,PetscInt*);
295: PetscErrorCode (*getmerged)(VecScatter,PetscBool *);
296: };
298: struct _p_VecScatter {
299: PETSCHEADER(struct _VecScatterOps);
300: PetscInt to_n,from_n;
301: PetscBool inuse; /* prevents corruption from mixing two scatters */
302: PetscBool beginandendtogether; /* indicates that the scatter begin and end function are called together, VecScatterEnd()
303: is then treated as a nop */
304: PetscBool packtogether; /* packs all the messages before sending, same with receive */
305: PetscBool reproduce; /* always receive the ghost points in the same order of processes */
306: void *fromdata,*todata;
307: void *spptr;
308: };
310: PETSC_INTERN PetscErrorCode VecStashCreate_Private(MPI_Comm,PetscInt,VecStash*);
311: PETSC_INTERN PetscErrorCode VecStashDestroy_Private(VecStash*);
312: PETSC_INTERN PetscErrorCode VecStashExpand_Private(VecStash*,PetscInt);
313: PETSC_INTERN PetscErrorCode VecStashScatterEnd_Private(VecStash*);
314: PETSC_INTERN PetscErrorCode VecStashSetInitialSize_Private(VecStash*,PetscInt);
315: PETSC_INTERN PetscErrorCode VecStashGetInfo_Private(VecStash*,PetscInt*,PetscInt*);
316: PETSC_INTERN PetscErrorCode VecStashScatterBegin_Private(VecStash*,PetscInt*);
317: PETSC_INTERN PetscErrorCode VecStashScatterGetMesg_Private(VecStash*,PetscMPIInt*,PetscInt**,PetscScalar**,PetscInt*);
318: PETSC_INTERN PetscErrorCode VecStashSortCompress_Private(VecStash*);
319: PETSC_INTERN PetscErrorCode VecStashGetOwnerList_Private(VecStash*,PetscLayout,PetscMPIInt*,PetscMPIInt**);
321: /*
322: VecStashValue_Private - inserts a single value into the stash.
324: Input Parameters:
325: stash - the stash
326: idx - the global of the inserted value
327: values - the value inserted
328: */
329: PETSC_STATIC_INLINE PetscErrorCode VecStashValue_Private(VecStash *stash,PetscInt row,PetscScalar value)
330: {
332: /* Check and see if we have sufficient memory */
333: if (((stash)->n + 1) > (stash)->nmax) {
334: VecStashExpand_Private(stash,1);
335: }
336: (stash)->idx[(stash)->n] = row;
337: (stash)->array[(stash)->n] = value;
338: (stash)->n++;
339: return 0;
340: }
342: /*
343: VecStashValuesBlocked_Private - inserts 1 block of values into the stash.
345: Input Parameters:
346: stash - the stash
347: idx - the global block index
348: values - the values inserted
349: */
350: PETSC_STATIC_INLINE PetscErrorCode VecStashValuesBlocked_Private(VecStash *stash,PetscInt row,PetscScalar *values)
351: {
352: PetscInt jj,stash_bs=(stash)->bs;
353: PetscScalar *array;
355: if (((stash)->n+1) > (stash)->nmax) {
356: VecStashExpand_Private(stash,1);
357: }
358: array = (stash)->array + stash_bs*(stash)->n;
359: (stash)->idx[(stash)->n] = row;
360: for (jj=0; jj<stash_bs; jj++) { array[jj] = values[jj];}
361: (stash)->n++;
362: return 0;
363: }
365: PETSC_INTERN PetscErrorCode VecStrideGather_Default(Vec,PetscInt,Vec,InsertMode);
366: PETSC_INTERN PetscErrorCode VecStrideScatter_Default(Vec,PetscInt,Vec,InsertMode);
367: PETSC_INTERN PetscErrorCode VecReciprocal_Default(Vec);
368: PETSC_INTERN PetscErrorCode VecStrideSubSetGather_Default(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);
369: PETSC_INTERN PetscErrorCode VecStrideSubSetScatter_Default(Vec,PetscInt,const PetscInt[],const PetscInt[],Vec,InsertMode);
371: #if defined(PETSC_HAVE_MATLAB_ENGINE)
372: PETSC_EXTERN PetscErrorCode VecMatlabEnginePut_Default(PetscObject,void*);
373: PETSC_EXTERN PetscErrorCode VecMatlabEngineGet_Default(PetscObject,void*);
374: #endif
376: PETSC_EXTERN PetscErrorCode PetscSectionGetField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);
377: PETSC_EXTERN PetscErrorCode PetscSectionRestoreField_Internal(PetscSection, PetscSection, Vec, PetscInt, PetscInt, PetscInt, IS *, Vec *);
379: #define VecCheckSameLocalSize(x,ar1,y,ar2) \
380: 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);
382: #define VecCheckSameSize(x,ar1,y,ar2) \
383: if ((x)->map->N != (y)->map->N) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector global lengths parameter # %d %D != paramter # %d %D", ar1,(x)->map->N, ar2,(y)->map->N); \
384: if ((x)->map->n != (y)->map->n) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible vector local lengths parameter # %d %D != parameter # %d %D", ar1,(x)->map->n, ar2,(y)->map->n);
386: typedef struct _VecTaggerOps *VecTaggerOps;
387: struct _VecTaggerOps {
388: PetscErrorCode (*create) (VecTagger);
389: PetscErrorCode (*destroy) (VecTagger);
390: PetscErrorCode (*setfromoptions) (PetscOptionItems*,VecTagger);
391: PetscErrorCode (*setup) (VecTagger);
392: PetscErrorCode (*view) (VecTagger,PetscViewer);
393: PetscErrorCode (*computeboxes) (VecTagger,Vec,PetscInt *,VecTaggerBox **);
394: PetscErrorCode (*computeis) (VecTagger,Vec,IS *);
395: };
396: struct _p_VecTagger {
397: PETSCHEADER(struct _VecTaggerOps);
398: void *data;
399: PetscInt blocksize;
400: PetscBool invert;
401: PetscBool setupcalled;
402: };
404: PETSC_EXTERN PetscBool VecTaggerRegisterAllCalled;
405: PETSC_EXTERN PetscErrorCode VecTaggerRegisterAll(void);
406: PETSC_EXTERN PetscErrorCode VecTaggerComputeIS_FromBoxes(VecTagger,Vec,IS*);
407: PETSC_EXTERN PetscMPIInt Petsc_Reduction_keyval;
409: #endif