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