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