xref: /petsc/src/dm/impls/moab/dmmbfield.cxx (revision 9daf19fd0d22011a060902ed53673b48ac0c4e7d)
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   PetscFunctionReturn(0);
195 }
196 
197 
198 /*@C
199   DMMoabGetFieldName - Gets the names of individual field components in multicomponent
200   vectors associated with a DMDA.
201 
202   Not Collective
203 
204   Input Parameter:
205 + dm     - the discretization manager object
206 . field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the
207         number of degrees of freedom per node within the DMMoab
208 
209   Output Parameter:
210 . fieldName - the name of the field (component)
211 
212   Level: intermediate
213 
214 .keywords: discretization manager, get, component name
215 
216 .seealso: DMMoabSetFieldName(), DMMoabSetFields()
217 @*/
218 PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName)
219 {
220   DM_Moab        *dmmoab;
221 
222   PetscFunctionBegin;
223   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
224   dmmoab = (DM_Moab*)(dm)->data;
225   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);
226 
227   *fieldName = dmmoab->fieldNames[field];
228   PetscFunctionReturn(0);
229 }
230 
231 
232 /*@C
233   DMMoabSetFieldName - Sets the name of a field (component) managed by the DM
234 
235   Not Collective
236 
237   Input Parameters:
238 + dm     - the discretization manager object
239 . field - the field number
240 . fieldName - the field (component) name
241 
242   Level: intermediate
243   Notes: Can only be called after DMMoabSetFields supplied with correct numFields
244 
245 .keywords: discretization manager, set, component name
246 
247 .seealso: DMMoabGetFieldName(), DMMoabSetFields()
248 @*/
249 PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName)
250 {
251   PetscErrorCode ierr;
252   DM_Moab        *dmmoab;
253 
254   PetscFunctionBegin;
255   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
256   PetscValidCharPointer(fieldName,3);
257 
258   dmmoab = (DM_Moab*)(dm)->data;
259   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);
260 
261   if (dmmoab->fieldNames[field]) {
262     ierr = PetscFree(dmmoab->fieldNames[field]);CHKERRQ(ierr);
263   }
264   ierr = PetscStrallocpy(fieldName, (char**) &dmmoab->fieldNames[field]);CHKERRQ(ierr);
265   PetscFunctionReturn(0);
266 }
267 
268 
269 /*@C
270   DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a
271   particular MOAB EntityHandle.
272 
273   Not Collective
274 
275   Input Parameters:
276 + dm     - the discretization manager object
277 . point - the MOAB EntityHandle container which holds the field degree-of-freedom values
278 . field - the field (component) index
279 
280   Output Parameter:
281 + dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat)
282 
283   Level: beginner
284 
285 .keywords: discretization manager, get, global degree of freedom
286 
287 .seealso: DMMoabGetFieldDofs(), DMMoabGetFieldDofsLocal()
288 @*/
289 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
290 {
291   DM_Moab        *dmmoab;
292 
293   PetscFunctionBegin;
294   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
295   dmmoab = (DM_Moab*)(dm)->data;
296 
297   *dof=dmmoab->gidmap[((PetscInt)point-dmmoab->seqstart)]*dmmoab->numFields+field;
298   PetscFunctionReturn(0);
299 }
300 
301 
302 /*@C
303   DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an
304   array of MOAB EntityHandles.
305 
306   Not Collective
307 
308   Input Parameters:
309 + dm     - the discretization manager object
310 . npoints - the total number of Entities in the points array
311 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
312 . field - the field (component) index
313 
314   Output Parameter:
315 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
316 
317   Level: intermediate
318 
319 .keywords: discretization manager, get, global degrees of freedom
320 
321 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofsLocal()
322 @*/
323 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
324 {
325   PetscInt        i;
326   PetscErrorCode  ierr;
327   DM_Moab        *dmmoab;
328 
329   PetscFunctionBegin;
330   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
331   PetscValidPointer(points,3);
332   dmmoab = (DM_Moab*)(dm)->data;
333 
334   if (!dof) {
335     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
336   }
337 
338   /* compute the DOF based on local blocking in the fields */
339   /* We also assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
340   /* TODO: eliminate the limitation using PetscSection to manage DOFs */
341   for (i=0; i<npoints; ++i)
342     dof[i] = (dmmoab->bs == 1 ? dmmoab->gidmap[(PetscInt)points[i]-dmmoab->seqstart]+field*dmmoab->n :
343                                 dmmoab->gidmap[(PetscInt)points[i]-dmmoab->seqstart]*dmmoab->numFields+field);
344   PetscFunctionReturn(0);
345 }
346 
347 
348 /*@C
349   DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an
350   array of MOAB EntityHandles.
351 
352   Not Collective
353 
354   Input Parameters:
355 + dm     - the discretization manager object
356 . npoints - the total number of Entities in the points array
357 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
358 . field - the field (component) index
359 
360   Output Parameter:
361 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
362 
363   Level: intermediate
364 
365 .keywords: discretization manager, get, local degrees of freedom
366 
367 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofs()
368 @*/
369 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
370 {
371   PetscInt i;
372   PetscErrorCode  ierr;
373   DM_Moab        *dmmoab;
374 
375   PetscFunctionBegin;
376   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
377   PetscValidPointer(points,3);
378   dmmoab = (DM_Moab*)(dm)->data;
379 
380   if (!dof) {
381     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
382   }
383 
384   /* compute the DOF based on local blocking in the fields */
385   /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
386   /* TODO: eliminate the limitation using PetscSection to manage DOFs */
387   for (i=0; i<npoints; ++i) {
388     dof[i] = (dmmoab->bs > 1 ? dmmoab->lidmap[(PetscInt)points[i]-dmmoab->seqstart]*dmmoab->numFields+field :
389                                dmmoab->lidmap[(PetscInt)points[i]-dmmoab->seqstart]+field*dmmoab->n);
390   }
391   PetscFunctionReturn(0);
392 }
393 
394 
395 /*@C
396   DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an
397   array of MOAB EntityHandles.
398 
399   Not Collective
400 
401   Input Parameters:
402 + dm     - the discretization manager object
403 . npoints - the total number of Entities in the points array
404 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
405 
406   Output Parameter:
407 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
408 
409   Level: intermediate
410 
411 .keywords: discretization manager, get, global degrees of freedom
412 
413 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofsLocal(), DMMoabGetDofsBlocked()
414 @*/
415 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
416 {
417   PetscInt        i,field,offset;
418   PetscErrorCode  ierr;
419   DM_Moab        *dmmoab;
420 
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
423   PetscValidPointer(points,3);
424   dmmoab = (DM_Moab*)(dm)->data;
425 
426   if (!dof) {
427     ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr);
428   }
429 
430   /* compute the DOF based on local blocking in the fields */
431   /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
432   /* TODO: eliminate the limitation using PetscSection to manage DOFs */
433   for (field=0; field<dmmoab->numFields; ++field) {
434     offset = field*dmmoab->n;
435     for (i=0; i<npoints; ++i)
436       dof[i*dmmoab->numFields+field] = (dmmoab->bs > 1 ? dmmoab->gidmap[(PetscInt)points[i]-dmmoab->seqstart]*dmmoab->numFields+field :
437                                                          dmmoab->gidmap[(PetscInt)points[i]-dmmoab->seqstart]+offset);
438   }
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "DMMoabGetDofsLocal"
444 /*@C
445   DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an
446   array of MOAB EntityHandles.
447 
448   Not Collective
449 
450   Input Parameters:
451 + dm     - the discretization manager object
452 . npoints - the total number of Entities in the points array
453 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
454 
455   Output Parameter:
456 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
457 
458   Level: intermediate
459 
460 .keywords: discretization manager, get, global degrees of freedom
461 
462 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofs(), DMMoabGetDofsBlocked()
463 @*/
464 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
465 {
466   PetscInt        i,field,offset;
467   PetscErrorCode  ierr;
468   DM_Moab        *dmmoab;
469 
470   PetscFunctionBegin;
471   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
472   PetscValidPointer(points,3);
473   dmmoab = (DM_Moab*)(dm)->data;
474 
475   if (!dof) {
476     ierr = PetscMalloc1(dmmoab->numFields*npoints,&dof);CHKERRQ(ierr);
477   }
478 
479   /* compute the DOF based on local blocking in the fields */
480   /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
481   /* TODO: eliminate the limitation using PetscSection to manage DOFs */
482   for (field=0; field<dmmoab->numFields; ++field) {
483     offset = field*dmmoab->n;
484     for (i=0; i<npoints; ++i)
485       dof[i*dmmoab->numFields+field] = (dmmoab->bs > 1 ? dmmoab->lidmap[(PetscInt)points[i]-dmmoab->seqstart]*dmmoab->numFields+field :
486                                                          dmmoab->lidmap[(PetscInt)points[i]-dmmoab->seqstart]+offset);
487   }
488   PetscFunctionReturn(0);
489 }
490 
491 
492 /*@C
493   DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
494   array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation
495   of element residuals and assembly of the discrete systems when all fields are co-located.
496 
497   Not Collective
498 
499   Input Parameters:
500 + dm     - the discretization manager object
501 . npoints - the total number of Entities in the points array
502 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
503 
504   Output Parameter:
505 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat)
506 
507   Level: intermediate
508 
509 .keywords: discretization manager, get, global degrees of freedom
510 
511 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal()
512 @*/
513 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
514 {
515   PetscInt        i;
516   DM_Moab        *dmmoab;
517   PetscErrorCode  ierr;
518 
519   PetscFunctionBegin;
520   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
521   PetscValidPointer(points,3);
522   dmmoab = (DM_Moab*)(dm)->data;
523 
524   if (!dof) {
525     ierr = PetscMalloc1(npoints,&dof);CHKERRQ(ierr);
526   }
527 
528   for (i=0; i<npoints; ++i) {
529     dof[i]=dmmoab->gidmap[(PetscInt)points[i]-dmmoab->seqstart];
530   }
531   PetscFunctionReturn(0);
532 }
533 
534 
535 /*@C
536   DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
537   array of MOAB EntityHandles. It is useful when performing local 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 local 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 DMMoabGetDofsBlockedLocal(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,3);
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->lidmap[(PetscInt)points[i]-dmmoab->seqstart];
573   PetscFunctionReturn(0);
574 }
575 
576 
577 /*@C
578   DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
579   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
580   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
581 
582   Not Collective
583 
584   Input Parameters:
585 + dm     - the discretization manager object
586 
587   Output Parameter:
588 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
589 
590   Level: intermediate
591 
592 .keywords: discretization manager, get, blocked degrees of freedom
593 
594 .seealso: DMMoabGetVertexDofsBlockedLocal(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
595 @*/
596 PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof)
597 {
598   DM_Moab        *dmmoab;
599 
600   PetscFunctionBegin;
601   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
602   dmmoab = (DM_Moab*)(dm)->data;
603 
604   *dof = dmmoab->gidmap;
605   PetscFunctionReturn(0);
606 }
607 
608 
609 /*@C
610   DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
611   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
612   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
613 
614   Not Collective
615 
616   Input Parameters:
617 + dm     - the discretization manager object
618 
619   Output Parameter:
620 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
621 
622   Level: intermediate
623 
624 .keywords: discretization manager, get, blocked degrees of freedom
625 
626 .seealso: DMMoabGetVertexDofsBlocked(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
627 @*/
628 PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof)
629 {
630   DM_Moab        *dmmoab;
631 
632   PetscFunctionBegin;
633   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
634   PetscValidPointer(dof,2);
635   dmmoab = (DM_Moab*)(dm)->data;
636 
637   *dof = dmmoab->lidmap;
638   PetscFunctionReturn(0);
639 }
640 
641