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