1*fa58d8a1SBlaise Bourdin static char help[] = "Test FEM layout with constraints\n\n"; 2*fa58d8a1SBlaise Bourdin 3*fa58d8a1SBlaise Bourdin #include <petsc.h> 4*fa58d8a1SBlaise Bourdin 5*fa58d8a1SBlaise Bourdin int main(int argc,char **argv) { 6*fa58d8a1SBlaise Bourdin DM dm,pdm; 7*fa58d8a1SBlaise Bourdin PetscSection section; 8*fa58d8a1SBlaise Bourdin const PetscInt field = 0; 9*fa58d8a1SBlaise Bourdin char ifilename[PETSC_MAX_PATH_LEN]; 10*fa58d8a1SBlaise Bourdin PetscInt sdim,s,pStart,pEnd,p,numVS,numPoints; 11*fa58d8a1SBlaise Bourdin PetscInt constraints[1]; 12*fa58d8a1SBlaise Bourdin IS setIS,pointIS; 13*fa58d8a1SBlaise Bourdin const PetscInt *setID,*pointID; 14*fa58d8a1SBlaise Bourdin 15*fa58d8a1SBlaise Bourdin PetscFunctionBeginUser; 16*fa58d8a1SBlaise Bourdin PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 17*fa58d8a1SBlaise Bourdin PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"FEM Layout Options","ex98"); 18*fa58d8a1SBlaise Bourdin PetscCall(PetscOptionsString("-i","Filename to read","ex98",ifilename,ifilename,sizeof(ifilename),NULL)); 19*fa58d8a1SBlaise Bourdin PetscOptionsEnd(); 20*fa58d8a1SBlaise Bourdin 21*fa58d8a1SBlaise Bourdin PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD,ifilename,NULL,PETSC_TRUE,&dm)); 22*fa58d8a1SBlaise Bourdin PetscCall(DMPlexDistributeSetDefault(dm,PETSC_FALSE)); 23*fa58d8a1SBlaise Bourdin PetscCall(DMSetFromOptions(dm)); 24*fa58d8a1SBlaise Bourdin 25*fa58d8a1SBlaise Bourdin PetscCall(DMPlexDistribute(dm,0,NULL,&pdm)); 26*fa58d8a1SBlaise Bourdin if (pdm) { 27*fa58d8a1SBlaise Bourdin PetscCall(DMDestroy(&dm)); 28*fa58d8a1SBlaise Bourdin dm = pdm; 29*fa58d8a1SBlaise Bourdin } 30*fa58d8a1SBlaise Bourdin PetscCall(DMViewFromOptions(dm,NULL,"-dm_view")); 31*fa58d8a1SBlaise Bourdin 32*fa58d8a1SBlaise Bourdin /* create a section */ 33*fa58d8a1SBlaise Bourdin PetscCall(DMGetDimension(dm,&sdim)); 34*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) dm),§ion)); 35*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetNumFields(section,1)); 36*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetFieldName(section,field,"U")); 37*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetFieldComponents(section,field,sdim)); 38*fa58d8a1SBlaise Bourdin PetscCall(DMPlexGetChart(dm,&pStart,&pEnd)); 39*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetChart(section,pStart,pEnd)); 40*fa58d8a1SBlaise Bourdin 41*fa58d8a1SBlaise Bourdin /* initialize the section storage for a P1 field */ 42*fa58d8a1SBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm,0,&pStart,&pEnd)); 43*fa58d8a1SBlaise Bourdin for (p = pStart; p < pEnd; ++p) { 44*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetDof(section,p,sdim)); 45*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetFieldDof(section,p,0,sdim)); 46*fa58d8a1SBlaise Bourdin } 47*fa58d8a1SBlaise Bourdin 48*fa58d8a1SBlaise Bourdin /* add constraints at all vertices belonging to a vertex set */ 49*fa58d8a1SBlaise Bourdin /* first pass is to reserve space */ 50*fa58d8a1SBlaise Bourdin PetscCall(DMGetLabelSize(dm,"Vertex Sets",&numVS)); 51*fa58d8a1SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD,"# Vertex set: %" PetscInt_FMT "\n",numVS)); 52*fa58d8a1SBlaise Bourdin PetscCall(DMGetLabelIdIS(dm,"Vertex Sets",&setIS)); 53*fa58d8a1SBlaise Bourdin PetscCall(ISGetIndices(setIS,&setID)); 54*fa58d8a1SBlaise Bourdin for (s = 0; s < numVS; ++s){ 55*fa58d8a1SBlaise Bourdin PetscCall(DMGetStratumIS(dm,"Vertex Sets",setID[s],&pointIS)); 56*fa58d8a1SBlaise Bourdin PetscCall(DMGetStratumSize(dm,"Vertex Sets",setID[s],&numPoints)); 57*fa58d8a1SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD,"set %" PetscInt_FMT " size: %" PetscInt_FMT "\n",s,numPoints)); 58*fa58d8a1SBlaise Bourdin PetscCall(ISGetIndices(pointIS,&pointID)); 59*fa58d8a1SBlaise Bourdin for (p = 0; p < numPoints; ++p){ 60*fa58d8a1SBlaise Bourdin PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\t point %" PetscInt_FMT "\n",pointID[p])); 61*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetConstraintDof(section,pointID[p],1)); 62*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetFieldConstraintDof(section,pointID[p],field,1)); 63*fa58d8a1SBlaise Bourdin } 64*fa58d8a1SBlaise Bourdin PetscCall(ISRestoreIndices(pointIS,&pointID)); 65*fa58d8a1SBlaise Bourdin PetscCall(ISDestroy(&pointIS)); 66*fa58d8a1SBlaise Bourdin } 67*fa58d8a1SBlaise Bourdin 68*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetUp(section)); 69*fa58d8a1SBlaise Bourdin 70*fa58d8a1SBlaise Bourdin /* add constraints at all vertices belonging to a vertex set */ 71*fa58d8a1SBlaise Bourdin /* second pass is to assign constraints to a specific component / dof */ 72*fa58d8a1SBlaise Bourdin for (s = 0; s < numVS; ++s){ 73*fa58d8a1SBlaise Bourdin PetscCall(DMGetStratumIS(dm,"Vertex Sets",setID[s],&pointIS)); 74*fa58d8a1SBlaise Bourdin PetscCall(DMGetStratumSize(dm,"Vertex Sets",setID[s],&numPoints)); 75*fa58d8a1SBlaise Bourdin PetscCall(ISGetIndices(pointIS,&pointID)); 76*fa58d8a1SBlaise Bourdin for (p = 0; p < numPoints; ++p){ 77*fa58d8a1SBlaise Bourdin constraints[0] = setID[s]%sdim; 78*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetConstraintIndices(section,pointID[p],constraints)); 79*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionSetFieldConstraintIndices(section,pointID[p],field,constraints)); 80*fa58d8a1SBlaise Bourdin } 81*fa58d8a1SBlaise Bourdin PetscCall(ISRestoreIndices(pointIS,&pointID)); 82*fa58d8a1SBlaise Bourdin PetscCall(ISDestroy(&pointIS)); 83*fa58d8a1SBlaise Bourdin } 84*fa58d8a1SBlaise Bourdin PetscCall(ISRestoreIndices(setIS,&setID)); 85*fa58d8a1SBlaise Bourdin PetscCall(ISDestroy(&setIS)); 86*fa58d8a1SBlaise Bourdin PetscCall(PetscObjectViewFromOptions((PetscObject) section,NULL,"-dm_section_view")); 87*fa58d8a1SBlaise Bourdin 88*fa58d8a1SBlaise Bourdin PetscCall(PetscSectionDestroy(§ion)); 89*fa58d8a1SBlaise Bourdin PetscCall(DMDestroy(&dm)); 90*fa58d8a1SBlaise Bourdin 91*fa58d8a1SBlaise Bourdin PetscCall(PetscFinalize()); 92*fa58d8a1SBlaise Bourdin return 0; 93*fa58d8a1SBlaise Bourdin } 94*fa58d8a1SBlaise Bourdin 95*fa58d8a1SBlaise Bourdin /*TEST 96*fa58d8a1SBlaise Bourdin build: 97*fa58d8a1SBlaise Bourdin requires: exodusii pnetcdf !complex 98*fa58d8a1SBlaise Bourdin testset: 99*fa58d8a1SBlaise Bourdin args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/SquareFaceSet.exo -dm_view -dm_section_view 100*fa58d8a1SBlaise Bourdin nsize: 1 101*fa58d8a1SBlaise Bourdin 102*fa58d8a1SBlaise Bourdin test: 103*fa58d8a1SBlaise Bourdin suffix: 0 104*fa58d8a1SBlaise Bourdin args: 105*fa58d8a1SBlaise Bourdin 106*fa58d8a1SBlaise Bourdin TEST*/ 107