1*b8ecf6d3SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 28cbae1a6SVijay Mahadevan 38cbae1a6SVijay Mahadevan #include <petscdmmoab.h> 48cbae1a6SVijay Mahadevan 58cbae1a6SVijay Mahadevan #undef __FUNCT__ 68cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabSetFieldVector" 799fa7e03SVijay Mahadevan /*@C 899fa7e03SVijay Mahadevan DMMoabSetFieldVector - Sets the vector reference that represents the solution associated 999fa7e03SVijay Mahadevan with a particular field component. 1099fa7e03SVijay Mahadevan 1199fa7e03SVijay Mahadevan Not Collective 1299fa7e03SVijay Mahadevan 1399fa7e03SVijay Mahadevan Input Parameters: 1499fa7e03SVijay Mahadevan + dm - the discretization manager object 1599fa7e03SVijay Mahadevan . ifield - the index of the field as set before via DMMoabSetFieldName. 16*b8ecf6d3SVijay Mahadevan . fvec - the Vector solution corresponding to the field (component) 1799fa7e03SVijay Mahadevan 1899fa7e03SVijay Mahadevan Level: intermediate 1999fa7e03SVijay Mahadevan 2099fa7e03SVijay Mahadevan .keywords: discretization manager, set, component solution 2199fa7e03SVijay Mahadevan 2299fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldName(), DMMoabSetGlobalFieldVector() 2399fa7e03SVijay Mahadevan @*/ 248cbae1a6SVijay Mahadevan PetscErrorCode DMMoabSetFieldVector(DM dm, PetscInt ifield, Vec fvec) 258cbae1a6SVijay Mahadevan { 268cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 278cbae1a6SVijay Mahadevan moab::Tag vtag,ntag; 288cbae1a6SVijay Mahadevan const PetscScalar *varray; 298cbae1a6SVijay Mahadevan PetscScalar *farray; 308cbae1a6SVijay Mahadevan moab::ErrorCode merr; 318cbae1a6SVijay Mahadevan PetscErrorCode ierr; 328cbae1a6SVijay Mahadevan std::string tag_name; 338cbae1a6SVijay Mahadevan 348cbae1a6SVijay Mahadevan PetscFunctionBegin; 358cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 368cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 378cbae1a6SVijay Mahadevan 38addae81cSVijay Mahadevan 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); 39addae81cSVijay Mahadevan 408cbae1a6SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 41addae81cSVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 428cbae1a6SVijay Mahadevan moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 438cbae1a6SVijay Mahadevan 448cbae1a6SVijay Mahadevan ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 458cbae1a6SVijay Mahadevan 468cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 478cbae1a6SVijay Mahadevan if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 488cbae1a6SVijay Mahadevan ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 498cbae1a6SVijay Mahadevan /* use the entity handle and the Dof index to set the right value */ 508cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)varray);MBERRNM(merr); 518cbae1a6SVijay Mahadevan ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 528cbae1a6SVijay Mahadevan } 538cbae1a6SVijay Mahadevan else { 548cbae1a6SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 558cbae1a6SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 568cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)farray);MBERRNM(merr); 578cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 588cbae1a6SVijay Mahadevan /* make sure the parallel exchange for ghosts are done appropriately */ 598cbae1a6SVijay Mahadevan ierr = PetscFree(farray);CHKERRQ(ierr); 608cbae1a6SVijay Mahadevan } 618cbae1a6SVijay Mahadevan merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vowned);MBERRNM(merr); 628cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 638cbae1a6SVijay Mahadevan } 648cbae1a6SVijay Mahadevan 658cbae1a6SVijay Mahadevan 668cbae1a6SVijay Mahadevan #undef __FUNCT__ 678cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabSetGlobalFieldVector" 6899fa7e03SVijay Mahadevan /*@C 6999fa7e03SVijay Mahadevan DMMoabSetGlobalFieldVector - Sets the vector reference that represents the global solution associated 7099fa7e03SVijay Mahadevan with all fields (components) managed by DM. 7199fa7e03SVijay Mahadevan A useful utility when updating the DM solution after a solve, to be serialized with the mesh for 7299fa7e03SVijay Mahadevan checkpointing purposes. 7399fa7e03SVijay Mahadevan 7499fa7e03SVijay Mahadevan Not Collective 7599fa7e03SVijay Mahadevan 7699fa7e03SVijay Mahadevan Input Parameters: 7799fa7e03SVijay Mahadevan + dm - the discretization manager object 78*b8ecf6d3SVijay Mahadevan . fvec - the global Vector solution corresponding to all the fields managed by DM 7999fa7e03SVijay Mahadevan 8099fa7e03SVijay Mahadevan Level: intermediate 8199fa7e03SVijay Mahadevan 8299fa7e03SVijay Mahadevan .keywords: discretization manager, set, component solution 8399fa7e03SVijay Mahadevan 8499fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldName(), DMMoabSetFieldVector() 8599fa7e03SVijay Mahadevan @*/ 868cbae1a6SVijay Mahadevan PetscErrorCode DMMoabSetGlobalFieldVector(DM dm, Vec fvec) 878cbae1a6SVijay Mahadevan { 888cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 898cbae1a6SVijay Mahadevan moab::Tag vtag,ntag; 908cbae1a6SVijay Mahadevan const PetscScalar *varray; 918cbae1a6SVijay Mahadevan PetscScalar *farray; 928cbae1a6SVijay Mahadevan moab::ErrorCode merr; 938cbae1a6SVijay Mahadevan PetscErrorCode ierr; 948cbae1a6SVijay Mahadevan PetscInt i,ifield; 958cbae1a6SVijay Mahadevan std::string tag_name; 968cbae1a6SVijay Mahadevan moab::Range::iterator iter; 978cbae1a6SVijay Mahadevan 988cbae1a6SVijay Mahadevan PetscFunctionBegin; 998cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1008cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1018cbae1a6SVijay Mahadevan 1028cbae1a6SVijay Mahadevan /* get the Tag corresponding to the global vector - possible that there is no tag associated.. */ 1038cbae1a6SVijay Mahadevan ierr = DMMoabGetVecTag(fvec,&vtag);CHKERRQ(ierr); 1048cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_name(vtag, tag_name); 1058cbae1a6SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*sizeof(PetscScalar),&farray);CHKERRQ(ierr); 1068cbae1a6SVijay Mahadevan if (!tag_name.length() && merr !=moab::MB_SUCCESS) { 1078cbae1a6SVijay Mahadevan /* not a MOAB vector - use VecGetSubVector to get the parts as needed */ 1088cbae1a6SVijay Mahadevan 1098cbae1a6SVijay Mahadevan ierr = VecGetArrayRead(fvec,&varray);CHKERRQ(ierr); 110addae81cSVijay Mahadevan for (ifield=0; ifield<dmmoab->numFields; ++ifield) { 1118cbae1a6SVijay Mahadevan 1128cbae1a6SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 113addae81cSVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 1148cbae1a6SVijay Mahadevan moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 1158cbae1a6SVijay Mahadevan 1168cbae1a6SVijay Mahadevan for(i=0;i<dmmoab->nloc;i++) { 1178cbae1a6SVijay Mahadevan if (dmmoab->bs == 1) 1188cbae1a6SVijay Mahadevan farray[i]=varray[ifield*dmmoab->nloc+i]; 1198cbae1a6SVijay Mahadevan else 120addae81cSVijay Mahadevan farray[i]=varray[i*dmmoab->numFields+ifield]; 1218cbae1a6SVijay Mahadevan } 1228cbae1a6SVijay Mahadevan 1238cbae1a6SVijay Mahadevan /* use the entity handle and the Dof index to set the right value */ 1248cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 1258cbae1a6SVijay Mahadevan } 1268cbae1a6SVijay Mahadevan ierr = VecRestoreArrayRead(fvec,&varray);CHKERRQ(ierr); 1278cbae1a6SVijay Mahadevan } 1288cbae1a6SVijay Mahadevan else { 1298cbae1a6SVijay Mahadevan ierr = PetscMalloc(dmmoab->nloc*dmmoab->bs*sizeof(PetscScalar),&varray);CHKERRQ(ierr); 1308cbae1a6SVijay Mahadevan 1318cbae1a6SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 1328cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(vtag, *dmmoab->vowned, (void*)varray);MBERRNM(merr); 133addae81cSVijay Mahadevan for (ifield=0; ifield<dmmoab->numFields; ++ifield) { 1348cbae1a6SVijay Mahadevan 1358cbae1a6SVijay Mahadevan /* Create a tag in MOAB mesh to index and keep track of number of Petsc vec tags */ 136addae81cSVijay Mahadevan merr = dmmoab->mbiface->tag_get_handle(dmmoab->fieldNames[ifield],1,moab::MB_TYPE_DOUBLE,ntag, 1378cbae1a6SVijay Mahadevan moab::MB_TAG_DENSE|moab::MB_TAG_CREAT);MBERRNM(merr); 1388cbae1a6SVijay Mahadevan 1398cbae1a6SVijay Mahadevan /* we are using a MOAB Vec - directly copy the tag data to new one */ 1408cbae1a6SVijay Mahadevan for(i=0; i < dmmoab->nloc; i++) { 1418cbae1a6SVijay Mahadevan farray[i] = varray[i*dmmoab->bs+ifield]; 1428cbae1a6SVijay Mahadevan } 1438cbae1a6SVijay Mahadevan 1448cbae1a6SVijay Mahadevan merr = dmmoab->mbiface->tag_set_data(ntag, *dmmoab->vowned, (const void*)farray);MBERRNM(merr); 1458cbae1a6SVijay Mahadevan /* make sure the parallel exchange for ghosts are done appropriately */ 1468cbae1a6SVijay Mahadevan merr = dmmoab->pcomm->exchange_tags(ntag, *dmmoab->vlocal);MBERRNM(merr); 1478cbae1a6SVijay Mahadevan } 1488cbae1a6SVijay Mahadevan ierr = PetscFree(varray);CHKERRQ(ierr); 1498cbae1a6SVijay Mahadevan } 1508cbae1a6SVijay Mahadevan ierr = PetscFree(farray);CHKERRQ(ierr); 1518cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 1528cbae1a6SVijay Mahadevan } 1538cbae1a6SVijay Mahadevan 1548cbae1a6SVijay Mahadevan 1558cbae1a6SVijay Mahadevan #undef __FUNCT__ 15699fa7e03SVijay Mahadevan #define __FUNCT__ "DMMoabSetFieldNames" 15799fa7e03SVijay Mahadevan /*@C 15899fa7e03SVijay Mahadevan DMMoabSetFieldNames - Sets the number of fields and their names to be managed by the DM 15999fa7e03SVijay Mahadevan 16099fa7e03SVijay Mahadevan Not Collective 16199fa7e03SVijay Mahadevan 16299fa7e03SVijay Mahadevan Input Parameters: 16399fa7e03SVijay Mahadevan + dm - the discretization manager object 16499fa7e03SVijay Mahadevan . numFields - the total number of fields 165*b8ecf6d3SVijay Mahadevan . fields - the array containing the names of each field (component); Can be NULL. 16699fa7e03SVijay Mahadevan 16799fa7e03SVijay Mahadevan Level: intermediate 16899fa7e03SVijay Mahadevan 16999fa7e03SVijay Mahadevan .keywords: discretization manager, set, component name 17099fa7e03SVijay Mahadevan 17199fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldName(), DMMoabSetFieldName() 17299fa7e03SVijay Mahadevan @*/ 17399fa7e03SVijay Mahadevan PetscErrorCode DMMoabSetFieldNames(DM dm,PetscInt numFields,const char** fields) 1748cbae1a6SVijay Mahadevan { 175addae81cSVijay Mahadevan PetscErrorCode ierr; 176addae81cSVijay Mahadevan PetscInt i; 1778cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 1788cbae1a6SVijay Mahadevan 1798cbae1a6SVijay Mahadevan PetscFunctionBegin; 1808cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1818cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 1828cbae1a6SVijay Mahadevan 183addae81cSVijay Mahadevan /* first deallocate the existing field structure */ 184addae81cSVijay Mahadevan if (dmmoab->fieldNames) { 185addae81cSVijay Mahadevan for(i=0; i<dmmoab->numFields; i++) { 186addae81cSVijay Mahadevan ierr = PetscFree(dmmoab->fieldNames[i]);CHKERRQ(ierr); 187addae81cSVijay Mahadevan } 188addae81cSVijay Mahadevan ierr = PetscFree(dmmoab->fieldNames);CHKERRQ(ierr); 189addae81cSVijay Mahadevan } 190addae81cSVijay Mahadevan 191addae81cSVijay Mahadevan /* now re-allocate and assign field names */ 192addae81cSVijay Mahadevan dmmoab->numFields = numFields; 193addae81cSVijay Mahadevan ierr = PetscMalloc(sizeof(char*)*numFields,&dmmoab->fieldNames);CHKERRQ(ierr); 194addae81cSVijay Mahadevan if (fields) { 195addae81cSVijay Mahadevan for(i=0; i<dmmoab->numFields; i++) { 196addae81cSVijay Mahadevan ierr = PetscStrallocpy(fields[i], (char**) &dmmoab->fieldNames[i]);CHKERRQ(ierr); 197addae81cSVijay Mahadevan } 198addae81cSVijay Mahadevan } 1998cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 2008cbae1a6SVijay Mahadevan } 2018cbae1a6SVijay Mahadevan 20299fa7e03SVijay Mahadevan 20399fa7e03SVijay Mahadevan #undef __FUNCT__ 20499fa7e03SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldName" 20599fa7e03SVijay Mahadevan /*@C 20699fa7e03SVijay Mahadevan DMMoabGetFieldName - Gets the names of individual field components in multicomponent 20799fa7e03SVijay Mahadevan vectors associated with a DMDA. 20899fa7e03SVijay Mahadevan 20999fa7e03SVijay Mahadevan Not Collective 21099fa7e03SVijay Mahadevan 21199fa7e03SVijay Mahadevan Input Parameter: 21299fa7e03SVijay Mahadevan + dm - the discretization manager object 213*b8ecf6d3SVijay Mahadevan . field - field number for the DMMoab (0, 1, ... dof-1), where dof indicates the 21499fa7e03SVijay Mahadevan number of degrees of freedom per node within the DMMoab 21599fa7e03SVijay Mahadevan 21699fa7e03SVijay Mahadevan Output Parameter: 21799fa7e03SVijay Mahadevan . fieldName - the name of the field (component) 21899fa7e03SVijay Mahadevan 21999fa7e03SVijay Mahadevan Level: intermediate 22099fa7e03SVijay Mahadevan 22199fa7e03SVijay Mahadevan .keywords: discretization manager, get, component name 22299fa7e03SVijay Mahadevan 22399fa7e03SVijay Mahadevan .seealso: DMMoabSetFieldName(), DMMoabSetFields() 22499fa7e03SVijay Mahadevan @*/ 22599fa7e03SVijay Mahadevan PetscErrorCode DMMoabGetFieldName(DM dm, PetscInt field, const char **fieldName) 22699fa7e03SVijay Mahadevan { 22799fa7e03SVijay Mahadevan DM_Moab *dmmoab; 22899fa7e03SVijay Mahadevan 22999fa7e03SVijay Mahadevan PetscFunctionBegin; 23099fa7e03SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 23199fa7e03SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 23299fa7e03SVijay Mahadevan 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); 23399fa7e03SVijay Mahadevan 23499fa7e03SVijay Mahadevan *fieldName = dmmoab->fieldNames[field]; 23599fa7e03SVijay Mahadevan PetscFunctionReturn(0); 23699fa7e03SVijay Mahadevan } 23799fa7e03SVijay Mahadevan 23899fa7e03SVijay Mahadevan 239c68c6878SVijay Mahadevan #undef __FUNCT__ 240c68c6878SVijay Mahadevan #define __FUNCT__ "DMMoabSetFieldName" 241c68c6878SVijay Mahadevan /*@C 24299fa7e03SVijay Mahadevan DMMoabSetFieldName - Sets the name of a field (component) managed by the DM 243c68c6878SVijay Mahadevan 244c68c6878SVijay Mahadevan Not Collective 245c68c6878SVijay Mahadevan 246c68c6878SVijay Mahadevan Input Parameters: 24799fa7e03SVijay Mahadevan + dm - the discretization manager object 248c68c6878SVijay Mahadevan . field - the field number 249*b8ecf6d3SVijay Mahadevan . fieldName - the field (component) name 250c68c6878SVijay Mahadevan 25199fa7e03SVijay Mahadevan Level: intermediate 25299fa7e03SVijay Mahadevan Notes: Can only be called after DMMoabSetFields supplied with correct numFields 253c68c6878SVijay Mahadevan 25499fa7e03SVijay Mahadevan .keywords: discretization manager, set, component name 25599fa7e03SVijay Mahadevan 25699fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldName(), DMMoabSetFields() 257c68c6878SVijay Mahadevan @*/ 2585f7ae369SVijay Mahadevan PetscErrorCode DMMoabSetFieldName(DM dm, PetscInt field, const char *fieldName) 259c68c6878SVijay Mahadevan { 260c68c6878SVijay Mahadevan PetscErrorCode ierr; 261c68c6878SVijay Mahadevan DM_Moab *dmmoab; 262c68c6878SVijay Mahadevan 263c68c6878SVijay Mahadevan PetscFunctionBegin; 264c68c6878SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2655f7ae369SVijay Mahadevan PetscValidCharPointer(fieldName,3); 266c68c6878SVijay Mahadevan 26799fa7e03SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 268c68c6878SVijay Mahadevan 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); 26999fa7e03SVijay Mahadevan 27099fa7e03SVijay Mahadevan if (dmmoab->fieldNames[field]) { 271c68c6878SVijay Mahadevan ierr = PetscFree(dmmoab->fieldNames[field]);CHKERRQ(ierr); 27299fa7e03SVijay Mahadevan } 2735f7ae369SVijay Mahadevan ierr = PetscStrallocpy(fieldName, (char**) &dmmoab->fieldNames[field]);CHKERRQ(ierr); 274c68c6878SVijay Mahadevan PetscFunctionReturn(0); 275c68c6878SVijay Mahadevan } 276c68c6878SVijay Mahadevan 2778cbae1a6SVijay Mahadevan 2788cbae1a6SVijay Mahadevan #undef __FUNCT__ 2798cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDof" 28099fa7e03SVijay Mahadevan /*@C 28199fa7e03SVijay Mahadevan DMMoabGetFieldDof - Gets the global degree-of-freedom of a field (component) defined on a 28299fa7e03SVijay Mahadevan particular MOAB EntityHandle. 28399fa7e03SVijay Mahadevan 28499fa7e03SVijay Mahadevan Not Collective 28599fa7e03SVijay Mahadevan 28699fa7e03SVijay Mahadevan Input Parameters: 28799fa7e03SVijay Mahadevan + dm - the discretization manager object 28899fa7e03SVijay Mahadevan . point - the MOAB EntityHandle container which holds the field degree-of-freedom values 289*b8ecf6d3SVijay Mahadevan . field - the field (component) index 29099fa7e03SVijay Mahadevan 29199fa7e03SVijay Mahadevan Output Parameter: 29299fa7e03SVijay Mahadevan + dof - the global degree-of-freedom index corresponding to the field in the discrete representation (Vec, Mat) 29399fa7e03SVijay Mahadevan 29499fa7e03SVijay Mahadevan Level: beginner 29599fa7e03SVijay Mahadevan 29699fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degree of freedom 29799fa7e03SVijay Mahadevan 29899fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldDofs(), DMMoabGetFieldDofsLocal() 29999fa7e03SVijay Mahadevan @*/ 3008cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetFieldDof(DM dm,moab::EntityHandle point,PetscInt field,PetscInt* dof) 3018cbae1a6SVijay Mahadevan { 3028cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 3038cbae1a6SVijay Mahadevan 3048cbae1a6SVijay Mahadevan PetscFunctionBegin; 3058cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3068cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 3078cbae1a6SVijay Mahadevan 3088cbae1a6SVijay Mahadevan *dof=dmmoab->gidmap[(PetscInt)point]; 3098cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 3108cbae1a6SVijay Mahadevan } 3118cbae1a6SVijay Mahadevan 3128cbae1a6SVijay Mahadevan 3138cbae1a6SVijay Mahadevan #undef __FUNCT__ 3148cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDofs" 31599fa7e03SVijay Mahadevan /*@C 31699fa7e03SVijay Mahadevan DMMoabGetFieldDofs - Gets the global degree-of-freedom of a field (component) defined on an 31799fa7e03SVijay Mahadevan array of MOAB EntityHandles. 31899fa7e03SVijay Mahadevan 31999fa7e03SVijay Mahadevan Not Collective 32099fa7e03SVijay Mahadevan 32199fa7e03SVijay Mahadevan Input Parameters: 32299fa7e03SVijay Mahadevan + dm - the discretization manager object 32399fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 32499fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 325*b8ecf6d3SVijay Mahadevan . field - the field (component) index 32699fa7e03SVijay Mahadevan 32799fa7e03SVijay Mahadevan Output Parameter: 32899fa7e03SVijay Mahadevan + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 32999fa7e03SVijay Mahadevan 33099fa7e03SVijay Mahadevan Level: intermediate 33199fa7e03SVijay Mahadevan 33299fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degrees of freedom 33399fa7e03SVijay Mahadevan 33499fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofsLocal() 33599fa7e03SVijay Mahadevan @*/ 3368cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetFieldDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 3378cbae1a6SVijay Mahadevan { 3388cbae1a6SVijay Mahadevan PetscInt i; 3398cbae1a6SVijay Mahadevan PetscErrorCode ierr; 3408cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 3418cbae1a6SVijay Mahadevan 3428cbae1a6SVijay Mahadevan PetscFunctionBegin; 3438cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 34499fa7e03SVijay Mahadevan PetscValidPointer(points,3); 3458cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 3468cbae1a6SVijay Mahadevan 3478cbae1a6SVijay Mahadevan if (!dof) { 3488cbae1a6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 3498cbae1a6SVijay Mahadevan } 3508cbae1a6SVijay Mahadevan 35199fa7e03SVijay Mahadevan if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 3528cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 353addae81cSVijay Mahadevan dof[i] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 3548cbae1a6SVijay Mahadevan } 3558cbae1a6SVijay Mahadevan else { 35699fa7e03SVijay Mahadevan /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 35799fa7e03SVijay Mahadevan /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 3588cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 3598cbae1a6SVijay Mahadevan dof[i] = dmmoab->gidmap[(PetscInt)points[i]]+field*dmmoab->n; 3608cbae1a6SVijay Mahadevan } 3618cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 3628cbae1a6SVijay Mahadevan } 3638cbae1a6SVijay Mahadevan 3648cbae1a6SVijay Mahadevan 3658cbae1a6SVijay Mahadevan #undef __FUNCT__ 3668cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetFieldDofsLocal" 36799fa7e03SVijay Mahadevan /*@C 36899fa7e03SVijay Mahadevan DMMoabGetFieldDofsLocal - Gets the local degrees-of-freedom of a field (component) defined on an 36999fa7e03SVijay Mahadevan array of MOAB EntityHandles. 37099fa7e03SVijay Mahadevan 37199fa7e03SVijay Mahadevan Not Collective 37299fa7e03SVijay Mahadevan 37399fa7e03SVijay Mahadevan Input Parameters: 37499fa7e03SVijay Mahadevan + dm - the discretization manager object 37599fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 37699fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 377*b8ecf6d3SVijay Mahadevan . field - the field (component) index 37899fa7e03SVijay Mahadevan 37999fa7e03SVijay Mahadevan Output Parameter: 38099fa7e03SVijay Mahadevan + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 38199fa7e03SVijay Mahadevan 38299fa7e03SVijay Mahadevan Level: intermediate 38399fa7e03SVijay Mahadevan 38499fa7e03SVijay Mahadevan .keywords: discretization manager, get, local degrees of freedom 38599fa7e03SVijay Mahadevan 38699fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldDof(), DMMoabGetFieldDofs() 38799fa7e03SVijay Mahadevan @*/ 3888cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetFieldDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt field,PetscInt* dof) 3898cbae1a6SVijay Mahadevan { 3908cbae1a6SVijay Mahadevan PetscInt i,offset; 3918cbae1a6SVijay Mahadevan PetscErrorCode ierr; 3928cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 3938cbae1a6SVijay Mahadevan 3948cbae1a6SVijay Mahadevan PetscFunctionBegin; 3958cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 39699fa7e03SVijay Mahadevan PetscValidPointer(points,3); 3978cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 3988cbae1a6SVijay Mahadevan 3998cbae1a6SVijay Mahadevan if (!dof) { 4008cbae1a6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 4018cbae1a6SVijay Mahadevan } 4028cbae1a6SVijay Mahadevan 40399fa7e03SVijay Mahadevan if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 4048cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 405addae81cSVijay Mahadevan dof[i] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 4068cbae1a6SVijay Mahadevan } 4078cbae1a6SVijay Mahadevan else { 40899fa7e03SVijay Mahadevan /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 40999fa7e03SVijay Mahadevan /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 41099fa7e03SVijay Mahadevan offset = field*dmmoab->n; 4118cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 4128cbae1a6SVijay Mahadevan dof[i] = dmmoab->lidmap[(PetscInt)points[i]]+offset; 4138cbae1a6SVijay Mahadevan } 4148cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 4158cbae1a6SVijay Mahadevan } 4168cbae1a6SVijay Mahadevan 4178cbae1a6SVijay Mahadevan 4188cbae1a6SVijay Mahadevan #undef __FUNCT__ 4198cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofs" 42099fa7e03SVijay Mahadevan /*@C 42199fa7e03SVijay Mahadevan DMMoabGetDofs - Gets the global degree-of-freedom for all fields (components) defined on an 42299fa7e03SVijay Mahadevan array of MOAB EntityHandles. 42399fa7e03SVijay Mahadevan 42499fa7e03SVijay Mahadevan Not Collective 42599fa7e03SVijay Mahadevan 42699fa7e03SVijay Mahadevan Input Parameters: 42799fa7e03SVijay Mahadevan + dm - the discretization manager object 42899fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 42999fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 43099fa7e03SVijay Mahadevan 43199fa7e03SVijay Mahadevan Output Parameter: 43299fa7e03SVijay Mahadevan + dof - the global degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 43399fa7e03SVijay Mahadevan 43499fa7e03SVijay Mahadevan Level: intermediate 43599fa7e03SVijay Mahadevan 43699fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degrees of freedom 43799fa7e03SVijay Mahadevan 43899fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldDofs(), DMMoabGetDofsLocal(), DMMoabGetDofsBlocked() 43999fa7e03SVijay Mahadevan @*/ 4408cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetDofs(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 4418cbae1a6SVijay Mahadevan { 4428cbae1a6SVijay Mahadevan PetscInt i,field,offset; 4438cbae1a6SVijay Mahadevan PetscErrorCode ierr; 4448cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 4458cbae1a6SVijay Mahadevan 4468cbae1a6SVijay Mahadevan PetscFunctionBegin; 4478cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 44899fa7e03SVijay Mahadevan PetscValidPointer(points,3); 4498cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 4508cbae1a6SVijay Mahadevan 4518cbae1a6SVijay Mahadevan if (!dof) { 452addae81cSVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->numFields*npoints, &dof);CHKERRQ(ierr); 4538cbae1a6SVijay Mahadevan } 4548cbae1a6SVijay Mahadevan 45599fa7e03SVijay Mahadevan if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 456addae81cSVijay Mahadevan for (field=0; field<dmmoab->numFields; ++field) { 4578cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 458addae81cSVijay Mahadevan dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 4598cbae1a6SVijay Mahadevan } 4608cbae1a6SVijay Mahadevan } 4618cbae1a6SVijay Mahadevan else { 46299fa7e03SVijay Mahadevan /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 46399fa7e03SVijay Mahadevan /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 464addae81cSVijay Mahadevan for (field=0; field<dmmoab->numFields; ++field) { 46599fa7e03SVijay Mahadevan offset = field*dmmoab->n; 4668cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 467addae81cSVijay Mahadevan dof[i*dmmoab->numFields+field] = dmmoab->gidmap[(PetscInt)points[i]]+offset; 4688cbae1a6SVijay Mahadevan } 4698cbae1a6SVijay Mahadevan } 4708cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 4718cbae1a6SVijay Mahadevan } 4728cbae1a6SVijay Mahadevan 4738cbae1a6SVijay Mahadevan 4748cbae1a6SVijay Mahadevan #undef __FUNCT__ 4758cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsLocal" 47699fa7e03SVijay Mahadevan /*@C 47799fa7e03SVijay Mahadevan DMMoabGetDofsLocal - Gets the local degree-of-freedom for all fields (components) defined on an 47899fa7e03SVijay Mahadevan array of MOAB EntityHandles. 47999fa7e03SVijay Mahadevan 48099fa7e03SVijay Mahadevan Not Collective 48199fa7e03SVijay Mahadevan 48299fa7e03SVijay Mahadevan Input Parameters: 48399fa7e03SVijay Mahadevan + dm - the discretization manager object 48499fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 48599fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 48699fa7e03SVijay Mahadevan 48799fa7e03SVijay Mahadevan Output Parameter: 48899fa7e03SVijay Mahadevan + dof - the local degree-of-freedom index array corresponding to the field in the discrete representation (Vec, Mat) 48999fa7e03SVijay Mahadevan 49099fa7e03SVijay Mahadevan Level: intermediate 49199fa7e03SVijay Mahadevan 49299fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degrees of freedom 49399fa7e03SVijay Mahadevan 49499fa7e03SVijay Mahadevan .seealso: DMMoabGetFieldDofs(), DMMoabGetDofs(), DMMoabGetDofsBlocked() 49599fa7e03SVijay Mahadevan @*/ 4968cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetDofsLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 4978cbae1a6SVijay Mahadevan { 4988cbae1a6SVijay Mahadevan PetscInt i,field,offset; 4998cbae1a6SVijay Mahadevan PetscErrorCode ierr; 5008cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 5018cbae1a6SVijay Mahadevan 5028cbae1a6SVijay Mahadevan PetscFunctionBegin; 5038cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 50499fa7e03SVijay Mahadevan PetscValidPointer(points,3); 5058cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 5068cbae1a6SVijay Mahadevan 5078cbae1a6SVijay Mahadevan if (!dof) { 508addae81cSVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*dmmoab->numFields*npoints, &dof);CHKERRQ(ierr); 5098cbae1a6SVijay Mahadevan } 5108cbae1a6SVijay Mahadevan 51199fa7e03SVijay Mahadevan if (dmmoab->bs > 1) { /* compute the DOF based on local blocking in the fields */ 512addae81cSVijay Mahadevan for (field=0; field<dmmoab->numFields; ++field) { 5138cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 514addae81cSVijay Mahadevan dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]*dmmoab->numFields+field; 5158cbae1a6SVijay Mahadevan } 5168cbae1a6SVijay Mahadevan } 5178cbae1a6SVijay Mahadevan else { 51899fa7e03SVijay Mahadevan /* assume all fields have equal distribution; i.e., all fields are either defined on vertices or elements and not on a mixture */ 51999fa7e03SVijay Mahadevan /* TODO: eliminate the limitation using PetscSection to manage DOFs */ 520addae81cSVijay Mahadevan for (field=0; field<dmmoab->numFields; ++field) { 52199fa7e03SVijay Mahadevan offset = field*dmmoab->n; 5228cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 523addae81cSVijay Mahadevan dof[i*dmmoab->numFields+field] = dmmoab->lidmap[(PetscInt)points[i]]+offset; 5248cbae1a6SVijay Mahadevan } 5258cbae1a6SVijay Mahadevan } 5268cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 5278cbae1a6SVijay Mahadevan } 5288cbae1a6SVijay Mahadevan 5298cbae1a6SVijay Mahadevan 5308cbae1a6SVijay Mahadevan #undef __FUNCT__ 5318cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsBlocked" 53299fa7e03SVijay Mahadevan /*@C 53399fa7e03SVijay Mahadevan DMMoabGetDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an 53499fa7e03SVijay Mahadevan array of MOAB EntityHandles. It is useful when performing Blocked(Get/Set) methods in computation 53599fa7e03SVijay Mahadevan of element residuals and assembly of the discrete systems when all fields are co-located. 53699fa7e03SVijay Mahadevan 53799fa7e03SVijay Mahadevan Not Collective 53899fa7e03SVijay Mahadevan 53999fa7e03SVijay Mahadevan Input Parameters: 54099fa7e03SVijay Mahadevan + dm - the discretization manager object 54199fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 54299fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 54399fa7e03SVijay Mahadevan 54499fa7e03SVijay Mahadevan Output Parameter: 54599fa7e03SVijay Mahadevan + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) 54699fa7e03SVijay Mahadevan 54799fa7e03SVijay Mahadevan Level: intermediate 54899fa7e03SVijay Mahadevan 54999fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degrees of freedom 55099fa7e03SVijay Mahadevan 55199fa7e03SVijay Mahadevan .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal() 55299fa7e03SVijay Mahadevan @*/ 5538cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetDofsBlocked(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 5548cbae1a6SVijay Mahadevan { 5558cbae1a6SVijay Mahadevan PetscInt i; 5568cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 5578cbae1a6SVijay Mahadevan PetscErrorCode ierr; 5588cbae1a6SVijay Mahadevan 5598cbae1a6SVijay Mahadevan PetscFunctionBegin; 5608cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5618cbae1a6SVijay Mahadevan PetscValidPointer(points,2); 5628cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 5638cbae1a6SVijay Mahadevan 5648cbae1a6SVijay Mahadevan if (!dof) { 5658cbae1a6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 5668cbae1a6SVijay Mahadevan } 5678cbae1a6SVijay Mahadevan 5688cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) { 5698cbae1a6SVijay Mahadevan dof[i]=dmmoab->gidmap[(PetscInt)points[i]]; 5708cbae1a6SVijay Mahadevan } 5718cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 5728cbae1a6SVijay Mahadevan } 5738cbae1a6SVijay Mahadevan 5748cbae1a6SVijay Mahadevan 5758cbae1a6SVijay Mahadevan #undef __FUNCT__ 5768cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetDofsBlockedLocal" 57799fa7e03SVijay Mahadevan /*@C 57899fa7e03SVijay Mahadevan DMMoabGetDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an 57999fa7e03SVijay Mahadevan array of MOAB EntityHandles. It is useful when performing local Blocked(Get/Set) methods in computation 58099fa7e03SVijay Mahadevan of element residuals and assembly of the discrete systems when all fields are co-located. 58199fa7e03SVijay Mahadevan 58299fa7e03SVijay Mahadevan Not Collective 58399fa7e03SVijay Mahadevan 58499fa7e03SVijay Mahadevan Input Parameters: 58599fa7e03SVijay Mahadevan + dm - the discretization manager object 58699fa7e03SVijay Mahadevan . npoints - the total number of Entities in the points array 58799fa7e03SVijay Mahadevan . points - the MOAB EntityHandle container array which holds the field degree-of-freedom values 58899fa7e03SVijay Mahadevan 58999fa7e03SVijay Mahadevan Output Parameter: 59099fa7e03SVijay Mahadevan + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) 59199fa7e03SVijay Mahadevan 59299fa7e03SVijay Mahadevan Level: intermediate 59399fa7e03SVijay Mahadevan 59499fa7e03SVijay Mahadevan .keywords: discretization manager, get, global degrees of freedom 59599fa7e03SVijay Mahadevan 59699fa7e03SVijay Mahadevan .seealso: DMMoabGetDofsLocal(), DMMoabGetDofs(), DMMoabGetDofsBlockedLocal() 59799fa7e03SVijay Mahadevan @*/ 5988cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetDofsBlockedLocal(DM dm,PetscInt npoints,const moab::EntityHandle* points,PetscInt* dof) 5998cbae1a6SVijay Mahadevan { 6008cbae1a6SVijay Mahadevan PetscInt i; 6018cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 6028cbae1a6SVijay Mahadevan PetscErrorCode ierr; 6038cbae1a6SVijay Mahadevan 6048cbae1a6SVijay Mahadevan PetscFunctionBegin; 6058cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6068cbae1a6SVijay Mahadevan PetscValidPointer(points,2); 6078cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 6088cbae1a6SVijay Mahadevan 6098cbae1a6SVijay Mahadevan if (!dof) { 6108cbae1a6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*npoints, &dof);CHKERRQ(ierr); 6118cbae1a6SVijay Mahadevan } 6128cbae1a6SVijay Mahadevan 6138cbae1a6SVijay Mahadevan for (i=0; i<npoints; ++i) 6148cbae1a6SVijay Mahadevan dof[i] = dmmoab->lidmap[(PetscInt)points[i]]; 6158cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 6168cbae1a6SVijay Mahadevan } 6178cbae1a6SVijay Mahadevan 6188cbae1a6SVijay Mahadevan 6198cbae1a6SVijay Mahadevan #undef __FUNCT__ 6208cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexDofsBlocked" 62199fa7e03SVijay Mahadevan /*@C 62299fa7e03SVijay Mahadevan DMMoabGetVertexDofsBlocked - Gets the global degree-of-freedom for the first field (component) defined on an 62399fa7e03SVijay Mahadevan array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations 62499fa7e03SVijay Mahadevan where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations. 62599fa7e03SVijay Mahadevan 62699fa7e03SVijay Mahadevan Not Collective 62799fa7e03SVijay Mahadevan 62899fa7e03SVijay Mahadevan Input Parameters: 62999fa7e03SVijay Mahadevan + dm - the discretization manager object 63099fa7e03SVijay Mahadevan 63199fa7e03SVijay Mahadevan Output Parameter: 63299fa7e03SVijay Mahadevan + dof - the blocked global degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering 63399fa7e03SVijay Mahadevan 63499fa7e03SVijay Mahadevan Level: intermediate 63599fa7e03SVijay Mahadevan 63699fa7e03SVijay Mahadevan .keywords: discretization manager, get, blocked degrees of freedom 63799fa7e03SVijay Mahadevan 63899fa7e03SVijay Mahadevan .seealso: DMMoabGetVertexDofsBlockedLocal(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal() 63999fa7e03SVijay Mahadevan @*/ 6408cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetVertexDofsBlocked(DM dm,PetscInt** dof) 6418cbae1a6SVijay Mahadevan { 6428cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 6438cbae1a6SVijay Mahadevan 6448cbae1a6SVijay Mahadevan PetscFunctionBegin; 6458cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6468cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 6478cbae1a6SVijay Mahadevan 6488cbae1a6SVijay Mahadevan *dof = dmmoab->gidmap; 6498cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 6508cbae1a6SVijay Mahadevan } 6518cbae1a6SVijay Mahadevan 6528cbae1a6SVijay Mahadevan 6538cbae1a6SVijay Mahadevan #undef __FUNCT__ 6548cbae1a6SVijay Mahadevan #define __FUNCT__ "DMMoabGetVertexDofsBlockedLocal" 65599fa7e03SVijay Mahadevan /*@C 65699fa7e03SVijay Mahadevan DMMoabGetVertexDofsBlockedLocal - Gets the local degree-of-freedom for the first field (component) defined on an 65799fa7e03SVijay Mahadevan array of locally owned MOAB mesh vertices. It's utility is when performing Finite-Difference type calculations 65899fa7e03SVijay Mahadevan where vertex traversal is faster than element-wise assembly that is typically done in FEM calculations. 65999fa7e03SVijay Mahadevan 66099fa7e03SVijay Mahadevan Not Collective 66199fa7e03SVijay Mahadevan 66299fa7e03SVijay Mahadevan Input Parameters: 66399fa7e03SVijay Mahadevan + dm - the discretization manager object 66499fa7e03SVijay Mahadevan 66599fa7e03SVijay Mahadevan Output Parameter: 66699fa7e03SVijay Mahadevan + dof - the blocked local degree-of-freedom index array in the discrete representation (Vec, Mat) that is vertex-based based on local numbering 66799fa7e03SVijay Mahadevan 66899fa7e03SVijay Mahadevan Level: intermediate 66999fa7e03SVijay Mahadevan 67099fa7e03SVijay Mahadevan .keywords: discretization manager, get, blocked degrees of freedom 67199fa7e03SVijay Mahadevan 67299fa7e03SVijay Mahadevan .seealso: DMMoabGetVertexDofsBlocked(), DMMoabGetDofsBlocked(), DMMoabGetDofsBlockedLocal() 67399fa7e03SVijay Mahadevan @*/ 6748cbae1a6SVijay Mahadevan PetscErrorCode DMMoabGetVertexDofsBlockedLocal(DM dm,PetscInt** dof) 6758cbae1a6SVijay Mahadevan { 6768cbae1a6SVijay Mahadevan DM_Moab *dmmoab; 6778cbae1a6SVijay Mahadevan 6788cbae1a6SVijay Mahadevan PetscFunctionBegin; 6798cbae1a6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 6808cbae1a6SVijay Mahadevan PetscValidPointer(dof,2); 6818cbae1a6SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 6828cbae1a6SVijay Mahadevan 6838cbae1a6SVijay Mahadevan *dof = dmmoab->lidmap; 6848cbae1a6SVijay Mahadevan PetscFunctionReturn(0); 6858cbae1a6SVijay Mahadevan } 6868cbae1a6SVijay Mahadevan 687