xref: /petsc/src/dm/impls/plex/tests/ex10.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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 
14c4762a1bSJed Brown PetscErrorCode ProcessOptions(AppCtx *options)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown   PetscInt       len;
17c4762a1bSJed Brown   PetscBool      flg;
18c4762a1bSJed Brown   PetscErrorCode ierr;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBegin;
21c4762a1bSJed Brown   options->numFields     = 1;
22c4762a1bSJed Brown   options->numComponents = NULL;
23c4762a1bSJed Brown   options->numDof        = NULL;
24c4762a1bSJed Brown   options->numGroups     = 0;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex10.c", options->numFields, &options->numFields, NULL,1));
28c4762a1bSJed Brown   if (options->numFields) {
29c4762a1bSJed Brown     len  = options->numFields;
305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(len, &options->numComponents));
315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsIntArray("-num_components", "The number of components per field", "ex10.c", options->numComponents, &len, &flg));
322c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg && (len != options->numFields),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %D should be %D", len, options->numFields);
33c4762a1bSJed Brown   }
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-num_groups", "Group permutation by this many label values", "ex10.c", options->numGroups, &options->numGroups, NULL,0));
35c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
36c4762a1bSJed Brown   PetscFunctionReturn(0);
37c4762a1bSJed Brown }
38c4762a1bSJed Brown 
39c4762a1bSJed Brown PetscErrorCode CleanupContext(AppCtx *user)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   PetscFunctionBegin;
425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->numComponents));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(user->numDof));
44c4762a1bSJed Brown   PetscFunctionReturn(0);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /* This mesh comes from~\cite{saad2003}, Fig. 2.10, p. 70. */
48c4762a1bSJed Brown PetscErrorCode CreateTestMesh(MPI_Comm comm, DM *dm, AppCtx *options)
49c4762a1bSJed Brown {
50a4a685f2SJacob Faibussowitsch   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,
51c4762a1bSJed Brown                                     2, 9, 7,  10, 9,  4,   1, 10, 12,  10,  4, 11,  12, 11, 3,  3, 11, 14,  11, 4, 13,  14, 13,  5};
52a4a685f2SJacob Faibussowitsch   const PetscReal   coords[15*2] = {0, -3,  0, -1,  2, -1,  0,  1,  2, 1,
53c4762a1bSJed Brown                                     0,  3,  1, -2,  1, -1,  0, -2,  2, 0,
54c4762a1bSJed Brown                                     1,  0,  1,  1,  0,  0,  1,  2,  0, 2};
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   PetscFunctionBegin;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateFromCellListPetsc(comm, 2, 16, 15, 3, PETSC_FALSE, cells, 2, coords, dm));
58c4762a1bSJed Brown   PetscFunctionReturn(0);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown PetscErrorCode TestReordering(DM dm, AppCtx *user)
62c4762a1bSJed Brown {
63c4762a1bSJed Brown   DM              pdm;
64c4762a1bSJed Brown   IS              perm;
65c4762a1bSJed Brown   Mat             A, pA;
66c4762a1bSJed Brown   PetscInt        bw, pbw;
67c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   PetscFunctionBegin;
705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetOrdering(dm, order, NULL, &perm));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexPermute(dm, perm, &pdm));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_"));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(pdm));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm,  NULL, "-orig_dm_view"));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(pdm, NULL, "-dm_view"));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &A));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(pdm, &pA));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeBandwidth(A, 0.0, &bw));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeBandwidth(pA, 0.0, &pbw));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,  NULL, "-orig_mat_view"));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&pA));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&pdm));
86c4762a1bSJed Brown   if (pbw > bw) {
875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s increased bandwidth from %D to %D\n", order, bw, pbw));
88c4762a1bSJed Brown   } else {
895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s reduced bandwidth from %D to %D\n", order, bw, pbw));
90c4762a1bSJed Brown   }
91c4762a1bSJed Brown   PetscFunctionReturn(0);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown PetscErrorCode CreateGroupLabel(DM dm, PetscInt numGroups, DMLabel *label, AppCtx *options)
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   const PetscInt groupA[10] = {15, 3, 13, 12, 2, 10, 7, 6, 0, 4};
97c4762a1bSJed Brown   const PetscInt groupB[6]  = {14, 11, 9, 1, 8, 5};
98c4762a1bSJed Brown   PetscInt       c;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBegin;
101c4762a1bSJed Brown   if (numGroups < 2) {*label = NULL; PetscFunctionReturn(0);}
1022c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numGroups != 2,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Test only coded for 2 groups, not %D", numGroups);
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "groups", label));
1045f80ce2aSJacob Faibussowitsch   for (c = 0; c < 10; ++c) CHKERRQ(DMLabelSetValue(*label, groupA[c], 101));
1055f80ce2aSJacob Faibussowitsch   for (c = 0; c < 6;  ++c) CHKERRQ(DMLabelSetValue(*label, groupB[c], 1001));
106c4762a1bSJed Brown   PetscFunctionReturn(0);
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109c4762a1bSJed Brown PetscErrorCode TestReorderingByGroup(DM dm, AppCtx *user)
110c4762a1bSJed Brown {
111c4762a1bSJed Brown   DM              pdm;
112c4762a1bSJed Brown   DMLabel         label;
113c4762a1bSJed Brown   Mat             A, pA;
114c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
115c4762a1bSJed Brown   IS              perm;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBegin;
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateGroupLabel(dm, user->numGroups, &label, user));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetOrdering(dm, order, label, &perm));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&label));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexPermute(dm, perm, &pdm));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_"));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(pdm));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm,  NULL, "-orig_dm_view"));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(pdm, NULL, "-perm_dm_view"));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(dm, &A));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(pdm, &pA));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(A,  NULL, "-orig_mat_view"));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&pA));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&pdm));
134c4762a1bSJed Brown   PetscFunctionReturn(0);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown int main(int argc, char **argv)
138c4762a1bSJed Brown {
139c4762a1bSJed Brown   DM             dm;
140c4762a1bSJed Brown   PetscSection   s;
141c4762a1bSJed Brown   AppCtx         user;
14230602db0SMatthew G. Knepley   PetscInt       dim;
143c4762a1bSJed Brown   PetscErrorCode ierr;
144c4762a1bSJed Brown 
145*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(&user));
147c4762a1bSJed Brown   if (user.numGroups < 1) {
1485f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
1495f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetType(dm, DMPLEX));
15030602db0SMatthew G. Knepley   } else {
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(CreateTestMesh(PETSC_COMM_WORLD, &dm, &user));
1520fdc7489SMatthew Knepley   }
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
15630602db0SMatthew G. Knepley   {
15730602db0SMatthew G. Knepley     PetscInt  len = (dim+1) * PetscMax(1, user.numFields);
15830602db0SMatthew G. Knepley     PetscBool flg;
15930602db0SMatthew G. Knepley 
1605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(len, &user.numDof));
16130602db0SMatthew G. Knepley     ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
1625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex10.c", user.numDof, &len, &flg));
1632c71b3e2SJacob Faibussowitsch     PetscCheckFalse(flg && (len != (dim+1) * PetscMax(1, user.numFields)),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %D should be %D", len, (dim+1) * PetscMax(1, user.numFields));
16430602db0SMatthew G. Knepley     ierr = PetscOptionsEnd();CHKERRQ(ierr);
16530602db0SMatthew G. Knepley   }
16630602db0SMatthew G. Knepley   if (user.numGroups < 1) {
1675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetNumFields(dm, user.numFields));
1685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateDS(dm));
1695f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s));
1705f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetLocalSection(dm, s));
1715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionDestroy(&s));
1725f80ce2aSJacob Faibussowitsch     CHKERRQ(TestReordering(dm, &user));
173c4762a1bSJed Brown   } else {
1745f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetNumFields(dm, user.numFields));
1755f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateDS(dm));
1765f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s));
1775f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetLocalSection(dm, s));
1785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionDestroy(&s));
1795f80ce2aSJacob Faibussowitsch     CHKERRQ(TestReorderingByGroup(dm, &user));
180c4762a1bSJed Brown   }
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(CleanupContext(&user));
183*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
184*b122ec5aSJacob Faibussowitsch   return 0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown /*TEST
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   # Two cell tests 0-3
190c4762a1bSJed Brown   test:
191c4762a1bSJed Brown     suffix: 0
1920fdc7489SMatthew Knepley     requires: triangle
19330602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -num_dof 1,0,0 -mat_view -dm_coord_space 0
194c4762a1bSJed Brown   test:
195c4762a1bSJed Brown     suffix: 1
19630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -num_dof 1,0,0 -mat_view -dm_coord_space 0
197c4762a1bSJed Brown   test:
198c4762a1bSJed Brown     suffix: 2
1990fdc7489SMatthew Knepley     requires: ctetgen
20030602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
201c4762a1bSJed Brown   test:
202c4762a1bSJed Brown     suffix: 3
20330602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
204c4762a1bSJed Brown   # Refined tests 4-7
205c4762a1bSJed Brown   test:
206c4762a1bSJed Brown     suffix: 4
207c4762a1bSJed Brown     requires: triangle
20830602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0
209c4762a1bSJed Brown   test:
210c4762a1bSJed Brown     suffix: 5
21130602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0
212c4762a1bSJed Brown   test:
213c4762a1bSJed Brown     suffix: 6
214c4762a1bSJed Brown     requires: ctetgen
21530602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0,0
216c4762a1bSJed Brown   test:
217c4762a1bSJed Brown     suffix: 7
21830602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0,0
219c4762a1bSJed Brown   # Parallel tests
220c4762a1bSJed Brown   # Grouping tests
221c4762a1bSJed Brown   test:
222c4762a1bSJed Brown     suffix: group_1
223c4762a1bSJed Brown     args: -num_groups 1 -num_dof 1,0,0 -is_view -orig_mat_view -perm_mat_view
224c4762a1bSJed Brown   test:
225c4762a1bSJed Brown     suffix: group_2
226c4762a1bSJed Brown     args: -num_groups 2 -num_dof 1,0,0 -is_view -perm_mat_view
227c4762a1bSJed Brown 
228c4762a1bSJed Brown TEST*/
229