xref: /petsc/src/dm/impls/plex/tests/ex10.c (revision a4a685f2b656679e5de30bfff53536422329f553)
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   PetscBool cellSimplex;       /* Flag for simplices */
8c4762a1bSJed Brown   PetscBool interpolate;       /* Flag for mesh interpolation */
9c4762a1bSJed Brown   PetscBool refinementUniform; /* Uniformly refine the mesh */
10c4762a1bSJed Brown   PetscReal refinementLimit;   /* Maximum volume of a refined cell */
11c4762a1bSJed Brown   PetscInt  numFields;         /* The number of section fields */
12c4762a1bSJed Brown   PetscInt *numComponents;     /* The number of field components */
13c4762a1bSJed Brown   PetscInt *numDof;            /* The dof signature for the section */
14c4762a1bSJed Brown   PetscInt  numGroups;         /* If greater than 1, use grouping in test */
15c4762a1bSJed Brown } AppCtx;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown PetscErrorCode ProcessOptions(AppCtx *options)
18c4762a1bSJed Brown {
19c4762a1bSJed Brown   PetscInt       len;
20c4762a1bSJed Brown   PetscBool      flg;
21c4762a1bSJed Brown   PetscErrorCode ierr;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   PetscFunctionBegin;
24c4762a1bSJed Brown   options->dim               = 2;
25c4762a1bSJed Brown   options->cellSimplex       = PETSC_TRUE;
26c4762a1bSJed Brown   options->interpolate       = PETSC_FALSE;
27c4762a1bSJed Brown   options->refinementUniform = PETSC_FALSE;
28c4762a1bSJed Brown   options->refinementLimit   = 0.0;
29c4762a1bSJed Brown   options->numFields         = 1;
30c4762a1bSJed Brown   options->numComponents     = NULL;
31c4762a1bSJed Brown   options->numDof            = NULL;
32c4762a1bSJed Brown   options->numGroups         = 0;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_SELF, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex10.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Flag for simplices", "ex10.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = PetscOptionsBool("-interpolate", "Flag for mesh interpolation", "ex10.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscOptionsBool("-refinement_uniform", "Uniformly refine the mesh", "ex10.c", options->refinementUniform, &options->refinementUniform, NULL);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = PetscOptionsReal("-refinement_limit", "The maximum volume of a refined cell", "ex10.c", options->refinementLimit, &options->refinementLimit, NULL);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_fields", "The number of section fields", "ex10.c", options->numFields, &options->numFields, NULL,1);CHKERRQ(ierr);
41c4762a1bSJed Brown   if (options->numFields) {
42c4762a1bSJed Brown     len  = options->numFields;
43c4762a1bSJed Brown     ierr = PetscCalloc1(len, &options->numComponents);CHKERRQ(ierr);
44c4762a1bSJed Brown     ierr = PetscOptionsIntArray("-num_components", "The number of components per field", "ex10.c", options->numComponents, &len, &flg);CHKERRQ(ierr);
45c4762a1bSJed Brown     if (flg && (len != options->numFields)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of components array is %D should be %D", len, options->numFields);
46c4762a1bSJed Brown   }
47c4762a1bSJed Brown   len  = (options->dim+1) * PetscMax(1, options->numFields);
48c4762a1bSJed Brown   ierr = PetscMalloc1(len, &options->numDof);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-num_dof", "The dof signature for the section", "ex10.c", options->numDof, &len, &flg);CHKERRQ(ierr);
50c4762a1bSJed Brown   if (flg && (len != (options->dim+1) * PetscMax(1, options->numFields))) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Length of dof array is %D should be %D", len, (options->dim+1) * PetscMax(1, options->numFields));
51c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_groups", "Group permutation by this many label values", "ex10.c", options->numGroups, &options->numGroups, NULL,0);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown PetscErrorCode CleanupContext(AppCtx *user)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   PetscErrorCode ierr;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   PetscFunctionBegin;
61c4762a1bSJed Brown   ierr = PetscFree(user->numComponents);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = PetscFree(user->numDof);CHKERRQ(ierr);
63c4762a1bSJed Brown   PetscFunctionReturn(0);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /* This mesh comes from~\cite{saad2003}, Fig. 2.10, p. 70. */
67c4762a1bSJed Brown PetscErrorCode CreateTestMesh(MPI_Comm comm, DM *dm, AppCtx *options)
68c4762a1bSJed Brown {
69*a4a685f2SJacob 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,
70c4762a1bSJed 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};
71*a4a685f2SJacob Faibussowitsch   const PetscReal   coords[15*2] = {0, -3,  0, -1,  2, -1,  0,  1,  2, 1,
72c4762a1bSJed Brown                                     0,  3,  1, -2,  1, -1,  0, -2,  2, 0,
73c4762a1bSJed Brown                                     1,  0,  1,  1,  0,  0,  1,  2,  0, 2};
74c4762a1bSJed Brown   PetscErrorCode ierr;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   PetscFunctionBegin;
77*a4a685f2SJacob Faibussowitsch   ierr = DMPlexCreateFromCellListPetsc(comm, 2, 16, 15, 3, options->interpolate, cells, 2, coords, dm);CHKERRQ(ierr);
78c4762a1bSJed Brown   PetscFunctionReturn(0);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown PetscErrorCode TestReordering(DM dm, AppCtx *user)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   DM              pdm;
84c4762a1bSJed Brown   IS              perm;
85c4762a1bSJed Brown   Mat             A, pA;
86c4762a1bSJed Brown   PetscInt        bw, pbw;
87c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
88c4762a1bSJed Brown   PetscErrorCode  ierr;
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   PetscFunctionBegin;
91c4762a1bSJed Brown   ierr = DMPlexGetOrdering(dm, order, NULL, &perm);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = DMPlexPermute(dm, perm, &pdm);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = DMSetFromOptions(pdm);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
95c4762a1bSJed Brown   ierr = DMCreateMatrix(dm, &A);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = DMCreateMatrix(pdm, &pA);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = MatComputeBandwidth(A, 0.0, &bw);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = MatComputeBandwidth(pA, 0.0, &pbw);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = MatDestroy(&pA);CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = DMDestroy(&pdm);CHKERRQ(ierr);
102c4762a1bSJed Brown   if (pbw > bw) {
103c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s increased bandwidth from %D to %D\n", order, bw, pbw);CHKERRQ(ierr);
104c4762a1bSJed Brown   } else {
105c4762a1bSJed Brown     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Ordering method %s reduced bandwidth from %D to %D\n", order, bw, pbw);CHKERRQ(ierr);
106c4762a1bSJed Brown   }
107c4762a1bSJed Brown   PetscFunctionReturn(0);
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown PetscErrorCode CreateGroupLabel(DM dm, PetscInt numGroups, DMLabel *label, AppCtx *options)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   const PetscInt groupA[10] = {15, 3, 13, 12, 2, 10, 7, 6, 0, 4};
113c4762a1bSJed Brown   const PetscInt groupB[6]  = {14, 11, 9, 1, 8, 5};
114c4762a1bSJed Brown   PetscInt       c;
115c4762a1bSJed Brown   PetscErrorCode ierr;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBegin;
118c4762a1bSJed Brown   if (numGroups < 2) {*label = NULL; PetscFunctionReturn(0);}
119c4762a1bSJed Brown   if (numGroups != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Test only coded for 2 groups, not %D", numGroups);
120c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, "groups", label);CHKERRQ(ierr);
121c4762a1bSJed Brown   for (c = 0; c < 10; ++c) {ierr = DMLabelSetValue(*label, groupA[c], 101);CHKERRQ(ierr);}
122c4762a1bSJed Brown   for (c = 0; c < 6;  ++c) {ierr = DMLabelSetValue(*label, groupB[c], 1001);CHKERRQ(ierr);}
123c4762a1bSJed Brown   PetscFunctionReturn(0);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown PetscErrorCode TestReorderingByGroup(DM dm, AppCtx *user)
127c4762a1bSJed Brown {
128c4762a1bSJed Brown   DM              pdm;
129c4762a1bSJed Brown   DMLabel         label;
130c4762a1bSJed Brown   Mat             A, pA;
131c4762a1bSJed Brown   MatOrderingType order = MATORDERINGRCM;
132c4762a1bSJed Brown   IS              perm;
133c4762a1bSJed Brown   PetscErrorCode  ierr;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   PetscFunctionBegin;
136c4762a1bSJed Brown   ierr = CreateGroupLabel(dm, user->numGroups, &label, user);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = DMPlexGetOrdering(dm, order, label, &perm);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
139c4762a1bSJed Brown   ierr = DMPlexPermute(dm, perm, &pdm);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = DMSetFromOptions(pdm);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = DMCreateMatrix(dm, &A);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = DMCreateMatrix(pdm, &pA);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = MatViewFromOptions(A,  NULL, "-orig_mat_view");CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = MatViewFromOptions(pA, NULL, "-perm_mat_view");CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = MatDestroy(&pA);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = DMDestroy(&pdm);CHKERRQ(ierr);
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown int main(int argc, char **argv)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   DM             dm;
155c4762a1bSJed Brown   PetscSection   s;
156c4762a1bSJed Brown   AppCtx         user;
157c4762a1bSJed Brown   PetscErrorCode ierr;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
160c4762a1bSJed Brown   ierr = ProcessOptions(&user);CHKERRQ(ierr);
161c4762a1bSJed Brown   if (user.numGroups < 1) {
162c4762a1bSJed Brown     ierr = DMPlexCreateDoublet(PETSC_COMM_WORLD, user.dim, user.cellSimplex, user.interpolate, user.refinementUniform, user.refinementLimit, &dm);CHKERRQ(ierr);
163c4762a1bSJed Brown     ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
164c4762a1bSJed Brown     ierr = DMSetNumFields(dm, user.numFields);CHKERRQ(ierr);
165c4762a1bSJed Brown     ierr = DMCreateDS(dm);CHKERRQ(ierr);
166c4762a1bSJed Brown     ierr = DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr);
167c4762a1bSJed Brown     ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr);
168c4762a1bSJed Brown     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
169c4762a1bSJed Brown     ierr = TestReordering(dm, &user);CHKERRQ(ierr);
170c4762a1bSJed Brown   } else {
171c4762a1bSJed Brown     ierr = CreateTestMesh(PETSC_COMM_WORLD, &dm, &user);CHKERRQ(ierr);
172c4762a1bSJed Brown     ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
173c4762a1bSJed Brown     ierr = DMSetNumFields(dm, user.numFields);CHKERRQ(ierr);
174c4762a1bSJed Brown     ierr = DMCreateDS(dm);CHKERRQ(ierr);
175c4762a1bSJed Brown     ierr = DMPlexCreateSection(dm, NULL, user.numComponents, user.numDof, 0, NULL, NULL, NULL, NULL, &s);CHKERRQ(ierr);
176c4762a1bSJed Brown     ierr = DMSetLocalSection(dm, s);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = PetscSectionDestroy(&s);CHKERRQ(ierr);
178c4762a1bSJed Brown     ierr = TestReorderingByGroup(dm, &user);CHKERRQ(ierr);
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
181c4762a1bSJed Brown   ierr = CleanupContext(&user);CHKERRQ(ierr);
182c4762a1bSJed Brown   ierr = PetscFinalize();
183c4762a1bSJed Brown   return ierr;
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /*TEST
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   # Two cell tests 0-3
189c4762a1bSJed Brown   test:
190c4762a1bSJed Brown     suffix: 0
191c4762a1bSJed Brown     args: -dim 2 -interpolate 1 -cell_simplex 1 -num_dof 1,0,0 -mat_view
192c4762a1bSJed Brown   test:
193c4762a1bSJed Brown     suffix: 1
194c4762a1bSJed Brown     args: -dim 2 -interpolate 1 -cell_simplex 0 -num_dof 1,0,0 -mat_view
195c4762a1bSJed Brown   test:
196c4762a1bSJed Brown     suffix: 2
197c4762a1bSJed Brown     args: -dim 3 -interpolate 1 -cell_simplex 1 -num_dof 1,0,0,0 -mat_view
198c4762a1bSJed Brown   test:
199c4762a1bSJed Brown     suffix: 3
200c4762a1bSJed Brown     args: -dim 3 -interpolate 1 -cell_simplex 0 -num_dof 1,0,0,0 -mat_view
201c4762a1bSJed Brown   # Refined tests 4-7
202c4762a1bSJed Brown   test:
203c4762a1bSJed Brown     suffix: 4
204c4762a1bSJed Brown     requires: triangle
205c4762a1bSJed Brown     args: -dim 2 -interpolate 1 -cell_simplex 1 -refinement_limit 0.00625 -num_dof 1,0,0
206c4762a1bSJed Brown   test:
207c4762a1bSJed Brown     suffix: 5
208c4762a1bSJed Brown     args: -dim 2 -interpolate 1 -cell_simplex 0 -refinement_uniform       -num_dof 1,0,0
209c4762a1bSJed Brown   test:
210c4762a1bSJed Brown     suffix: 6
211c4762a1bSJed Brown     requires: ctetgen
212c4762a1bSJed Brown     args: -dim 3 -interpolate 1 -cell_simplex 1 -refinement_limit 0.00625 -num_dof 1,0,0,0
213c4762a1bSJed Brown   test:
214c4762a1bSJed Brown     suffix: 7
215c4762a1bSJed Brown     args: -dim 3 -interpolate 1 -cell_simplex 0 -refinement_uniform       -num_dof 1,0,0,0
216c4762a1bSJed Brown   # Parallel tests
217c4762a1bSJed Brown   # Grouping tests
218c4762a1bSJed Brown   test:
219c4762a1bSJed Brown     suffix: group_1
220c4762a1bSJed Brown     args: -num_groups 1 -num_dof 1,0,0 -is_view -orig_mat_view -perm_mat_view
221c4762a1bSJed Brown   test:
222c4762a1bSJed Brown     suffix: group_2
223c4762a1bSJed Brown     args: -num_groups 2 -num_dof 1,0,0 -is_view -perm_mat_view
224c4762a1bSJed Brown 
225c4762a1bSJed Brown TEST*/
226