xref: /petsc/src/dm/impls/moab/dmmbfield.cxx (revision 20e9df35cdcd4867fd9770ed14d3194493387312)
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 = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&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   *varray;
91   PetscScalar   *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 = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&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,&varray);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]=varray[ifield*dmmoab->nloc+i];
119         else
120           farray[i]=varray[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,&varray);CHKERRQ(ierr);
127   }
128   else {
129     ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&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         farray[i] = varray[i*dmmoab->bs+ifield];
142       }
143 
144       merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr);
145       /* make sure the parallel exchange for ghosts are done appropriately */
146       merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr);
147     }
148     ierr = PetscFree(varray);CHKERRQ(ierr);
149   }
150   ierr = PetscFree(farray);CHKERRQ(ierr);
151   PetscFunctionReturn(0);
152 }
153 
154 
155 #undef __FUNCT__
156 #define __FUNCT__ "DMMoabSetFieldNames"
157 /*@C
158   DMMoabSetFieldNames - Sets the number of fields and their names to be managed by the DM
159 
160   Not Collective
161 
162   Input Parameters:
163 + dm     - the discretization manager object
164 . numFields - the total number of fields
165 . fields - the array containing the names of each field (component); Can be NULL.
166 
167   Level: intermediate
168 
169 .keywords: discretization manager, set, component name
170 
171 .seealso: DMMoabGetFieldName(), DMMoabSetFieldName()
172 @*/
173 PetscErrorCode DMMoabSetFieldNames(DM dm,PetscInt numFields,const char** fields)
174 {
175   PetscErrorCode ierr;
176   PetscInt       i;
177   DM_Moab        *dmmoab;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
181   dmmoab = (DM_Moab*)(dm)->data;
182 
183   /* first deallocate the existing field structure */
184   if (dmmoab->fieldNames) {
185     for(i=0; i<dmmoab->numFields; i++) {
186       ierr = PetscFree(dmmoab->fieldNames[i]);CHKERRQ(ierr);
187     }
188     ierr = PetscFree(dmmoab->fieldNames);CHKERRQ(ierr);
189   }
190 
191   /* now re-allocate and assign field names  */
192   dmmoab->numFields = numFields;
193   ierr = PetscMalloc(sizeof(char*)*numFields,&dmmoab->fieldNames);CHKERRQ(ierr);
194   if (fields) {
195     for(i=0; i<dmmoab->numFields; i++) {
196       ierr = PetscStrallocpy(fields[i], (char**) &dmmoab->fieldNames[i]);CHKERRQ(ierr);
197     }
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "DMMoabGetFieldName"
205 /*@C
206   DMMoabGetFieldName - Gets the names of individual field components in multicomponent
207   vectors associated with a DMDA.
208 
209   Not Collective
210 
211   Input Parameter:
212 + dm     - the discretization manager object
213 . field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the
214         number of degrees of freedom per node within the DMMoab
215 
216   Output Parameter:
217 . fieldName - the name of the field (component)
218 
219   Level: intermediate
220 
221 .keywords: discretization manager, get, component name
222 
223 .seealso: DMMoabSetFieldName(), DMMoabSetFields()
224 @*/
225 PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName)
226 {
227   DM_Moab        *dmmoab;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
231   dmmoab = (DM_Moab*)(dm)->data;
232   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);
233 
234   *fieldName = dmmoab->fieldNames[field];
235   PetscFunctionReturn(0);
236 }
237 
238 
239 #undef __FUNCT__
240 #define __FUNCT__ "DMMoabSetFieldName"
241 /*@C
242   DMMoabSetFieldName - Sets the name of a field (component) managed by the DM
243 
244   Not Collective
245 
246   Input Parameters:
247 + dm     - the discretization manager object
248 . field - the field number
249 . fieldName - the field (component) name
250 
251   Level: intermediate
252   Notes: Can only be called after DMMoabSetFields supplied with correct numFields
253 
254 .keywords: discretization manager, set, component name
255 
256 .seealso: DMMoabGetFieldName(), DMMoabSetFields()
257 @*/
258 PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName)
259 {
260   PetscErrorCode ierr;
261   DM_Moab        *dmmoab;
262 
263   PetscFunctionBegin;
264   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
265   PetscValidCharPointer(fieldName,3);
266 
267   dmmoab = (DM_Moab*)(dm)->data;
268   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);
269 
270   if (dmmoab->fieldNames[field]) {
271     ierr = PetscFree(dmmoab->fieldNames[field]);CHKERRQ(ierr);
272   }
273   ierr = PetscStrallocpy(fieldName, (char**) &dmmoab->fieldNames[field]);CHKERRQ(ierr);
274   PetscFunctionReturn(0);
275 }
276 
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "DMMoabGetFieldDof"
280 /*@C
281   DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a
282   particular MOAB EntityHandle.
283 
284   Not Collective
285 
286   Input Parameters:
287 + dm     - the discretization manager object
288 . point - the MOAB EntityHandle container which holds the field degree-of-freedom values
289 . field - the field (component) index
290 
291   Output Parameter:
292 + dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat)
293 
294   Level: beginner
295 
296 .keywords: discretization manager, get, global degree of freedom
297 
298 .seealso: DMMoabGetFieldDofs(), DMMoabGetFieldDofsLocal()
299 @*/
300 PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof)
301 {
302   DM_Moab        *dmmoab;
303 
304   PetscFunctionBegin;
305   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
306   dmmoab = (DM_Moab*)(dm)->data;
307 
308   *dof=dmmoab->gidmap[(PetscInt)point];
309   PetscFunctionReturn(0);
310 }
311 
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "DMMoabGetFieldDofs"
315 /*@C
316   DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an
317   array of MOAB EntityHandles.
318 
319   Not Collective
320 
321   Input Parameters:
322 + dm     - the discretization manager object
323 . npoints - the total number of Entities in the points array
324 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
325 . field - the field (component) index
326 
327   Output Parameter:
328 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
329 
330   Level: intermediate
331 
332 .keywords: discretization manager, get, global degrees of freedom
333 
334 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofsLocal()
335 @*/
336 PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
337 {
338   PetscInt        i;
339   PetscErrorCode  ierr;
340   DM_Moab        *dmmoab;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
344   PetscValidPointer(points,3);
345   dmmoab = (DM_Moab*)(dm)->data;
346 
347   if (!dof) {
348     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
349   }
350 
351   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
352     for (i=0; i<npoints; ++i)
353       dof[i] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
354   }
355   else {
356     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
357     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
358     for (i=0; i<npoints; ++i)
359       dof[i] = dmmoab->gidmap[(PetscInt)points[i]]+field*dmmoab->n;
360   }
361   PetscFunctionReturn(0);
362 }
363 
364 
365 #undef __FUNCT__
366 #define __FUNCT__ "DMMoabGetFieldDofsLocal"
367 /*@C
368   DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an
369   array of MOAB EntityHandles.
370 
371   Not Collective
372 
373   Input Parameters:
374 + dm     - the discretization manager object
375 . npoints - the total number of Entities in the points array
376 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
377 . field - the field (component) index
378 
379   Output Parameter:
380 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
381 
382   Level: intermediate
383 
384 .keywords: discretization manager, get, local degrees of freedom
385 
386 .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofs()
387 @*/
388 PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof)
389 {
390   PetscInt i,offset;
391   PetscErrorCode  ierr;
392   DM_Moab        *dmmoab;
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
396   PetscValidPointer(points,3);
397   dmmoab = (DM_Moab*)(dm)->data;
398 
399   if (!dof) {
400     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
401   }
402 
403   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
404     for (i=0; i<npoints; ++i)
405       dof[i] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
406   }
407   else {
408     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
409     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
410     offset = field*dmmoab->n;
411     for (i=0; i<npoints; ++i)
412       dof[i] = dmmoab->lidmap[(PetscInt)points[i]]+offset;
413   }
414   PetscFunctionReturn(0);
415 }
416 
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "DMMoabGetDofs"
420 /*@C
421   DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an
422   array of MOAB EntityHandles.
423 
424   Not Collective
425 
426   Input Parameters:
427 + dm     - the discretization manager object
428 . npoints - the total number of Entities in the points array
429 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
430 
431   Output Parameter:
432 + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
433 
434   Level: intermediate
435 
436 .keywords: discretization manager, get, global degrees of freedom
437 
438 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofsLocal(), DMMoabGetDofsBlocked()
439 @*/
440 PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
441 {
442   PetscInt        i,field,offset;
443   PetscErrorCode  ierr;
444   DM_Moab        *dmmoab;
445 
446   PetscFunctionBegin;
447   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
448   PetscValidPointer(points,3);
449   dmmoab = (DM_Moab*)(dm)->data;
450 
451   if (!dof) {
452     ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->numFields*npoints, &dof);CHKERRQ(ierr);
453   }
454 
455   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
456     for (field=0; field<dmmoab->numFields; ++field) {
457       for (i=0; i<npoints; ++i)
458         dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
459     }
460   }
461   else {
462     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
463     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
464     for (field=0; field<dmmoab->numFields; ++field) {
465       offset = field*dmmoab->n;
466       for (i=0; i<npoints; ++i)
467         dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]+offset;
468     }
469   }
470   PetscFunctionReturn(0);
471 }
472 
473 
474 #undef __FUNCT__
475 #define __FUNCT__ "DMMoabGetDofsLocal"
476 /*@C
477   DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an
478   array of MOAB EntityHandles.
479 
480   Not Collective
481 
482   Input Parameters:
483 + dm     - the discretization manager object
484 . npoints - the total number of Entities in the points array
485 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
486 
487   Output Parameter:
488 + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat)
489 
490   Level: intermediate
491 
492 .keywords: discretization manager, get, global degrees of freedom
493 
494 .seealso: DMMoabGetFieldDofs(), DMMoabGetDofs(), DMMoabGetDofsBlocked()
495 @*/
496 PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
497 {
498   PetscInt        i,field,offset;
499   PetscErrorCode  ierr;
500   DM_Moab        *dmmoab;
501 
502   PetscFunctionBegin;
503   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
504   PetscValidPointer(points,3);
505   dmmoab = (DM_Moab*)(dm)->data;
506 
507   if (!dof) {
508     ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->numFields*npoints, &dof);CHKERRQ(ierr);
509   }
510 
511   if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */
512     for (field=0; field<dmmoab->numFields; ++field) {
513       for (i=0; i<npoints; ++i)
514         dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field;
515     }
516   }
517   else {
518     /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */
519     /* TODO: eliminate the limitation using PetscSection to manage DOFs */
520     for (field=0; field<dmmoab->numFields; ++field) {
521       offset = field*dmmoab->n;
522       for (i=0; i<npoints; ++i)
523         dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]+offset;
524     }
525   }
526   PetscFunctionReturn(0);
527 }
528 
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "DMMoabGetDofsBlocked"
532 /*@C
533   DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
534   array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation
535   of element residuals and assembly of the discrete systems when all fields are co-located.
536 
537   Not Collective
538 
539   Input Parameters:
540 + dm     - the discretization manager object
541 . npoints - the total number of Entities in the points array
542 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
543 
544   Output Parameter:
545 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat)
546 
547   Level: intermediate
548 
549 .keywords: discretization manager, get, global degrees of freedom
550 
551 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal()
552 @*/
553 PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
554 {
555   PetscInt        i;
556   DM_Moab        *dmmoab;
557   PetscErrorCode  ierr;
558 
559   PetscFunctionBegin;
560   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
561   PetscValidPointer(points,2);
562   dmmoab = (DM_Moab*)(dm)->data;
563 
564   if (!dof) {
565     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
566   }
567 
568   for (i=0; i<npoints; ++i) {
569     dof[i]=dmmoab->gidmap[(PetscInt)points[i]];
570   }
571   PetscFunctionReturn(0);
572 }
573 
574 
575 #undef __FUNCT__
576 #define __FUNCT__ "DMMoabGetDofsBlockedLocal"
577 /*@C
578   DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
579   array of MOAB EntityHandles. It is useful when performing local Blocked(Get/Set) methods in computation
580   of element residuals and assembly of the discrete systems when all fields are co-located.
581 
582   Not Collective
583 
584   Input Parameters:
585 + dm     - the discretization manager object
586 . npoints - the total number of Entities in the points array
587 . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values
588 
589   Output Parameter:
590 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat)
591 
592   Level: intermediate
593 
594 .keywords: discretization manager, get, global degrees of freedom
595 
596 .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal()
597 @*/
598 PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof)
599 {
600   PetscInt        i;
601   DM_Moab        *dmmoab;
602   PetscErrorCode  ierr;
603 
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
606   PetscValidPointer(points,2);
607   dmmoab = (DM_Moab*)(dm)->data;
608 
609   if (!dof) {
610     ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr);
611   }
612 
613   for (i=0; i<npoints; ++i)
614     dof[i] = dmmoab->lidmap[(PetscInt)points[i]];
615   PetscFunctionReturn(0);
616 }
617 
618 
619 #undef __FUNCT__
620 #define __FUNCT__ "DMMoabGetVertexDofsBlocked"
621 /*@C
622   DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an
623   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
624   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
625 
626   Not Collective
627 
628   Input Parameters:
629 + dm     - the discretization manager object
630 
631   Output Parameter:
632 + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
633 
634   Level: intermediate
635 
636 .keywords: discretization manager, get, blocked degrees of freedom
637 
638 .seealso: DMMoabGetVertexDofsBlockedLocal(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
639 @*/
640 PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof)
641 {
642   DM_Moab        *dmmoab;
643 
644   PetscFunctionBegin;
645   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
646   dmmoab = (DM_Moab*)(dm)->data;
647 
648   *dof = dmmoab->gidmap;
649   PetscFunctionReturn(0);
650 }
651 
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "DMMoabGetVertexDofsBlockedLocal"
655 /*@C
656   DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an
657   array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations
658   where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations.
659 
660   Not Collective
661 
662   Input Parameters:
663 + dm     - the discretization manager object
664 
665   Output Parameter:
666 + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering
667 
668   Level: intermediate
669 
670 .keywords: discretization manager, get, blocked degrees of freedom
671 
672 .seealso: DMMoabGetVertexDofsBlocked(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal()
673 @*/
674 PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof)
675 {
676   DM_Moab        *dmmoab;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
680   PetscValidPointer(dof,2);
681   dmmoab = (DM_Moab*)(dm)->data;
682 
683   *dof = dmmoab->lidmap;
684   PetscFunctionReturn(0);
685 }
686 
687