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