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