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