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