xref: /petsc/src/dm/impls/moab/dmmbfield.cxx (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2 
3 #include <petscdmmoab.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMMoabSetFieldVector"
7 /*@C
8   DMMoabSetFieldVector - Sets the vector reference that represents the solution associated
9   with a particular field component.
10 
11   Not Collective
12 
13   Input Parameters:
14 + dm     - the discretization manager object
15 . ifield - the index of the field as set before via DMMoabSetFieldName.
16 . fvec - the Vector solution corresponding to the field (component)
17 
18   Level: intermediate
19 
20 .keywords: discretization manager, set, component solution
21 
22 .seealso: DMMoabGetFieldName(), DMMoabSetGlobalFieldVector()
23 @*/
24 PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
25 {
26   DM_Moab        *dmmoab;
27   moab::Tag     vtag,ntag;
28   const PetscScalar *varray;
29   PetscScalar *farray;
30   moab::ErrorCode merr;
31   PetscErrorCode  ierr;
32   std::string tag_name;
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36   dmmoab = (DM_Moab*)(dm)->data;
37 
38   if ((ifield < 0) || (ifield >= dmmoab->numFields)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The field %d should be positive and less than %d.", ifield, dmmoab->numFields);
39 
40   /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
41   merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
42                                           moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
43 
44   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
45 
46   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
47   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
48     ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr);
49     /* use the entity handle and the Dof index to set the right value */
50     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr);
51     ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr);
52   }
53   else {
54     ierr = PetscMalloc1(dmmoab->nloc,&farray);CHKERRQ(ierr);
55     /* we are using a MOAB Vec - directly copy the tag data to new one */
56     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr);
57     merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
58     /* make sure the parallel exchange for ghosts are done appropriately */
59     ierr = PetscFree(farray);CHKERRQ(ierr);
60   }
61   merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vowned);MBERRNM(merr);
62   PetscFunctionReturn(0);
63 }
64 
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "DMMoabSetGlobalFieldVector"
68 /*@C
69   DMMoabSetGlobalFieldVector - Sets the vector reference that represents the global solution associated
70   with all fields (components) managed by DM.
71   A useful utility when updating the DM solution after a solve, to be serialized with the mesh for
72   checkpointing purposes.
73 
74   Not Collective
75 
76   Input Parameters:
77 + dm     - the discretization manager object
78 . fvec - the global Vector solution corresponding to all the fields managed by DM
79 
80   Level: intermediate
81 
82 .keywords: discretization manager, set, component solution
83 
84 .seealso: DMMoabGetFieldName(), DMMoabSetFieldVector()
85 @*/
86 PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec)
87 {
88   DM_Moab        *dmmoab;
89   moab::Tag     vtag,ntag;
90   const PetscScalar   *rarray;
91   PetscScalar   *varray,*farray;
92   moab::ErrorCode merr;
93   PetscErrorCode  ierr;
94   PetscInt i,ifield;
95   std::string tag_name;
96   moab::Range::iterator iter;
97 
98   PetscFunctionBegin;
99   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
100   dmmoab = (DM_Moab*)(dm)->data;
101 
102   /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */
103   ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr);
104   merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
105   ierr = PetscMalloc1(dmmoab->nloc,&farray);CHKERRQ(ierr);
106   if (!tag_name.length() && merr !=moab::MB_SUCCESS) {
107     /* not a MOAB vector - use VecGetSubVector to get the parts as needed */
108 
109     ierr = VecGetArrayRead(fvec,&rarray);CHKERRQ(ierr);
110     for (ifield=0; ifield<dmmoab->numFields; ++ifield) {
111 
112       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
113       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
114                                             moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
115 
116       for(i=0;i<dmmoab->nloc;i++) {
117         if (dmmoab->bs == 1)
118           farray[i]=rarray[ifield*dmmoab->nloc+i];
119         else
120           farray[i]=rarray[i*dmmoab->numFields+ifield];
121       }
122 
123       /* use the entity handle and the Dof index to set the right value */
124       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
125     }
126     ierr = VecRestoreArrayRead(fvec,&rarray);CHKERRQ(ierr);
127   }
128   else {
129     ierr = PetscMalloc1(dmmoab->nloc*dmmoab->numFields,&varray);CHKERRQ(ierr);
130 
131     /* we are using a MOAB Vec - directly copy the tag data to new one */
132     merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr);
133     for (ifield=0; ifield<dmmoab->numFields; ++ifield) {
134 
135       /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
136       merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag,
137                                             moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr);
138 
139       /* we are using a MOAB Vec - directly copy the tag data to new one */
140       for(i=0; i < dmmoab->nloc; i++) {
141         if (dmmoab->bs == 1)
142           farray[i]=varray[ifield*dmmoab->nloc+i];
143         else
144           farray[i]=varray[i*dmmoab->numFields+ifield];
145       }
146 
147       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
148       /* make sure the parallel exchange for ghosts are done appropriately */
149       merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
150     }
151     ierr = PetscFree(varray);CHKERRQ(ierr);
152   }
153   ierr = PetscFree(farray);CHKERRQ(ierr);
154   PetscFunctionReturn(0);
155 }
156 
157 
158 #undef __FUNCT__
159 #define __FUNCT__ "DMMoabSetFieldNames"
160 /*@C
161   DMMoabSetFieldNames - Sets the number of fields and their names to be managed by the DM
162 
163   Not Collective
164 
165   Input Parameters:
166 + dm     - the discretization manager object
167 . numFields - the total number of fields
168 . fields - the array containing the names of each field (component); Can be NULL.
169 
170   Level: intermediate
171 
172 .keywords: discretization manager, set, component name
173 
174 .seealso: DMMoabGetFieldName(), DMMoabSetFieldName()
175 @*/
176 PetscErrorCode DMMoabSetFieldNames(DM dm,PetscInt numFields,const char* fields[])
177 {
178   PetscErrorCode ierr;
179   PetscInt       i;
180   DM_Moab        *dmmoab;
181 
182   PetscFunctionBegin;
183   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
184   dmmoab = (DM_Moab*)(dm)->data;
185 
186   /* first deallocate the existing field structure */
187   if (dmmoab->fieldNames) {
188     for(i=0; i<dmmoab->numFields; i++) {
189       ierr = PetscFree(dmmoab->fieldNames[i]);CHKERRQ(ierr);
190     }
191     ierr = PetscFree(dmmoab->fieldNames);CHKERRQ(ierr);
192   }
193 
194   /* now re-allocate and assign field names  */
195   dmmoab->numFields = numFields;
196   ierr = PetscMalloc(numFields*sizeof(char*),&dmmoab->fieldNames);CHKERRQ(ierr);
197   if (fields) {
198     for(i=0; i<dmmoab->numFields; i++) {
199       ierr = PetscStrallocpy(fields[i], (char**) &dmmoab->fieldNames[i]);CHKERRQ(ierr);
200     }
201   }
202   PetscFunctionReturn(0);
203 }
204 
205 
206 #undef __FUNCT__
207 #define __FUNCT__ "DMMoabGetFieldName"
208 /*@C
209   DMMoabGetFieldName - Gets the names of individual field components in multicomponent
210   vectors associated with a DMDA.
211 
212   Not Collective
213 
214   Input Parameter:
215 + dm     - the discretization manager object
216 . field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the
217         number of degrees of freedom per node within the DMMoab
218 
219   Output Parameter:
220 . fieldName - the name of the field (component)
221 
222   Level: intermediate
223 
224 .keywords: discretization manager, get, component name
225 
226 .seealso: DMMoabSetFieldName(), DMMoabSetFields()
227 @*/
228 PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName)
229 {
230   DM_Moab        *dmmoab;
231 
232   PetscFunctionBegin;
233   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
234   dmmoab = (DM_Moab*)(dm)->data;
235   if ((field < 0) || (field >= dmmoab->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields);
236 
237   *fieldName = dmmoab->fieldNames[field];
238   PetscFunctionReturn(0);
239 }
240 
241 
242 #undef __FUNCT__
243 #define __FUNCT__ "DMMoabSetFieldName"
244 /*@C
245   DMMoabSetFieldName - Sets the name of a field (component) managed by the DM
246 
247   Not Collective
248 
249   Input Parameters:
250 + dm     - the discretization manager object
251 . field - the field number
252 . fieldName - the field (component) name
253 
254   Level: intermediate
255   Notes: Can only be called after DMMoabSetFields supplied with correct numFields
256 
257 .keywords: discretization manager, set, component name
258 
259 .seealso: DMMoabGetFieldName(), DMMoabSetFields()
260 @*/
261 PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName)
262 {
263   PetscErrorCode ierr;
264   DM_Moab        *dmmoab;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
268   PetscValidCharPointer(fieldName,3);
269 
270   dmmoab = (DM_Moab*)(dm)->data;
271   if ((field < 0) || (field >= dmmoab->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields);
272 
273   if (dmmoab->fieldNames[field]) {
274     ierr = PetscFree(dmmoab->fieldNames[field]);CHKERRQ(ierr);
275   }
276   ierr = PetscStrallocpy(fieldName, (char**) &dmmoab->fieldNames[field]);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 
281 #undef __FUNCT__
282 #define __FUNCT__ "DMMoabGetFieldDof"
283 /*@C
284   DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a
285   particular MOAB EntityHandle.
286 
287   Not Collective
288 
289   Input Parameters:
290 + dm     - the discretization manager object
291 . point - the MOAB EntityHandle container which holds the field degree-of-freedom values
292 . field - the field (component) index
293 
294   Output Parameter:
295 + dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat)
296 
297   Level: beginner
298 
299 .keywords: discretization manager, get, global degree of freedom
300 
301 .seealso: DMMoabGetFieldDofs(), DMMoabGetFieldDofsLocal()
302 @*/
303 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
304 {
305   DM_Moab        *dmmoab;
306 
307   PetscFunctionBegin;
308   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
309   dmmoab = (DM_Moab*)(dm)->data;
310 
311   *dof=dmmoab->gidmap[(PetscInt)point];
312   PetscFunctionReturn(0);
313 }
314 
315 
316 #undef __FUNCT__
317 #define __FUNCT__ "DMMoabGetFieldDofs"
318 /*@C
319   DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an
320   array of MOAB EntityHandles.
321 
322   Not Collective
323 
324   Input Parameters:
325 + dm     - the discretization manager object
326 . npoints - the total number of Entities in the points array
327 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
328 . field - the field (component) index
329 
330   Output Parameter:
331 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
332 
333   Level: intermediate
334 
335 .keywords: discretization manager, get, global degrees of freedom
336 
337 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofsLocal()
338 @*/
339 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
340 {
341   PetscInt        i;
342   PetscErrorCode  ierr;
343   DM_Moab        *dmmoab;
344 
345   PetscFunctionBegin;
346   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
347   PetscValidPointer(points,3);
348   dmmoab = (DM_Moab*)(dm)->data;
349 
350   if (!dof) {
351     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
352   }
353 
354   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
355     for (i=0; i<npoints; ++i)
356       dof[i] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
357   }
358   else {
359     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
360     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
361     for (i=0; i<npoints; ++i)
362       dof[i] = dmmoab->gidmap[(PetscInt)points[i]]+field*dmmoab->n;
363   }
364   PetscFunctionReturn(0);
365 }
366 
367 
368 #undef __FUNCT__
369 #define __FUNCT__ "DMMoabGetFieldDofsLocal"
370 /*@C
371   DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an
372   array of MOAB EntityHandles.
373 
374   Not Collective
375 
376   Input Parameters:
377 + dm     - the discretization manager object
378 . npoints - the total number of Entities in the points array
379 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
380 . field - the field (component) index
381 
382   Output Parameter:
383 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
384 
385   Level: intermediate
386 
387 .keywords: discretization manager, get, local degrees of freedom
388 
389 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofs()
390 @*/
391 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
392 {
393   PetscInt i,offset;
394   PetscErrorCode  ierr;
395   DM_Moab        *dmmoab;
396 
397   PetscFunctionBegin;
398   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
399   PetscValidPointer(points,3);
400   dmmoab = (DM_Moab*)(dm)->data;
401 
402   if (!dof) {
403     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
404   }
405 
406   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
407     for (i=0; i<npoints; ++i)
408       dof[i] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
409   }
410   else {
411     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
412     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
413     offset = field*dmmoab->n;
414     for (i=0; i<npoints; ++i)
415       dof[i] = dmmoab->lidmap[(PetscInt)points[i]]+offset;
416   }
417   PetscFunctionReturn(0);
418 }
419 
420 
421 #undef __FUNCT__
422 #define __FUNCT__ "DMMoabGetDofs"
423 /*@C
424   DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an
425   array of MOAB EntityHandles.
426 
427   Not Collective
428 
429   Input Parameters:
430 + dm     - the discretization manager object
431 . npoints - the total number of Entities in the points array
432 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
433 
434   Output Parameter:
435 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
436 
437   Level: intermediate
438 
439 .keywords: discretization manager, get, global degrees of freedom
440 
441 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofsLocal(), DMMoabGetDofsBlocked()
442 @*/
443 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
444 {
445   PetscInt        i,field,offset;
446   PetscErrorCode  ierr;
447   DM_Moab        *dmmoab;
448 
449   PetscFunctionBegin;
450   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
451   PetscValidPointer(points,3);
452   dmmoab = (DM_Moab*)(dm)->data;
453 
454   if (!dof) {
455     ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr);
456   }
457 
458   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
459     for (field=0; field<dmmoab->numFields; ++field) {
460       for (i=0; i<npoints; ++i)
461         dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
462     }
463   }
464   else {
465     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
466     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
467     for (field=0; field<dmmoab->numFields; ++field) {
468       offset = field*dmmoab->n;
469       for (i=0; i<npoints; ++i)
470         dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]+offset;
471     }
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "DMMoabGetDofsLocal"
479 /*@C
480   DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an
481   array of MOAB EntityHandles.
482 
483   Not Collective
484 
485   Input Parameters:
486 + dm     - the discretization manager object
487 . npoints - the total number of Entities in the points array
488 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
489 
490   Output Parameter:
491 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
492 
493   Level: intermediate
494 
495 .keywords: discretization manager, get, global degrees of freedom
496 
497 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofs(), DMMoabGetDofsBlocked()
498 @*/
499 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
500 {
501   PetscInt        i,field,offset;
502   PetscErrorCode  ierr;
503   DM_Moab        *dmmoab;
504 
505   PetscFunctionBegin;
506   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
507   PetscValidPointer(points,3);
508   dmmoab = (DM_Moab*)(dm)->data;
509 
510   if (!dof) {
511     ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr);
512   }
513 
514   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
515     for (field=0; field<dmmoab->numFields; ++field) {
516       for (i=0; i<npoints; ++i)
517         dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
518     }
519   }
520   else {
521     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
522     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
523     for (field=0; field<dmmoab->numFields; ++field) {
524       offset = field*dmmoab->n;
525       for (i=0; i<npoints; ++i)
526         dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]+offset;
527     }
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "DMMoabGetDofsBlocked"
535 /*@C
536   DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
537   array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation
538   of element residuals and assembly of the discrete systems when all fields are co-located.
539 
540   Not Collective
541 
542   Input Parameters:
543 + dm     - the discretization manager object
544 . npoints - the total number of Entities in the points array
545 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
546 
547   Output Parameter:
548 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat)
549 
550   Level: intermediate
551 
552 .keywords: discretization manager, get, global degrees of freedom
553 
554 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal()
555 @*/
556 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
557 {
558   PetscInt        i;
559   DM_Moab        *dmmoab;
560   PetscErrorCode  ierr;
561 
562   PetscFunctionBegin;
563   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
564   PetscValidPointer(points,2);
565   dmmoab = (DM_Moab*)(dm)->data;
566 
567   if (!dof) {
568     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
569   }
570 
571   for (i=0; i<npoints; ++i) {
572     dof[i]=dmmoab->gidmap[(PetscInt)points[i]];
573   }
574   PetscFunctionReturn(0);
575 }
576 
577 
578 #undef __FUNCT__
579 #define __FUNCT__ "DMMoabGetDofsBlockedLocal"
580 /*@C
581   DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
582   array of MOAB EntityHandles. It is useful when performing local Blocked(Get/Set) methods in computation
583   of element residuals and assembly of the discrete systems when all fields are co-located.
584 
585   Not Collective
586 
587   Input Parameters:
588 + dm     - the discretization manager object
589 . npoints - the total number of Entities in the points array
590 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
591 
592   Output Parameter:
593 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat)
594 
595   Level: intermediate
596 
597 .keywords: discretization manager, get, global degrees of freedom
598 
599 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal()
600 @*/
601 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
602 {
603   PetscInt        i;
604   DM_Moab        *dmmoab;
605   PetscErrorCode  ierr;
606 
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
609   PetscValidPointer(points,2);
610   dmmoab = (DM_Moab*)(dm)->data;
611 
612   if (!dof) {
613     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
614   }
615 
616   for (i=0; i<npoints; ++i)
617     dof[i] = dmmoab->lidmap[(PetscInt)points[i]];
618   PetscFunctionReturn(0);
619 }
620 
621 
622 #undef __FUNCT__
623 #define __FUNCT__ "DMMoabGetVertexDofsBlocked"
624 /*@C
625   DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
626   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
627   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
628 
629   Not Collective
630 
631   Input Parameters:
632 + dm     - the discretization manager object
633 
634   Output Parameter:
635 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
636 
637   Level: intermediate
638 
639 .keywords: discretization manager, get, blocked degrees of freedom
640 
641 .seealso: DMMoabGetVertexDofsBlockedLocal(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
642 @*/
643 PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof)
644 {
645   DM_Moab        *dmmoab;
646 
647   PetscFunctionBegin;
648   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
649   dmmoab = (DM_Moab*)(dm)->data;
650 
651   *dof = dmmoab->gidmap;
652   PetscFunctionReturn(0);
653 }
654 
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "DMMoabGetVertexDofsBlockedLocal"
658 /*@C
659   DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
660   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
661   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
662 
663   Not Collective
664 
665   Input Parameters:
666 + dm     - the discretization manager object
667 
668   Output Parameter:
669 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
670 
671   Level: intermediate
672 
673 .keywords: discretization manager, get, blocked degrees of freedom
674 
675 .seealso: DMMoabGetVertexDofsBlocked(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
676 @*/
677 PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof)
678 {
679   DM_Moab        *dmmoab;
680 
681   PetscFunctionBegin;
682   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
683   PetscValidPointer(dof,2);
684   dmmoab = (DM_Moab*)(dm)->data;
685 
686   *dof = dmmoab->lidmap;
687   PetscFunctionReturn(0);
688 }
689 
690