xref: /petsc/src/dm/impls/plex/tests/ex10.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
19c4762a1bSJed Brown   PetscFunctionBegin;
20c4762a1bSJed Brown   options->numFields     = 1;
21c4762a1bSJed Brown   options->numComponents = NULL;
22c4762a1bSJed Brown   options->numDof        = NULL;
23c4762a1bSJed Brown   options->numGroups     = 0;
24c4762a1bSJed Brown 
25d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex10.c", options->numFields, &options->numFields, NULL,1));
27c4762a1bSJed Brown   if (options->numFields) {
28c4762a1bSJed Brown     len  = options->numFields;
299566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(len, &options->numComponents));
309566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_components", "The number of components per field", "ex10.c", options->numComponents, &len, &flg));
3163a3b9bcSJacob 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);
32c4762a1bSJed Brown   }
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-num_groups", "Group permutation by this many label values", "ex10.c", options->numGroups, &options->numGroups, NULL,0));
34d0609cedSBarry Smith   PetscOptionsEnd();
35c4762a1bSJed Brown   PetscFunctionReturn(0);
36c4762a1bSJed Brown }
37c4762a1bSJed Brown 
38c4762a1bSJed Brown PetscErrorCode CleanupContext(AppCtx *user)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   PetscFunctionBegin;
419566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numComponents));
429566063dSJacob Faibussowitsch   PetscCall(PetscFree(user->numDof));
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* This mesh comes from~\cite{saad2003}, Fig. 2.10, p. 70. */
47c4762a1bSJed Brown PetscErrorCode CreateTestMesh(MPI_Comm comm, DM *dm, AppCtx *options)
48c4762a1bSJed Brown {
49a4a685f2SJacob 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,
50c4762a1bSJed 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};
51a4a685f2SJacob Faibussowitsch   const PetscReal   coords[15*2] = {0, -3,  0, -1,  2, -1,  0,  1,  2, 1,
52c4762a1bSJed Brown                                     0,  3,  1, -2,  1, -1,  0, -2,  2, 0,
53c4762a1bSJed Brown                                     1,  0,  1,  1,  0,  0,  1,  2,  0, 2};
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListPetsc(comm, 2, 16, 15, 3, PETSC_FALSE, cells, 2, coords, dm));
57c4762a1bSJed Brown   PetscFunctionReturn(0);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown PetscErrorCode TestReordering(DM dm, AppCtx *user)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   DM              pdm;
63c4762a1bSJed Brown   IS              perm;
64c4762a1bSJed Brown   Mat             A, pA;
65c4762a1bSJed Brown   PetscInt        bw, pbw;
66c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(DMPlexGetOrdering(dm, order, NULL, &perm));
709566063dSJacob Faibussowitsch   PetscCall(DMPlexPermute(dm, perm, &pdm));
719566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_"));
729566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pdm));
739566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
749566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm,  NULL, "-orig_dm_view"));
759566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
769566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &A));
779566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pdm, &pA));
789566063dSJacob Faibussowitsch   PetscCall(MatComputeBandwidth(A, 0.0, &bw));
799566063dSJacob Faibussowitsch   PetscCall(MatComputeBandwidth(pA, 0.0, &pbw));
809566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A,  NULL, "-orig_mat_view"));
819566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pA));
849566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pdm));
85c4762a1bSJed Brown   if (pbw > bw) {
8663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s increased bandwidth from %" PetscInt_FMT " to %" PetscInt_FMT "\n", order, bw, pbw));
87c4762a1bSJed Brown   } else {
8863a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s reduced bandwidth from %" PetscInt_FMT " to %" PetscInt_FMT "\n", order, bw, pbw));
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown PetscErrorCode CreateGroupLabel(DM dm, PetscInt numGroups, DMLabel *label, AppCtx *options)
94c4762a1bSJed Brown {
95c4762a1bSJed Brown   const PetscInt groupA[10] = {15, 3, 13, 12, 2, 10, 7, 6, 0, 4};
96c4762a1bSJed Brown   const PetscInt groupB[6]  = {14, 11, 9, 1, 8, 5};
97c4762a1bSJed Brown   PetscInt       c;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   PetscFunctionBegin;
100c4762a1bSJed Brown   if (numGroups < 2) {*label = NULL; PetscFunctionReturn(0);}
10163a3b9bcSJacob Faibussowitsch   PetscCheck(numGroups == 2,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Test only coded for 2 groups, not %" PetscInt_FMT, numGroups);
1029566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "groups", label));
1039566063dSJacob Faibussowitsch   for (c = 0; c < 10; ++c) PetscCall(DMLabelSetValue(*label, groupA[c], 101));
1049566063dSJacob Faibussowitsch   for (c = 0; c < 6;  ++c) PetscCall(DMLabelSetValue(*label, groupB[c], 1001));
105c4762a1bSJed Brown   PetscFunctionReturn(0);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown PetscErrorCode TestReorderingByGroup(DM dm, AppCtx *user)
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   DM              pdm;
111c4762a1bSJed Brown   DMLabel         label;
112c4762a1bSJed Brown   Mat             A, pA;
113c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
114c4762a1bSJed Brown   IS              perm;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(CreateGroupLabel(dm, user->numGroups, &label, user));
1189566063dSJacob Faibussowitsch   PetscCall(DMPlexGetOrdering(dm, order, label, &perm));
1199566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
1209566063dSJacob Faibussowitsch   PetscCall(DMPlexPermute(dm, perm, &pdm));
1219566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) pdm, "perm_"));
1229566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pdm));
1239566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm,  NULL, "-orig_dm_view"));
1249566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(pdm, NULL, "-perm_dm_view"));
1259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1269566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(dm, &A));
1279566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pdm, &pA));
1289566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(A,  NULL, "-orig_mat_view"));
1299566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(pA, NULL, "-perm_mat_view"));
1309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&pA));
1329566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&pdm));
133c4762a1bSJed Brown   PetscFunctionReturn(0);
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown int main(int argc, char **argv)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   DM             dm;
139c4762a1bSJed Brown   PetscSection   s;
140c4762a1bSJed Brown   AppCtx         user;
14130602db0SMatthew G. Knepley   PetscInt       dim;
142c4762a1bSJed Brown 
143*327415f7SBarry Smith   PetscFunctionBeginUser;
1449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1459566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(&user));
146c4762a1bSJed Brown   if (user.numGroups < 1) {
1479566063dSJacob Faibussowitsch     PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
1489566063dSJacob Faibussowitsch     PetscCall(DMSetType(dm, DMPLEX));
14930602db0SMatthew G. Knepley   } else {
1509566063dSJacob Faibussowitsch     PetscCall(CreateTestMesh(PETSC_COMM_WORLD, &dm, &user));
1510fdc7489SMatthew Knepley   }
1529566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1539566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1549566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
15530602db0SMatthew G. Knepley   {
15630602db0SMatthew G. Knepley     PetscInt  len = (dim+1) * PetscMax(1, user.numFields);
15730602db0SMatthew G. Knepley     PetscBool flg;
15830602db0SMatthew G. Knepley 
1599566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(len, &user.numDof));
160d0609cedSBarry Smith     PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");
1619566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex10.c", user.numDof, &len, &flg));
16263a3b9bcSJacob 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));
163d0609cedSBarry Smith     PetscOptionsEnd();
16430602db0SMatthew G. Knepley   }
16530602db0SMatthew G. Knepley   if (user.numGroups < 1) {
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(TestReordering(dm, &user));
172c4762a1bSJed Brown   } else {
1739566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(dm, user.numFields));
1749566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(dm));
1759566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s));
1769566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(dm, s));
1779566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&s));
1789566063dSJacob Faibussowitsch     PetscCall(TestReorderingByGroup(dm, &user));
179c4762a1bSJed Brown   }
1809566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1819566063dSJacob Faibussowitsch   PetscCall(CleanupContext(&user));
1829566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
183b122ec5aSJacob Faibussowitsch   return 0;
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /*TEST
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   # Two cell tests 0-3
189c4762a1bSJed Brown   test:
190c4762a1bSJed Brown     suffix: 0
1910fdc7489SMatthew Knepley     requires: triangle
19230602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -num_dof 1,0,0 -mat_view -dm_coord_space 0
193c4762a1bSJed Brown   test:
194c4762a1bSJed Brown     suffix: 1
19530602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -num_dof 1,0,0 -mat_view -dm_coord_space 0
196c4762a1bSJed Brown   test:
197c4762a1bSJed Brown     suffix: 2
1980fdc7489SMatthew Knepley     requires: ctetgen
19930602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
200c4762a1bSJed Brown   test:
201c4762a1bSJed Brown     suffix: 3
20230602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -num_dof 1,0,0,0 -mat_view -dm_coord_space 0
203c4762a1bSJed Brown   # Refined tests 4-7
204c4762a1bSJed Brown   test:
205c4762a1bSJed Brown     suffix: 4
206c4762a1bSJed Brown     requires: triangle
20730602db0SMatthew G. Knepley     args: -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0
208c4762a1bSJed Brown   test:
209c4762a1bSJed Brown     suffix: 5
21030602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0
211c4762a1bSJed Brown   test:
212c4762a1bSJed Brown     suffix: 6
213c4762a1bSJed Brown     requires: ctetgen
21430602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 1 -dm_refine_volume_limit_pre 0.00625 -num_dof 1,0,0,0
215c4762a1bSJed Brown   test:
216c4762a1bSJed Brown     suffix: 7
21730602db0SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_refine 1 -num_dof 1,0,0,0
218c4762a1bSJed Brown   # Parallel tests
219c4762a1bSJed Brown   # Grouping tests
220c4762a1bSJed Brown   test:
221c4762a1bSJed Brown     suffix: group_1
222c4762a1bSJed Brown     args: -num_groups 1 -num_dof 1,0,0 -is_view -orig_mat_view -perm_mat_view
223c4762a1bSJed Brown   test:
224c4762a1bSJed Brown     suffix: group_2
225c4762a1bSJed Brown     args: -num_groups 2 -num_dof 1,0,0 -is_view -perm_mat_view
226c4762a1bSJed Brown 
227c4762a1bSJed Brown TEST*/
228