xref: /petsc/src/dm/impls/plex/tests/ex10.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Test for mesh reordering\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   PetscInt  dim;             /* The topological mesh dimension */
7c4762a1bSJed Brown   PetscReal refinementLimit; /* Maximum volume of a refined cell */
8c4762a1bSJed Brown   PetscInt  numFields;       /* The number of section fields */
9c4762a1bSJed Brown   PetscInt *numComponents;   /* The number of field components */
10c4762a1bSJed Brown   PetscInt *numDof;          /* The dof signature for the section */
11c4762a1bSJed Brown   PetscInt  numGroups;       /* If greater than 1, use grouping in test */
12c4762a1bSJed Brown } AppCtx;
13c4762a1bSJed Brown 
14*9371c9d4SSatish Balay PetscErrorCode ProcessOptions(AppCtx *options) {
15c4762a1bSJed Brown   PetscInt  len;
16c4762a1bSJed Brown   PetscBool flg;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   PetscFunctionBegin;
19c4762a1bSJed Brown   options->numFields     = 1;
20c4762a1bSJed Brown   options->numComponents = NULL;
21c4762a1bSJed Brown   options->numDof        = NULL;
22c4762a1bSJed Brown   options->numGroups     = 0;
23c4762a1bSJed Brown 
24d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex10.c", options->numFields, &options->numFields, NULL, 1));
26c4762a1bSJed Brown   if (options->numFields) {
27c4762a1bSJed Brown     len = options->numFields;
289566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(len, &options->numComponents));
299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex10.c", options->numComponents, &len, &flg));
3063a3b9bcSJacob Faibussowitsch     PetscCheck(!flg || !(len != options->numFields), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %" PetscInt_FMT " should be %" PetscInt_FMT, len, options->numFields);
31c4762a1bSJed Brown   }
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_groups", "Group permutation by this many label values", "ex10.c", options->numGroups, &options->numGroups, NULL, 0));
33d0609cedSBarry Smith   PetscOptionsEnd();
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37*9371c9d4SSatish Balay PetscErrorCode CleanupContext(AppCtx *user) {
38c4762a1bSJed Brown   PetscFunctionBegin;
399566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numComponents));
409566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numDof));
41c4762a1bSJed Brown   PetscFunctionReturn(0);
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /* This mesh comes from~\cite{saad2003}, Fig. 2.10, p. 70. */
45*9371c9d4SSatish Balay PetscErrorCode CreateTestMesh(MPI_Comm comm, DM *dm, AppCtx *options) {
46*9371c9d4SSatish Balay   const PetscInt  cells[16 * 3]  = {6, 7, 8, 7, 9, 10, 10, 11, 12, 11, 13, 14, 0, 6, 8, 6, 2, 7, 1, 8, 7, 1, 7, 10, 2, 9, 7, 10, 9, 4, 1, 10, 12, 10, 4, 11, 12, 11, 3, 3, 11, 14, 11, 4, 13, 14, 13, 5};
47*9371c9d4SSatish Balay   const PetscReal coords[15 * 2] = {0, -3, 0, -1, 2, -1, 0, 1, 2, 1, 0, 3, 1, -2, 1, -1, 0, -2, 2, 0, 1, 0, 1, 1, 0, 0, 1, 2, 0, 2};
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListPetsc(comm, 2, 16, 15, 3, PETSC_FALSE, cells, 2, coords, dm));
51c4762a1bSJed Brown   PetscFunctionReturn(0);
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54*9371c9d4SSatish Balay PetscErrorCode TestReordering(DM dm, AppCtx *user) {
55c4762a1bSJed Brown   DM              pdm;
56c4762a1bSJed Brown   IS              perm;
57c4762a1bSJed Brown   Mat             A, pA;
58c4762a1bSJed Brown   PetscInt        bw, pbw;
59c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscFunctionBegin;
629566063dSJacob Faibussowitsch   PetscCall(DMPlexGetOrdering(dm, order, NULL, &perm));
639566063dSJacob Faibussowitsch   PetscCall(DMPlexPermute(dm, perm, &pdm));
649566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pdm, "perm_"));
659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pdm));
669566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
679566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view"));
689566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
699566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &A));
709566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pdm, &pA));
719566063dSJacob Faibussowitsch   PetscCall(MatComputeBandwidth(A, 0.0, &bw));
729566063dSJacob Faibussowitsch   PetscCall(MatComputeBandwidth(pA, 0.0, &pbw));
739566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-orig_mat_view"));
749566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pA));
779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pdm));
78c4762a1bSJed Brown   if (pbw > bw) {
7963a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Ordering method %s increased bandwidth from %" PetscInt_FMT " to %" PetscInt_FMT "\n", order, bw, pbw));
80c4762a1bSJed Brown   } else {
8163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Ordering method %s reduced bandwidth from %" PetscInt_FMT " to %" PetscInt_FMT "\n", order, bw, pbw));
82c4762a1bSJed Brown   }
83c4762a1bSJed Brown   PetscFunctionReturn(0);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86*9371c9d4SSatish Balay PetscErrorCode CreateGroupLabel(DM dm, PetscInt numGroups, DMLabel *label, AppCtx *options) {
87c4762a1bSJed Brown   const PetscInt groupA[10] = {15, 3, 13, 12, 2, 10, 7, 6, 0, 4};
88c4762a1bSJed Brown   const PetscInt groupB[6]  = {14, 11, 9, 1, 8, 5};
89c4762a1bSJed Brown   PetscInt       c;
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   PetscFunctionBegin;
92*9371c9d4SSatish Balay   if (numGroups < 2) {
93*9371c9d4SSatish Balay     *label = NULL;
94*9371c9d4SSatish Balay     PetscFunctionReturn(0);
95*9371c9d4SSatish Balay   }
9663a3b9bcSJacob Faibussowitsch   PetscCheck(numGroups == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Test only coded for 2 groups, not %" PetscInt_FMT, numGroups);
979566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "groups", label));
989566063dSJacob Faibussowitsch   for (c = 0; c < 10; ++c) PetscCall(DMLabelSetValue(*label, groupA[c], 101));
999566063dSJacob Faibussowitsch   for (c = 0; c < 6; ++c) PetscCall(DMLabelSetValue(*label, groupB[c], 1001));
100c4762a1bSJed Brown   PetscFunctionReturn(0);
101c4762a1bSJed Brown }
102c4762a1bSJed Brown 
103*9371c9d4SSatish Balay PetscErrorCode TestReorderingByGroup(DM dm, AppCtx *user) {
104c4762a1bSJed Brown   DM              pdm;
105c4762a1bSJed Brown   DMLabel         label;
106c4762a1bSJed Brown   Mat             A, pA;
107c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
108c4762a1bSJed Brown   IS              perm;
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   PetscFunctionBegin;
1119566063dSJacob Faibussowitsch   PetscCall(CreateGroupLabel(dm, user->numGroups, &label, user));
1129566063dSJacob Faibussowitsch   PetscCall(DMPlexGetOrdering(dm, order, label, &perm));
1139566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
1149566063dSJacob Faibussowitsch   PetscCall(DMPlexPermute(dm, perm, &pdm));
1159566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pdm, "perm_"));
1169566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pdm));
1179566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-orig_dm_view"));
1189566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(pdm, NULL, "-perm_dm_view"));
1199566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1209566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &A));
1219566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pdm, &pA));
1229566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A, NULL, "-orig_mat_view"));
1239566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
1249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pA));
1269566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pdm));
127c4762a1bSJed Brown   PetscFunctionReturn(0);
128c4762a1bSJed Brown }
129c4762a1bSJed Brown 
130*9371c9d4SSatish Balay int main(int argc, char **argv) {
131c4762a1bSJed Brown   DM           dm;
132c4762a1bSJed Brown   PetscSection s;
133c4762a1bSJed Brown   AppCtx       user;
13430602db0SMatthew G. Knepley   PetscInt     dim;
135c4762a1bSJed Brown 
136327415f7SBarry Smith   PetscFunctionBeginUser;
1379566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1389566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
139c4762a1bSJed Brown   if (user.numGroups < 1) {
1409566063dSJacob Faibussowitsch     PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
1419566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
14230602db0SMatthew G. Knepley   } else {
1439566063dSJacob Faibussowitsch     PetscCall(CreateTestMesh(PETSC_COMM_WORLD, &dm, &user));
1440fdc7489SMatthew Knepley   }
1459566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1469566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1479566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
14830602db0SMatthew G. Knepley   {
14930602db0SMatthew G. Knepley     PetscInt  len = (dim + 1) * PetscMax(1, user.numFields);
15030602db0SMatthew G. Knepley     PetscBool flg;
15130602db0SMatthew G. Knepley 
1529566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(len, &user.numDof));
153d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
1549566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex10.c", user.numDof, &len, &flg));
15563a3b9bcSJacob Faibussowitsch     if (flg) PetscCheck(len == ((dim + 1) * PetscMax(1, user.numFields)), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %" PetscInt_FMT " should be %" PetscInt_FMT, len, (dim + 1) * PetscMax(1, user.numFields));
156d0609cedSBarry Smith     PetscOptionsEnd();
15730602db0SMatthew G. Knepley   }
15830602db0SMatthew G. Knepley   if (user.numGroups < 1) {
1599566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(dm, user.numFields));
1609566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(dm));
1619566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s));
1629566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(dm, s));
1639566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&s));
1649566063dSJacob Faibussowitsch     PetscCall(TestReordering(dm, &user));
165c4762a1bSJed Brown   } else {
1669566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(dm, user.numFields));
1679566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(dm));
1689566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s));
1699566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(dm, s));
1709566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&s));
1719566063dSJacob Faibussowitsch     PetscCall(TestReorderingByGroup(dm, &user));
172c4762a1bSJed Brown   }
1739566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1749566063dSJacob Faibussowitsch   PetscCall(CleanupContext(&user));
1759566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
176b122ec5aSJacob Faibussowitsch   return 0;
177c4762a1bSJed Brown }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown /*TEST
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   # Two cell tests 0-3
182c4762a1bSJed Brown   test:
183c4762a1bSJed Brown     suffix: 0
1840fdc7489SMatthew Knepley     requires: triangle
18530602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -num_dof 1,0,0 -mat_view -dm_coord_space 0
186c4762a1bSJed Brown   test:
187c4762a1bSJed Brown     suffix: 1
18830602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -num_dof 1,0,0 -mat_view -dm_coord_space 0
189c4762a1bSJed Brown   test:
190c4762a1bSJed Brown     suffix: 2
1910fdc7489SMatthew Knepley     requires: ctetgen
19230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
193c4762a1bSJed Brown   test:
194c4762a1bSJed Brown     suffix: 3
19530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
196c4762a1bSJed Brown   # Refined tests 4-7
197c4762a1bSJed Brown   test:
198c4762a1bSJed Brown     suffix: 4
199c4762a1bSJed Brown     requires: triangle
20030602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0
201c4762a1bSJed Brown   test:
202c4762a1bSJed Brown     suffix: 5
20330602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0
204c4762a1bSJed Brown   test:
205c4762a1bSJed Brown     suffix: 6
206c4762a1bSJed Brown     requires: ctetgen
20730602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0,0
208c4762a1bSJed Brown   test:
209c4762a1bSJed Brown     suffix: 7
21030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0,0
211c4762a1bSJed Brown   # Parallel tests
212c4762a1bSJed Brown   # Grouping tests
213c4762a1bSJed Brown   test:
214c4762a1bSJed Brown     suffix: group_1
215c4762a1bSJed Brown     args: -num_groups 1 -num_dof 1,0,0 -is_view -orig_mat_view -perm_mat_view
216c4762a1bSJed Brown   test:
217c4762a1bSJed Brown     suffix: group_2
218c4762a1bSJed Brown     args: -num_groups 2 -num_dof 1,0,0 -is_view -perm_mat_view
219c4762a1bSJed Brown 
220c4762a1bSJed Brown TEST*/
221