Actual source code: dmmbfield.cxx
1: #include <petsc/private/dmmbimpl.h>
3: #include <petscdmmoab.h>
5: /*@C
6: DMMoabSetFieldVector - Sets the vector reference that represents the solution associated
7: with a particular field component.
9: Not Collective
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)
16: Level: intermediate
18: .seealso: `DMMoabGetFieldName()`, `DMMoabSetGlobalFieldVector()`
19: @*/
20: PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec)
21: {
22: DM_Moab *dmmoab;
23: moab::Tag vtag, ntag;
24: const PetscScalar *varray;
25: PetscScalar *farray;
26: moab::ErrorCode merr;
27: std::string tag_name;
29: PetscFunctionBegin;
31: dmmoab = (DM_Moab *)(dm)->data;
33: PetscCheck(!(ifield < 0) && !(ifield >= dmmoab->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The field %d should be positive and less than %d.", ifield, dmmoab->numFields);
35: /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
36: merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield], 1, moab::MB_TYPE_DOUBLE, ntag, moab::MB_TAG_DENSE | moab::MB_TAG_CREAT);
37: MBERRNM(merr);
39: PetscCall(DMMoabGetVecTag(fvec, &vtag));
41: merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
42: if (!tag_name.length() && merr != moab::MB_SUCCESS) {
43: PetscCall(VecGetArrayRead(fvec, &varray));
44: /* use the entity handle and the Dof index to set the right value */
45: merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void *)varray);
46: MBERRNM(merr);
47: PetscCall(VecRestoreArrayRead(fvec, &varray));
48: } else {
49: PetscCall(PetscMalloc1(dmmoab->nloc, &farray));
50: /* we are using a MOAB Vec - directly copy the tag data to new one */
51: merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void *)farray);
52: MBERRNM(merr);
53: merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void *)farray);
54: MBERRNM(merr);
55: /* make sure the parallel exchange for ghosts are done appropriately */
56: PetscCall(PetscFree(farray));
57: }
58: #ifdef MOAB_HAVE_MPI
59: merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vowned);
60: MBERRNM(merr);
61: #endif
62: PetscFunctionReturn(PETSC_SUCCESS);
63: }
65: /*@C
66: DMMoabSetGlobalFieldVector - Sets the vector reference that represents the global solution associated
67: with all fields (components) managed by DM.
68: A useful utility when updating the DM solution after a solve, to be serialized with the mesh for
69: checkpointing purposes.
71: Not Collective
73: Input Parameters:
74: + dm - the discretization manager object
75: - fvec - the global Vector solution corresponding to all the fields managed by DM
77: Level: intermediate
79: .seealso: `DMMoabGetFieldName()`, `DMMoabSetFieldVector()`
80: @*/
81: PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec)
82: {
83: DM_Moab *dmmoab;
84: moab::Tag vtag, ntag;
85: const PetscScalar *rarray;
86: PetscScalar *varray, *farray;
87: moab::ErrorCode merr;
88: PetscInt i, ifield;
89: std::string tag_name;
90: moab::Range::iterator iter;
92: PetscFunctionBegin;
94: dmmoab = (DM_Moab *)(dm)->data;
96: /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */
97: PetscCall(DMMoabGetVecTag(fvec, &vtag));
98: merr = dmmoab->mbiface->tag_get_name(vtag, tag_name);
99: PetscCall(PetscMalloc1(dmmoab->nloc, &farray));
100: if (!tag_name.length() && merr != moab::MB_SUCCESS) {
101: /* not a MOAB vector - use VecGetSubVector to get the parts as needed */
102: PetscCall(VecGetArrayRead(fvec, &rarray));
103: for (ifield = 0; ifield < dmmoab->numFields; ++ifield) {
104: /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
105: merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield], 1, moab::MB_TYPE_DOUBLE, ntag, moab::MB_TAG_DENSE | moab::MB_TAG_CREAT);
106: MBERRNM(merr);
108: for (i = 0; i < dmmoab->nloc; i++) farray[i] = (dmmoab->bs == 1 ? rarray[ifield * dmmoab->nloc + i] : rarray[i * dmmoab->numFields + ifield]);
110: /* use the entity handle and the Dof index to set the right value */
111: merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void *)farray);
112: MBERRNM(merr);
113: }
114: PetscCall(VecRestoreArrayRead(fvec, &rarray));
115: } else {
116: PetscCall(PetscMalloc1(dmmoab->nloc * dmmoab->numFields, &varray));
118: /* we are using a MOAB Vec - directly copy the tag data to new one */
119: merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void *)varray);
120: MBERRNM(merr);
121: for (ifield = 0; ifield < dmmoab->numFields; ++ifield) {
122: /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */
123: merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield], 1, moab::MB_TYPE_DOUBLE, ntag, moab::MB_TAG_DENSE | moab::MB_TAG_CREAT);
124: MBERRNM(merr);
126: /* we are using a MOAB Vec - directly copy the tag data to new one */
127: for (i = 0; i < dmmoab->nloc; i++) farray[i] = (dmmoab->bs == 1 ? varray[ifield * dmmoab->nloc + i] : varray[i * dmmoab->numFields + ifield]);
129: merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void *)farray);
130: MBERRNM(merr);
132: #ifdef MOAB_HAVE_MPI
133: /* make sure the parallel exchange for ghosts are done appropriately */
134: merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);
135: MBERRNM(merr);
136: #endif
137: }
138: PetscCall(PetscFree(varray));
139: }
140: PetscCall(PetscFree(farray));
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }
144: /*@C
145: DMMoabSetFieldNames - Sets the number of fields and their names to be managed by the DM
147: Not Collective
149: Input Parameters:
150: + dm - the discretization manager object
151: . numFields - the total number of fields
152: - fields - the array containing the names of each field (component); Can be `NULL`.
154: Level: intermediate
156: .seealso: `DMMoabGetFieldName()`, `DMMoabSetFieldName()`
157: @*/
158: PetscErrorCode DMMoabSetFieldNames(DM dm, PetscInt numFields, const char *fields[])
159: {
160: PetscInt i;
161: DM_Moab *dmmoab;
163: PetscFunctionBegin;
165: dmmoab = (DM_Moab *)(dm)->data;
167: /* first deallocate the existing field structure */
168: if (dmmoab->fieldNames) {
169: for (i = 0; i < dmmoab->numFields; i++) PetscCall(PetscFree(dmmoab->fieldNames[i]));
170: PetscCall(PetscFree(dmmoab->fieldNames));
171: }
173: /* now re-allocate and assign field names */
174: dmmoab->numFields = numFields;
175: PetscCall(PetscMalloc1(numFields, &dmmoab->fieldNames));
176: if (fields) {
177: for (i = 0; i < dmmoab->numFields; i++) PetscCall(PetscStrallocpy(fields[i], (char **)&dmmoab->fieldNames[i]));
178: }
179: PetscCall(DMSetNumFields(dm, numFields));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: /*@C
184: DMMoabGetFieldName - Gets the names of individual field components in multicomponent
185: vectors associated with a DMDA.
187: Not Collective
189: Input Parameters:
190: + dm - the discretization manager object
191: - field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the
192: number of degrees of freedom per node within the DMMoab
194: Output Parameter:
195: . fieldName - the name of the field (component)
197: Level: intermediate
199: .seealso: `DMMoabSetFieldName()`, `DMMoabSetFields()`
200: @*/
201: PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName)
202: {
203: DM_Moab *dmmoab;
205: PetscFunctionBegin;
207: dmmoab = (DM_Moab *)(dm)->data;
208: PetscCheck(!(field < 0) && !(field >= dmmoab->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields);
210: *fieldName = dmmoab->fieldNames[field];
211: PetscFunctionReturn(PETSC_SUCCESS);
212: }
214: /*@C
215: DMMoabSetFieldName - Sets the name of a field (component) managed by the DM
217: Not Collective
219: Input Parameters:
220: + dm - the discretization manager object
221: . field - the field number
222: - fieldName - the field (component) name
224: Level: intermediate
226: Notes:
227: Can only be called after DMMoabSetFields supplied with correct numFields
229: .seealso: `DMMoabGetFieldName()`, `DMMoabSetFields()`
230: @*/
231: PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName)
232: {
233: DM_Moab *dmmoab;
235: PetscFunctionBegin;
237: PetscAssertPointer(fieldName, 3);
239: dmmoab = (DM_Moab *)(dm)->data;
240: PetscCheck(!(field < 0) && !(field >= dmmoab->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "DM field %d should be in [%d, %d)", field, 0, dmmoab->numFields);
242: if (dmmoab->fieldNames[field]) PetscCall(PetscFree(dmmoab->fieldNames[field]));
243: PetscCall(PetscStrallocpy(fieldName, (char **)&dmmoab->fieldNames[field]));
244: PetscFunctionReturn(PETSC_SUCCESS);
245: }
247: /*@C
248: DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a
249: particular MOAB EntityHandle.
251: Not Collective
253: Input Parameters:
254: + dm - the discretization manager object
255: . point - the MOAB EntityHandle container which holds the field degree-of-freedom values
256: - field - the field (component) index
258: Output Parameter:
259: . dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat)
261: Level: beginner
263: .seealso: `DMMoabGetFieldDofs()`, `DMMoabGetFieldDofsLocal()`
264: @*/
265: PetscErrorCode DMMoabGetFieldDof(DM dm, moab::EntityHandle point, PetscInt field, PetscInt *dof)
266: {
267: DM_Moab *dmmoab;
269: PetscFunctionBegin;
271: dmmoab = (DM_Moab *)(dm)->data;
273: *dof = (dmmoab->bs == 1 ? dmmoab->gidmap[dmmoab->mbiface->id_from_handle(point) - dmmoab->seqstart] + field * dmmoab->n : dmmoab->gidmap[dmmoab->mbiface->id_from_handle(point) - dmmoab->seqstart] * dmmoab->numFields + field);
274: PetscFunctionReturn(PETSC_SUCCESS);
275: }
277: /*@C
278: DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an
279: array of MOAB EntityHandles.
281: Not Collective
283: Input Parameters:
284: + dm - the discretization manager object
285: . npoints - the total number of Entities in the points array
286: . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
287: - field - the field (component) index
289: Output Parameter:
290: . dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
292: Level: intermediate
294: .seealso: `DMMoabGetFieldDof()`, `DMMoabGetFieldDofsLocal()`
295: @*/
296: PetscErrorCode DMMoabGetFieldDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof)
297: {
298: PetscInt i;
299: DM_Moab *dmmoab;
301: PetscFunctionBegin;
303: PetscAssertPointer(points, 3);
304: dmmoab = (DM_Moab *)(dm)->data;
306: if (!dof) PetscCall(PetscMalloc1(npoints, &dof));
308: /* compute the DOF based on local blocking in the fields */
309: /* We also assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
310: /* TODO: eliminate the limitation using PetscSection to manage DOFs */
311: for (i = 0; i < npoints; ++i)
312: dof[i] = (dmmoab->bs == 1 ? dmmoab->gidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] + field * dmmoab->n : dmmoab->gidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] * dmmoab->numFields + field);
313: PetscFunctionReturn(PETSC_SUCCESS);
314: }
316: /*@C
317: DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an
318: array of MOAB EntityHandles.
320: Not Collective
322: Input Parameters:
323: + dm - the discretization manager object
324: . npoints - the total number of Entities in the points array
325: . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
326: - field - the field (component) index
328: Output Parameter:
329: . dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
331: Level: intermediate
333: .seealso: `DMMoabGetFieldDof()`, `DMMoabGetFieldDofs()`
334: @*/
335: PetscErrorCode DMMoabGetFieldDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt field, PetscInt *dof)
336: {
337: PetscInt i;
338: DM_Moab *dmmoab;
340: PetscFunctionBegin;
342: PetscAssertPointer(points, 3);
343: dmmoab = (DM_Moab *)(dm)->data;
345: if (!dof) PetscCall(PetscMalloc1(npoints, &dof));
347: /* compute the DOF based on local blocking in the fields */
348: /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
349: /* TODO: eliminate the limitation using PetscSection to manage DOFs */
350: for (i = 0; i < npoints; ++i) {
351: dof[i] = (dmmoab->bs > 1 ? dmmoab->lidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] * dmmoab->numFields + field : dmmoab->lidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] + field * dmmoab->n);
352: }
353: PetscFunctionReturn(PETSC_SUCCESS);
354: }
356: /*@C
357: DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an
358: array of MOAB EntityHandles.
360: Not Collective
362: Input Parameters:
363: + dm - the discretization manager object
364: . npoints - the total number of Entities in the points array
365: - points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
367: Output Parameter:
368: . dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
370: Level: intermediate
372: .seealso: `DMMoabGetFieldDofs()`, `DMMoabGetDofsLocal()`, `DMMoabGetDofsBlocked()`
373: @*/
374: PetscErrorCode DMMoabGetDofs(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof)
375: {
376: PetscInt i, field, offset;
377: DM_Moab *dmmoab;
379: PetscFunctionBegin;
381: PetscAssertPointer(points, 3);
382: dmmoab = (DM_Moab *)(dm)->data;
384: if (!dof) PetscCall(PetscMalloc1(dmmoab->numFields * npoints, &dof));
386: /* compute the DOF based on local blocking in the fields */
387: /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
388: /* TODO: eliminate the limitation using PetscSection to manage DOFs */
389: for (field = 0; field < dmmoab->numFields; ++field) {
390: offset = field * dmmoab->n;
391: for (i = 0; i < npoints; ++i)
392: dof[i * dmmoab->numFields + field] = (dmmoab->bs > 1 ? dmmoab->gidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] * dmmoab->numFields + field : dmmoab->gidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] + offset);
393: }
394: PetscFunctionReturn(PETSC_SUCCESS);
395: }
397: /*@C
398: DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an
399: array of MOAB EntityHandles.
401: Not Collective
403: Input Parameters:
404: + dm - the discretization manager object
405: . npoints - the total number of Entities in the points array
406: - points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
408: Output Parameter:
409: . dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
411: Level: intermediate
413: .seealso: `DMMoabGetFieldDofs()`, `DMMoabGetDofs()`, `DMMoabGetDofsBlocked()`
414: @*/
415: PetscErrorCode DMMoabGetDofsLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof)
416: {
417: PetscInt i, field, offset;
418: DM_Moab *dmmoab;
420: PetscFunctionBegin;
422: PetscAssertPointer(points, 3);
423: dmmoab = (DM_Moab *)(dm)->data;
425: if (!dof) PetscCall(PetscMalloc1(dmmoab->numFields * npoints, &dof));
427: /* compute the DOF based on local blocking in the fields */
428: /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
429: /* TODO: eliminate the limitation using PetscSection to manage DOFs */
430: for (field = 0; field < dmmoab->numFields; ++field) {
431: offset = field * dmmoab->n;
432: for (i = 0; i < npoints; ++i)
433: dof[i * dmmoab->numFields + field] = (dmmoab->bs > 1 ? dmmoab->lidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] * dmmoab->numFields + field : dmmoab->lidmap[dmmoab->mbiface->id_from_handle(points[i]) - dmmoab->seqstart] + offset);
434: }
435: PetscFunctionReturn(PETSC_SUCCESS);
436: }
438: /*@C
439: DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
440: array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation
441: of element residuals and assembly of the discrete systems when all fields are co-located.
443: Not Collective
445: Input Parameters:
446: + dm - the discretization manager object
447: . npoints - the total number of Entities in the points array
448: - points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
450: Output Parameter:
451: . dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat)
453: Level: intermediate
455: .seealso: `DMMoabGetDofsLocal()`, `DMMoabGetDofs()`, `DMMoabGetDofsBlockedLocal()`
456: @*/
457: PetscErrorCode DMMoabGetDofsBlocked(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof)
458: {
459: PetscInt i;
460: DM_Moab *dmmoab;
462: PetscFunctionBegin;
464: PetscAssertPointer(points, 3);
465: dmmoab = (DM_Moab *)(dm)->data;
467: if (!dof) PetscCall(PetscMalloc1(npoints, &dof));
469: for (i = 0; i < npoints; ++i) dof[i] = dmmoab->gidmap[(PetscInt)points[i] - dmmoab->seqstart];
470: PetscFunctionReturn(PETSC_SUCCESS);
471: }
473: /*@C
474: DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
475: array of MOAB EntityHandles. It is useful when performing local Blocked(Get/Set) methods in computation
476: of element residuals and assembly of the discrete systems when all fields are co-located.
478: Not Collective
480: Input Parameters:
481: + dm - the discretization manager object
482: . npoints - the total number of Entities in the points array
483: - points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
485: Output Parameter:
486: . dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat)
488: Level: intermediate
490: .seealso: `DMMoabGetDofsLocal()`, `DMMoabGetDofs()`
491: @*/
492: PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm, PetscInt npoints, const moab::EntityHandle *points, PetscInt *dof)
493: {
494: PetscInt i;
495: DM_Moab *dmmoab;
497: PetscFunctionBegin;
499: PetscAssertPointer(points, 3);
500: dmmoab = (DM_Moab *)(dm)->data;
502: if (!dof) PetscCall(PetscMalloc1(npoints, &dof));
504: for (i = 0; i < npoints; ++i) dof[i] = dmmoab->lidmap[(PetscInt)points[i] - dmmoab->seqstart];
505: PetscFunctionReturn(PETSC_SUCCESS);
506: }
508: /*@C
509: DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
510: array of locally owned MOAB mesh vertices.
512: Not Collective
514: Input Parameters:
515: . dm - the discretization manager object
517: Output Parameter:
518: . dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
520: Level: intermediate
522: Note:
523: It's utility is when performing Finite-Difference type calculations where vertex traversal is
524: faster than element-wise assembly that is typically done in FEM calculations.
526: .seealso: `DMMoabGetVertexDofsBlockedLocal()`, `DMMoabGetDofsBlocked()`, `DMMoabGetDofsBlockedLocal()`
527: @*/
528: PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm, PetscInt **dof)
529: {
530: DM_Moab *dmmoab;
532: PetscFunctionBegin;
534: dmmoab = (DM_Moab *)(dm)->data;
536: *dof = dmmoab->gidmap;
537: PetscFunctionReturn(PETSC_SUCCESS);
538: }
540: /*@C
541: DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
542: array of locally owned MOAB mesh vertices.
544: Not Collective
546: Input Parameters:
547: . dm - the discretization manager object
549: Output Parameter:
550: . dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
552: Level: intermediate
554: Note:
555: It's utility is when performing Finite-Difference type calculations where vertex traversal is
556: faster than element-wise assembly that is typically done in FEM calculations.
558: .seealso: `DMMoabGetVertexDofsBlocked()`, `DMMoabGetDofsBlocked()`, `DMMoabGetDofsBlockedLocal()`
559: @*/
560: PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm, PetscInt **dof)
561: {
562: DM_Moab *dmmoab;
564: PetscFunctionBegin;
566: PetscAssertPointer(dof, 2);
567: dmmoab = (DM_Moab *)(dm)->data;
569: *dof = dmmoab->lidmap;
570: PetscFunctionReturn(PETSC_SUCCESS);
571: }