xref: /petsc/src/dm/impls/moab/dmmbvec.cxx (revision 0369a392ee39ea4d01e519e42e389892d4e1e2e4)
1 #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2 #include <petsc/private/vecimpl.h>
3 
4 #include <petscdmmoab.h>
5 #include <MBTagConventions.hpp>
6 
7 #define USE_NATIVE_PETSCVEC
8 
9 /* declare some private DMMoab specific overrides */
10 static PetscErrorCode DMCreateVector_Moab_Private(DM dm,moab::Tag tag,const moab::Range* userrange,PetscBool is_global_vec,PetscBool destroy_tag,Vec *vec);
11 static PetscErrorCode DMVecUserDestroy_Moab(void *user);
12 static PetscErrorCode DMVecDuplicate_Moab(Vec x,Vec *y);
13 #ifdef MOAB_HAVE_MPI
14 static PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface,moab::ParallelComm *pcomm,char** tag_name);
15 #else
16 static PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface,char** tag_name);
17 #endif
18 
19 /*@C
20   DMMoabCreateVector - Create a Vec from either an existing tag, or a specified tag size, and a range of entities
21 
22   Collective
23 
24   Input Parameters:
25 + dm              - The DMMoab object being set
26 . tag             - If non-zero, block size will be taken from the tag size
27 . range           - If non-empty, Vec corresponds to these entities, otherwise to the entities set on the DMMoab
28 . is_global_vec   - If true, this is a local representation of the Vec (including ghosts in parallel), otherwise a truly parallel one
29 - destroy_tag     - If true, MOAB tag is destroyed with Vec, otherwise it is left on MOAB
30 
31   Output Parameter:
32 . vec             - The created vector
33 
34   Level: beginner
35 
36 .seealso: VecCreate()
37 @*/
38 PetscErrorCode DMMoabCreateVector(DM dm, moab::Tag tag, const moab::Range* range, PetscBool is_global_vec, PetscBool destroy_tag, Vec *vec)
39 {
40   PetscFunctionBegin;
41   PetscCheck(tag || (range && !range->empty()),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Both tag and range cannot be null.");
42 
43   PetscCall(DMCreateVector_Moab_Private(dm, tag, range, is_global_vec, destroy_tag, vec));
44   PetscFunctionReturn(0);
45 }
46 
47 /*@C
48   DMMoabGetVecTag - Get the MOAB tag associated with this Vec
49 
50   Input Parameter:
51 . vec - Vec being queried
52 
53   Output Parameter:
54 . tag - Tag associated with this Vec. NULL if vec is a native PETSc Vec.
55 
56   Level: beginner
57 
58 .seealso: DMMoabCreateVector(), DMMoabGetVecRange()
59 @*/
60 PetscErrorCode DMMoabGetVecTag(Vec vec, moab::Tag *tag)
61 {
62   PetscContainer  moabdata;
63   Vec_MOAB        *vmoab;
64 
65   PetscFunctionBegin;
66   PetscValidPointer(tag, 2);
67 
68   /* Get the MOAB private data */
69   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
70   PetscCall(PetscContainerGetPointer(moabdata, (void**) &vmoab));
71 
72   *tag = vmoab->tag;
73   PetscFunctionReturn(0);
74 }
75 
76 /*@C
77   DMMoabGetVecRange - Get the MOAB entities associated with this Vec
78 
79   Input Parameter:
80 . vec   - Vec being queried
81 
82   Output Parameter:
83 . range - Entities associated with this Vec. NULL if vec is a native PETSc Vec.
84 
85   Level: beginner
86 
87 .seealso: DMMoabCreateVector(), DMMoabGetVecTag()
88 @*/
89 PetscErrorCode DMMoabGetVecRange(Vec vec, moab::Range *range)
90 {
91   PetscContainer  moabdata;
92   Vec_MOAB        *vmoab;
93 
94   PetscFunctionBegin;
95   PetscValidPointer(range, 2);
96 
97   /* Get the MOAB private data handle */
98   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
99   PetscCall(PetscContainerGetPointer(moabdata, (void**) &vmoab));
100 
101   *range = *vmoab->tag_range;
102   PetscFunctionReturn(0);
103 }
104 
105 /*@C
106   DMMoabVecGetArray - Returns the writable direct access array to the local representation of MOAB tag data for the underlying vector using locally owned+ghosted range of entities
107 
108   Collective
109 
110   Input Parameters:
111 + dm              - The DMMoab object being set
112 - vec             - The Vector whose underlying data is requested
113 
114   Output Parameter:
115 . array           - The local data array
116 
117   Level: intermediate
118 
119 .seealso: DMMoabVecRestoreArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead()
120 @*/
121 PetscErrorCode  DMMoabVecGetArray(DM dm, Vec vec, void* array)
122 {
123   DM_Moab        *dmmoab;
124   moab::ErrorCode merr;
125   PetscInt        count, i, f;
126   moab::Tag       vtag;
127   PetscScalar     **varray;
128   PetscScalar     *marray;
129   PetscContainer  moabdata;
130   Vec_MOAB        *vmoab, *xmoab;
131 
132   PetscFunctionBegin;
133   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
134   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
135   PetscValidPointer(array, 3);
136   dmmoab = (DM_Moab*)dm->data;
137 
138   /* Get the Vec_MOAB struct for the original vector */
139   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
140   PetscCall(PetscContainerGetPointer(moabdata, (void**)&vmoab));
141 
142   /* Get the real scalar array handle */
143   varray = reinterpret_cast<PetscScalar**>(array);
144 
145   if (vmoab->is_native_vec) {
146 
147     /* get the local representation of the arrays from Vectors */
148     PetscCall(DMGetLocalVector(dm, &vmoab->local));
149     PetscCall(DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vmoab->local));
150     PetscCall(DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vmoab->local));
151 
152     /* Get the Vec_MOAB struct for the original vector */
153     PetscCall(PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata));
154     PetscCall(PetscContainerGetPointer(moabdata, (void**)&xmoab));
155 
156     /* get the local representation of the arrays from Vectors */
157     PetscCall(VecGhostGetLocalForm(vmoab->local, &xmoab->local));
158     PetscCall(VecGhostUpdateBegin(vmoab->local, INSERT_VALUES, SCATTER_FORWARD));
159     PetscCall(VecGhostUpdateEnd(vmoab->local, INSERT_VALUES, SCATTER_FORWARD));
160 
161     PetscCall(VecGetArray(xmoab->local, varray));
162   }
163   else {
164 
165     /* Get the MOAB private data */
166     PetscCall(DMMoabGetVecTag(vec, &vtag));
167 
168 #ifdef MOAB_HAVE_MPI
169     /* exchange the data into ghost cells first */
170     merr = dmmoab->pcomm->exchange_tags(vtag, *dmmoab->vlocal); MBERRNM(merr);
171 #endif
172 
173     PetscCall(PetscMalloc1((dmmoab->nloc + dmmoab->nghost) * dmmoab->numFields, varray));
174 
175     /* Get the array data for local entities */
176     merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr);
177     PetscCheckFalse(count != dmmoab->nloc + dmmoab->nghost,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost);
178 
179     i = 0;
180     for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
181       for (f = 0; f < dmmoab->numFields; f++, i++)
182         (*varray)[dmmoab->lidmap[(PetscInt)*iter - dmmoab->seqstart]*dmmoab->numFields + f] = marray[i];
183       //(*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]]=marray[i];
184     }
185   }
186   PetscFunctionReturn(0);
187 }
188 
189 /*@C
190   DMMoabVecRestoreArray - Restores the writable direct access array obtained via DMMoabVecGetArray
191 
192   Collective
193 
194   Input Parameters:
195 + dm              - The DMMoab object being set
196 . vec             - The Vector whose underlying data is requested
197 - array           - The local data array
198 
199   Level: intermediate
200 
201 .seealso: DMMoabVecGetArray(), DMMoabVecGetArrayRead(), DMMoabVecRestoreArrayRead()
202 @*/
203 PetscErrorCode  DMMoabVecRestoreArray(DM dm, Vec vec, void* array)
204 {
205   DM_Moab        *dmmoab;
206   moab::ErrorCode merr;
207   moab::Tag       vtag;
208   PetscInt        count, i, f;
209   PetscScalar     **varray;
210   PetscScalar     *marray;
211   PetscContainer  moabdata;
212   Vec_MOAB        *vmoab, *xmoab;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
216   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
217   PetscValidPointer(array, 3);
218   dmmoab = (DM_Moab*)dm->data;
219 
220   /* Get the Vec_MOAB struct for the original vector */
221   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
222   PetscCall(PetscContainerGetPointer(moabdata, (void**)&vmoab));
223 
224   /* Get the real scalar array handle */
225   varray = reinterpret_cast<PetscScalar**>(array);
226 
227   if (vmoab->is_native_vec) {
228 
229     /* Get the Vec_MOAB struct for the original vector */
230     PetscCall(PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata));
231     PetscCall(PetscContainerGetPointer(moabdata, (void**)&xmoab));
232 
233     /* get the local representation of the arrays from Vectors */
234     PetscCall(VecRestoreArray(xmoab->local, varray));
235     PetscCall(VecGhostUpdateBegin(vmoab->local, ADD_VALUES, SCATTER_REVERSE));
236     PetscCall(VecGhostUpdateEnd(vmoab->local, ADD_VALUES, SCATTER_REVERSE));
237     PetscCall(VecGhostRestoreLocalForm(vmoab->local, &xmoab->local));
238 
239     /* restore local pieces */
240     PetscCall(DMLocalToGlobalBegin(dm, vmoab->local, INSERT_VALUES, vec));
241     PetscCall(DMLocalToGlobalEnd(dm, vmoab->local, INSERT_VALUES, vec));
242     PetscCall(DMRestoreLocalVector(dm, &vmoab->local));
243   }
244   else {
245 
246     /* Get the MOAB private data */
247     PetscCall(DMMoabGetVecTag(vec, &vtag));
248 
249     /* Get the array data for local entities */
250     merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr);
251     PetscCheckFalse(count != dmmoab->nloc + dmmoab->nghost,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost);
252 
253     i = 0;
254     for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
255       for (f = 0; f < dmmoab->numFields; f++, i++)
256         marray[i] = (*varray)[dmmoab->lidmap[(PetscInt) * iter - dmmoab->seqstart] * dmmoab->numFields + f];
257       //marray[i] = (*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]];
258     }
259 
260 #ifdef MOAB_HAVE_MPI
261     /* reduce the tags correctly -> should probably let the user choose how to reduce in the future
262       For all FEM residual based assembly calculations, MPI_SUM should serve well */
263     merr = dmmoab->pcomm->reduce_tags(vtag, MPI_SUM, *dmmoab->vlocal); MBERRV(dmmoab->mbiface, merr);
264 #endif
265     PetscCall(PetscFree(*varray));
266   }
267   PetscFunctionReturn(0);
268 }
269 
270 /*@C
271   DMMoabVecGetArrayRead - Returns the read-only direct access array to the local representation of MOAB tag data for the underlying vector using locally owned+ghosted range of entities
272 
273   Collective
274 
275   Input Parameters:
276 + dm              - The DMMoab object being set
277 - vec             - The Vector whose underlying data is requested
278 
279   Output Parameter:
280 . array           - The local data array
281 
282   Level: intermediate
283 
284 .seealso: DMMoabVecRestoreArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray()
285 @*/
286 PetscErrorCode  DMMoabVecGetArrayRead(DM dm, Vec vec, void* array)
287 {
288   DM_Moab        *dmmoab;
289   moab::ErrorCode merr;
290   PetscInt        count, i, f;
291   moab::Tag       vtag;
292   PetscScalar     **varray;
293   PetscScalar     *marray;
294   PetscContainer  moabdata;
295   Vec_MOAB        *vmoab, *xmoab;
296 
297   PetscFunctionBegin;
298   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
299   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
300   PetscValidPointer(array, 3);
301   dmmoab = (DM_Moab*)dm->data;
302 
303   /* Get the Vec_MOAB struct for the original vector */
304   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
305   PetscCall(PetscContainerGetPointer(moabdata, (void**)&vmoab));
306 
307   /* Get the real scalar array handle */
308   varray = reinterpret_cast<PetscScalar**>(array);
309 
310   if (vmoab->is_native_vec) {
311     /* get the local representation of the arrays from Vectors */
312     PetscCall(DMGetLocalVector(dm, &vmoab->local));
313     PetscCall(DMGlobalToLocalBegin(dm, vec, INSERT_VALUES, vmoab->local));
314     PetscCall(DMGlobalToLocalEnd(dm, vec, INSERT_VALUES, vmoab->local));
315 
316     /* Get the Vec_MOAB struct for the original vector */
317     PetscCall(PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata));
318     PetscCall(PetscContainerGetPointer(moabdata, (void**)&xmoab));
319 
320     /* get the local representation of the arrays from Vectors */
321     PetscCall(VecGhostGetLocalForm(vmoab->local, &xmoab->local));
322     PetscCall(VecGhostUpdateBegin(vmoab->local, INSERT_VALUES, SCATTER_FORWARD));
323     PetscCall(VecGhostUpdateEnd(vmoab->local, INSERT_VALUES, SCATTER_FORWARD));
324     PetscCall(VecGetArray(xmoab->local, varray));
325   }
326   else {
327     /* Get the MOAB private data */
328     PetscCall(DMMoabGetVecTag(vec, &vtag));
329 
330 #ifdef MOAB_HAVE_MPI
331     /* exchange the data into ghost cells first */
332     merr = dmmoab->pcomm->exchange_tags(vtag, *dmmoab->vlocal); MBERRNM(merr);
333 #endif
334     PetscCall(PetscMalloc1((dmmoab->nloc + dmmoab->nghost) * dmmoab->numFields, varray));
335 
336     /* Get the array data for local entities */
337     merr = dmmoab->mbiface->tag_iterate(vtag, dmmoab->vlocal->begin(), dmmoab->vlocal->end(), count, reinterpret_cast<void*&>(marray), false); MBERRNM(merr);
338     PetscCheckFalse(count != dmmoab->nloc + dmmoab->nghost,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mismatch between local vertices and tag partition for Vec. %D != %D.", count, dmmoab->nloc + dmmoab->nghost);
339 
340     i = 0;
341     for (moab::Range::iterator iter = dmmoab->vlocal->begin(); iter != dmmoab->vlocal->end(); iter++) {
342       for (f = 0; f < dmmoab->numFields; f++, i++)
343         (*varray)[dmmoab->lidmap[(PetscInt)*iter - dmmoab->seqstart]*dmmoab->numFields + f] = marray[i];
344       //(*varray)[dmmoab->llmap[dmmoab->lidmap[((PetscInt)*iter-dmmoab->seqstart)]*dmmoab->numFields+f]]=marray[i];
345     }
346   }
347   PetscFunctionReturn(0);
348 }
349 
350 /*@C
351   DMMoabVecRestoreArrayRead - Restores the read-only direct access array obtained via DMMoabVecGetArray
352 
353   Collective
354 
355   Input Parameters:
356 + dm              - The DMMoab object being set
357 . vec             - The Vector whose underlying data is requested
358 - array           - The local data array
359 
360   Level: intermediate
361 
362 .seealso: DMMoabVecGetArrayRead(), DMMoabVecGetArray(), DMMoabVecRestoreArray()
363 @*/
364 PetscErrorCode  DMMoabVecRestoreArrayRead(DM dm, Vec vec, void* array)
365 {
366   PetscScalar     **varray;
367   PetscContainer  moabdata;
368   Vec_MOAB        *vmoab, *xmoab;
369 
370   PetscFunctionBegin;
371   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
372   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
373   PetscValidPointer(array, 3);
374 
375   /* Get the Vec_MOAB struct for the original vector */
376   PetscCall(PetscObjectQuery((PetscObject)vec, "MOABData", (PetscObject*) &moabdata));
377   PetscCall(PetscContainerGetPointer(moabdata, (void**)&vmoab));
378 
379   /* Get the real scalar array handle */
380   varray = reinterpret_cast<PetscScalar**>(array);
381 
382   if (vmoab->is_native_vec) {
383     /* Get the Vec_MOAB struct for the original vector */
384     PetscCall(PetscObjectQuery((PetscObject)vmoab->local, "MOABData", (PetscObject*) &moabdata));
385     PetscCall(PetscContainerGetPointer(moabdata, (void**)&xmoab));
386 
387     /* restore the local representation of the arrays from Vectors */
388     PetscCall(VecRestoreArray(xmoab->local, varray));
389     PetscCall(VecGhostRestoreLocalForm(vmoab->local, &xmoab->local));
390 
391     /* restore local pieces */
392     PetscCall(DMRestoreLocalVector(dm, &vmoab->local));
393   }
394   else {
395     /* Nothing to do but just free the memory allocated before */
396     PetscCall(PetscFree(*varray));
397 
398   }
399   PetscFunctionReturn(0);
400 }
401 
402 PetscErrorCode DMCreateVector_Moab_Private(DM dm, moab::Tag tag, const moab::Range* userrange, PetscBool is_global_vec, PetscBool destroy_tag, Vec *vec)
403 {
404   moab::ErrorCode   merr;
405   PetscBool         is_newtag;
406   const moab::Range *range;
407   PetscInt          count, lnative_vec, gnative_vec;
408   std::string       ttname;
409   PetscScalar       *data_ptr, *defaultvals;
410 
411   Vec_MOAB *vmoab;
412   DM_Moab *dmmoab = (DM_Moab*)dm->data;
413 #ifdef MOAB_HAVE_MPI
414   moab::ParallelComm *pcomm = dmmoab->pcomm;
415 #endif
416   moab::Interface *mbiface = dmmoab->mbiface;
417 
418   PetscFunctionBegin;
419   PetscCheckFalse(sizeof(PetscReal) != sizeof(PetscScalar),PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "MOAB tags only support Real types (Complex-type unsupported)");
420   if (!userrange) range = dmmoab->vowned;
421   else range = userrange;
422   PetscCheck(range,PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Input range cannot be empty or call DMSetUp first.");
423 
424 #ifndef USE_NATIVE_PETSCVEC
425   /* If the tag data is in a single sequence, use PETSc native vector since tag_iterate isn't useful anymore */
426   lnative_vec = (range->psize() - 1);
427 #else
428   lnative_vec = 1; /* NOTE: Testing PETSc vector: will force to create native vector all the time */
429 //  lnative_vec=0; /* NOTE: Testing MOAB vector: will force to create MOAB tag_iterate based vector all the time */
430 #endif
431 
432 #ifdef MOAB_HAVE_MPI
433   PetscCall(MPIU_Allreduce(&lnative_vec, &gnative_vec, 1, MPI_INT, MPI_MAX, (((PetscObject)dm)->comm)));
434 #else
435   gnative_vec = lnative_vec;
436 #endif
437 
438   /* Create the MOAB internal data object */
439   PetscCall(PetscNew(&vmoab));
440   vmoab->is_native_vec = (gnative_vec > 0 ? PETSC_TRUE : PETSC_FALSE);
441 
442   if (!vmoab->is_native_vec) {
443     merr = moab::MB_SUCCESS;
444     if (tag != 0) merr = mbiface->tag_get_name(tag, ttname);
445     if (!ttname.length() || merr != moab::MB_SUCCESS) {
446       /* get the new name for the anonymous MOABVec -> the tag_name will be destroyed along with Tag */
447       char *tag_name = NULL;
448 #ifdef MOAB_HAVE_MPI
449       PetscCall(DMVecCreateTagName_Moab_Private(mbiface,pcomm,&tag_name));
450 #else
451       PetscCall(DMVecCreateTagName_Moab_Private(mbiface,&tag_name));
452 #endif
453       is_newtag = PETSC_TRUE;
454 
455       /* Create the default value for the tag (all zeros) */
456       PetscCall(PetscCalloc1(dmmoab->numFields, &defaultvals));
457 
458       /* Create the tag */
459       merr = mbiface->tag_get_handle(tag_name, dmmoab->numFields, moab::MB_TYPE_DOUBLE, tag,
460                                      moab::MB_TAG_DENSE | moab::MB_TAG_CREAT, defaultvals); MBERRNM(merr);
461       PetscCall(PetscFree(tag_name));
462       PetscCall(PetscFree(defaultvals));
463     }
464     else {
465       /* Make sure the tag data is of type "double" */
466       moab::DataType tag_type;
467       merr = mbiface->tag_get_data_type(tag, tag_type); MBERRNM(merr);
468       PetscCheckFalse(tag_type != moab::MB_TYPE_DOUBLE,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Tag data type must be MB_TYPE_DOUBLE");
469       is_newtag = destroy_tag;
470     }
471 
472     vmoab->tag = tag;
473     vmoab->new_tag = is_newtag;
474   }
475   vmoab->mbiface = mbiface;
476 #ifdef MOAB_HAVE_MPI
477   vmoab->pcomm = pcomm;
478 #endif
479   vmoab->is_global_vec = is_global_vec;
480   vmoab->tag_size = dmmoab->bs;
481 
482   if (vmoab->is_native_vec) {
483 
484     /* Create the PETSc Vector directly and attach our functions accordingly */
485     if (!is_global_vec) {
486       /* This is an MPI Vector with ghosted padding */
487       PetscCall(VecCreateGhostBlock((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * dmmoab->nloc,dmmoab->numFields * dmmoab->n, dmmoab->nghost, &dmmoab->gsindices[dmmoab->nloc], vec));
488     } else {
489       /* This is an MPI/SEQ Vector */
490       PetscCall(VecCreate((((PetscObject)dm)->comm), vec));
491       PetscCall(VecSetSizes(*vec, dmmoab->numFields * dmmoab->nloc, PETSC_DECIDE));
492       PetscCall(VecSetBlockSize(*vec, dmmoab->bs));
493       PetscCall(VecSetType(*vec, VECMPI));
494     }
495   } else {
496     /* Call tag_iterate. This will cause MOAB to allocate memory for the
497        tag data if it hasn't already happened */
498     merr = mbiface->tag_iterate(tag, range->begin(), range->end(), count, (void*&)data_ptr); MBERRNM(merr);
499 
500     /* set the reference for vector range */
501     vmoab->tag_range = new moab::Range(*range);
502     merr = mbiface->tag_get_length(tag, dmmoab->numFields); MBERRNM(merr);
503 
504     /* Create the PETSc Vector
505       Query MOAB mesh to check if there are any ghosted entities
506         -> if we do, create a ghosted vector to map correctly to the same layout
507         -> else, create a non-ghosted parallel vector */
508     if (!is_global_vec) {
509       /* This is an MPI Vector with ghosted padding */
510       PetscCall(VecCreateGhostBlockWithArray((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * dmmoab->nloc,dmmoab->numFields * dmmoab->n, dmmoab->nghost, &dmmoab->gsindices[dmmoab->nloc], data_ptr, vec));
511     } else {
512       /* This is an MPI Vector without ghosted padding */
513       PetscCall(VecCreateMPIWithArray((((PetscObject)dm)->comm), dmmoab->bs, dmmoab->numFields * range->size(),PETSC_DECIDE, data_ptr, vec));
514     }
515   }
516   PetscCall(VecSetFromOptions(*vec));
517 
518   /* create a container and store the internal MOAB data for faster access based on Entities etc */
519   PetscContainer moabdata;
520   PetscCall(PetscContainerCreate(PETSC_COMM_WORLD, &moabdata));
521   PetscCall(PetscContainerSetPointer(moabdata, vmoab));
522   PetscCall(PetscContainerSetUserDestroy(moabdata, DMVecUserDestroy_Moab));
523   PetscCall(PetscObjectCompose((PetscObject) * vec, "MOABData", (PetscObject)moabdata));
524   (*vec)->ops->duplicate = DMVecDuplicate_Moab;
525   PetscCall(PetscContainerDestroy(&moabdata));
526 
527   /* Vector created, manually set local to global mapping */
528   if (dmmoab->ltog_map) {
529     PetscCall(VecSetLocalToGlobalMapping(*vec, dmmoab->ltog_map));
530   }
531 
532   /* set the DM reference to the vector */
533   PetscCall(VecSetDM(*vec, dm));
534   PetscFunctionReturn(0);
535 }
536 
537 /*  DMVecCreateTagName_Moab_Private
538  *
539  *  Creates a unique tag name that will be shared across processes. If
540  *  pcomm is NULL, then this is a serial vector. A unique tag name
541  *  will be returned in tag_name in either case.
542  *
543  *  The tag names have the format _PETSC_VEC_N where N is some integer.
544  *
545  *  NOTE: The tag_name is allocated in this routine; The user needs to free
546  *        the character array.
547  */
548 #ifdef MOAB_HAVE_MPI
549 PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface, moab::ParallelComm *pcomm, char** tag_name)
550 #else
551 PetscErrorCode DMVecCreateTagName_Moab_Private(moab::Interface *mbiface, char** tag_name)
552 #endif
553 {
554   moab::ErrorCode mberr;
555   PetscInt        n, global_n;
556   moab::Tag indexTag;
557 
558   PetscFunctionBegin;
559   const char*       PVEC_PREFIX      = "__PETSC_VEC_";
560   PetscCall(PetscMalloc1(PETSC_MAX_PATH_LEN, tag_name));
561 
562   moab::EntityHandle rootset = mbiface->get_root_set();
563 
564   /* Check to see if there are any PETSc vectors defined */
565   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
566   mberr = mbiface->tag_get_handle("__PETSC_VECS__", 1, moab::MB_TYPE_INTEGER, indexTag,
567                                   moab::MB_TAG_SPARSE | moab::MB_TAG_CREAT, 0); MBERRNM(mberr);
568   mberr = mbiface->tag_get_data(indexTag, &rootset, 1, &n);
569   if (mberr == moab::MB_TAG_NOT_FOUND) n = 0; /* this is the first temporary vector */
570   else MBERRNM(mberr);
571 
572   /* increment the new value of n */
573   ++n;
574 
575 #ifdef MOAB_HAVE_MPI
576   /* Make sure that n is consistent across all processes */
577   PetscCall(MPIU_Allreduce(&n, &global_n, 1, MPI_INT, MPI_MAX, pcomm->comm()));
578 #else
579   global_n = n;
580 #endif
581 
582   /* Set the new name accordingly and return */
583   PetscCall(PetscSNPrintf(*tag_name, PETSC_MAX_PATH_LEN - 1, "%s_%D", PVEC_PREFIX, global_n));
584   mberr = mbiface->tag_set_data(indexTag, &rootset, 1, (const void*)&global_n); MBERRNM(mberr);
585   PetscFunctionReturn(0);
586 }
587 
588 PETSC_EXTERN PetscErrorCode DMCreateGlobalVector_Moab(DM dm, Vec *gvec)
589 {
590   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
591 
592   PetscFunctionBegin;
593   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
594   PetscValidPointer(gvec,2);
595   PetscCall(DMCreateVector_Moab_Private(dm,NULL,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,gvec));
596   PetscFunctionReturn(0);
597 }
598 
599 PETSC_EXTERN PetscErrorCode DMCreateLocalVector_Moab(DM dm, Vec *lvec)
600 {
601   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
602 
603   PetscFunctionBegin;
604   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
605   PetscValidPointer(lvec,2);
606   PetscCall(DMCreateVector_Moab_Private(dm,NULL,dmmoab->vlocal,PETSC_FALSE,PETSC_TRUE,lvec));
607   PetscFunctionReturn(0);
608 }
609 
610 PetscErrorCode DMVecDuplicate_Moab(Vec x, Vec *y)
611 {
612   DM             dm;
613   PetscContainer  moabdata;
614   Vec_MOAB        *vmoab;
615 
616   PetscFunctionBegin;
617   PetscValidHeaderSpecific(x, VEC_CLASSID, 1);
618   PetscValidPointer(y, 2);
619 
620   /* Get the Vec_MOAB struct for the original vector */
621   PetscCall(PetscObjectQuery((PetscObject)x, "MOABData", (PetscObject*) &moabdata));
622   PetscCall(PetscContainerGetPointer(moabdata, (void**)&vmoab));
623 
624   PetscCall(VecGetDM(x, &dm));
625   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
626 
627   PetscCall(DMCreateVector_Moab_Private(dm,NULL,vmoab->tag_range,vmoab->is_global_vec,PETSC_TRUE,y));
628   PetscCall(VecSetDM(*y, dm));
629   PetscFunctionReturn(0);
630 }
631 
632 PetscErrorCode DMVecUserDestroy_Moab(void *user)
633 {
634   Vec_MOAB        *vmoab = (Vec_MOAB*)user;
635   moab::ErrorCode merr;
636 
637   PetscFunctionBegin;
638   if (vmoab->new_tag && vmoab->tag) {
639     /* Tag was created via a call to VecDuplicate, delete the underlying tag in MOAB */
640     merr = vmoab->mbiface->tag_delete(vmoab->tag); MBERRNM(merr);
641   }
642   delete vmoab->tag_range;
643   vmoab->tag = NULL;
644   vmoab->mbiface = NULL;
645 #ifdef MOAB_HAVE_MPI
646   vmoab->pcomm = NULL;
647 #endif
648   PetscCall(PetscFree(vmoab));
649   PetscFunctionReturn(0);
650 }
651 
652 PETSC_EXTERN PetscErrorCode  DMGlobalToLocalBegin_Moab(DM dm, Vec g, InsertMode mode, Vec l)
653 {
654   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
655 
656   PetscFunctionBegin;
657   PetscCall(VecScatterBegin(dmmoab->ltog_sendrecv, g, l, mode, SCATTER_REVERSE));
658   PetscFunctionReturn(0);
659 }
660 
661 PETSC_EXTERN PetscErrorCode  DMGlobalToLocalEnd_Moab(DM dm, Vec g, InsertMode mode, Vec l)
662 {
663   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
664 
665   PetscFunctionBegin;
666   PetscCall(VecScatterEnd(dmmoab->ltog_sendrecv, g, l, mode, SCATTER_REVERSE));
667   PetscFunctionReturn(0);
668 }
669 
670 PETSC_EXTERN PetscErrorCode  DMLocalToGlobalBegin_Moab(DM dm, Vec l, InsertMode mode, Vec g)
671 {
672   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
673 
674   PetscFunctionBegin;
675   PetscCall(VecScatterBegin(dmmoab->ltog_sendrecv, l, g, mode, SCATTER_FORWARD));
676   PetscFunctionReturn(0);
677 }
678 
679 PETSC_EXTERN PetscErrorCode  DMLocalToGlobalEnd_Moab(DM dm, Vec l, InsertMode mode, Vec g)
680 {
681   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
682 
683   PetscFunctionBegin;
684   PetscCall(VecScatterEnd(dmmoab->ltog_sendrecv, l, g, mode, SCATTER_FORWARD));
685   PetscFunctionReturn(0);
686 }
687