xref: /petsc/src/dm/impls/plex/tests/ex16.c (revision a44be8dc6eeccc8754dbff0bc56030e29894a2e1)
1c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
59371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) {
6d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
79566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
89566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
99566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
109566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
11c4762a1bSJed Brown   PetscFunctionReturn(0);
12c4762a1bSJed Brown }
13c4762a1bSJed Brown 
14d3ef4daaSMatthew G. Knepley // Label half of the cells
159371c9d4SSatish Balay static PetscErrorCode CreateHalfCellsLabel(DM dm, PetscBool lower, DMLabel *label) {
16d3ef4daaSMatthew G. Knepley   PetscInt cStart, cEnd, cStartSub, cEndSub;
17d3ef4daaSMatthew G. Knepley 
18d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
19d3ef4daaSMatthew G. Knepley   PetscCall(DMCreateLabel(dm, "cells"));
20d3ef4daaSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "cells", label));
21d3ef4daaSMatthew G. Knepley   PetscCall(DMLabelClearStratum(*label, 1));
22d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
239371c9d4SSatish Balay   if (lower) {
249371c9d4SSatish Balay     cStartSub = cStart;
259371c9d4SSatish Balay     cEndSub   = cEnd / 2;
269371c9d4SSatish Balay   } else {
279371c9d4SSatish Balay     cStartSub = cEnd / 2;
289371c9d4SSatish Balay     cEndSub   = cEnd;
299371c9d4SSatish Balay   }
30d3ef4daaSMatthew G. Knepley   for (PetscInt c = cStartSub; c < cEndSub; ++c) PetscCall(DMLabelSetValue(*label, c, 1));
31d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexLabelComplete(dm, *label));
32d3ef4daaSMatthew G. Knepley   PetscFunctionReturn(0);
33d3ef4daaSMatthew G. Knepley }
34d3ef4daaSMatthew G. Knepley 
35d3ef4daaSMatthew G. Knepley // Label everything on the right half of the domain
36*a44be8dcSMatthew G. Knepley static PetscErrorCode CreateHalfDomainLabel(DM dm, PetscBool lower, PetscReal height, DMLabel *label) {
37d3ef4daaSMatthew G. Knepley   PetscReal centroid[3];
38d3ef4daaSMatthew G. Knepley   PetscInt  cStart, cEnd, cdim;
39d3ef4daaSMatthew G. Knepley 
40d3ef4daaSMatthew G. Knepley   PetscFunctionBeginUser;
41*a44be8dcSMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
42d3ef4daaSMatthew G. Knepley   PetscCall(DMCreateLabel(dm, "cells"));
43d3ef4daaSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "cells", label));
44d3ef4daaSMatthew G. Knepley   PetscCall(DMLabelClearStratum(*label, 1));
45d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
46d3ef4daaSMatthew G. Knepley   PetscCall(DMGetCoordinateDim(dm, &cdim));
47d3ef4daaSMatthew G. Knepley   for (PetscInt c = cStart; c < cEnd; ++c) {
48d3ef4daaSMatthew G. Knepley     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, NULL, centroid, NULL));
49*a44be8dcSMatthew G. Knepley     if (height > 0.0 && PetscAbsReal(centroid[1] - height) > PETSC_SMALL) continue;
509371c9d4SSatish Balay     if (lower) {
519371c9d4SSatish Balay       if (centroid[0] < 0.5) PetscCall(DMLabelSetValue(*label, c, 1));
529371c9d4SSatish Balay     } else {
539371c9d4SSatish Balay       if (centroid[0] > 0.5) PetscCall(DMLabelSetValue(*label, c, 1));
549371c9d4SSatish Balay     }
55d3ef4daaSMatthew G. Knepley   }
56d3ef4daaSMatthew G. Knepley   PetscCall(DMPlexLabelComplete(dm, *label));
57d3ef4daaSMatthew G. Knepley   PetscFunctionReturn(0);
58d3ef4daaSMatthew G. Knepley }
59d3ef4daaSMatthew G. Knepley 
60836a3779SMatthew G. Knepley // Create a line of faces at a given x value
61836a3779SMatthew G. Knepley static PetscErrorCode CreateLineLabel(DM dm, PetscReal x, DMLabel *label) {
62836a3779SMatthew G. Knepley   PetscReal centroid[3];
63836a3779SMatthew G. Knepley   PetscInt  fStart, fEnd;
64836a3779SMatthew G. Knepley 
65836a3779SMatthew G. Knepley   PetscFunctionBeginUser;
66836a3779SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
67836a3779SMatthew G. Knepley   PetscCall(DMCreateLabel(dm, "faces"));
68836a3779SMatthew G. Knepley   PetscCall(DMGetLabel(dm, "faces", label));
69836a3779SMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
70836a3779SMatthew G. Knepley   for (PetscInt f = fStart; f < fEnd; ++f) {
71836a3779SMatthew G. Knepley     PetscCall(DMPlexComputeCellGeometryFVM(dm, f, NULL, centroid, NULL));
72836a3779SMatthew G. Knepley     if (PetscAbsReal(centroid[0] - x) < PETSC_SMALL) PetscCall(DMLabelSetValue(*label, f, 1));
73836a3779SMatthew G. Knepley   }
74836a3779SMatthew G. Knepley   PetscCall(DMPlexLabelComplete(dm, *label));
75836a3779SMatthew G. Knepley   PetscFunctionReturn(0);
76836a3779SMatthew G. Knepley }
77836a3779SMatthew G. Knepley 
78*a44be8dcSMatthew G. Knepley static PetscErrorCode CreateVolumeSubmesh(DM dm, PetscBool domain, PetscBool lower, PetscReal height, DM *subdm) {
79c4762a1bSJed Brown   DMLabel label, map;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBegin;
82*a44be8dcSMatthew G. Knepley   if (domain) PetscCall(CreateHalfDomainLabel(dm, lower, height, &label));
83d3ef4daaSMatthew G. Knepley   else PetscCall(CreateHalfCellsLabel(dm, lower, &label));
849566063dSJacob Faibussowitsch   PetscCall(DMPlexFilter(dm, label, 1, subdm));
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*subdm, "Submesh"));
86d3ef4daaSMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*subdm, "sub_"));
879566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*subdm, NULL, "-dm_view"));
889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSubpointMap(*subdm, &map));
89d3ef4daaSMatthew G. Knepley   PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view"));
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93836a3779SMatthew G. Knepley static PetscErrorCode TestBoundaryField(DM dm) {
94836a3779SMatthew G. Knepley   DM       subdm;
95836a3779SMatthew G. Knepley   DMLabel  label, map;
96836a3779SMatthew G. Knepley   PetscFV  fvm;
97836a3779SMatthew G. Knepley   Vec      gv;
98836a3779SMatthew G. Knepley   PetscInt cdim;
99836a3779SMatthew G. Knepley 
100836a3779SMatthew G. Knepley   PetscFunctionBeginUser;
101836a3779SMatthew G. Knepley   PetscCall(CreateLineLabel(dm, 0.5, &label));
102836a3779SMatthew G. Knepley   PetscCall(DMPlexFilter(dm, label, 1, &subdm));
103836a3779SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)subdm, "Submesh"));
104836a3779SMatthew G. Knepley   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subdm, "sub_"));
105836a3779SMatthew G. Knepley   PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view"));
106836a3779SMatthew G. Knepley   PetscCall(DMPlexGetSubpointMap(subdm, &map));
107836a3779SMatthew G. Knepley   PetscCall(PetscObjectViewFromOptions((PetscObject)map, NULL, "-map_view"));
108836a3779SMatthew G. Knepley 
109836a3779SMatthew G. Knepley   PetscCall(PetscFVCreate(PetscObjectComm((PetscObject)subdm), &fvm));
110836a3779SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fvm, "testField"));
111836a3779SMatthew G. Knepley   PetscCall(PetscFVSetNumComponents(fvm, 1));
112836a3779SMatthew G. Knepley   PetscCall(DMGetCoordinateDim(subdm, &cdim));
113836a3779SMatthew G. Knepley   PetscCall(PetscFVSetSpatialDimension(fvm, cdim));
114836a3779SMatthew G. Knepley   PetscCall(PetscFVSetFromOptions(fvm));
115836a3779SMatthew G. Knepley 
116836a3779SMatthew G. Knepley   PetscCall(DMAddField(subdm, NULL, (PetscObject)fvm));
117836a3779SMatthew G. Knepley   PetscCall(PetscFVDestroy(&fvm));
118836a3779SMatthew G. Knepley   PetscCall(DMCreateDS(subdm));
119836a3779SMatthew G. Knepley 
120836a3779SMatthew G. Knepley   PetscCall(DMCreateGlobalVector(subdm, &gv));
121836a3779SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)gv, "potential"));
122836a3779SMatthew G. Knepley   PetscCall(VecSet(gv, 1.));
123836a3779SMatthew G. Knepley   PetscCall(VecViewFromOptions(gv, NULL, "-vec_view"));
124836a3779SMatthew G. Knepley   PetscCall(VecDestroy(&gv));
125836a3779SMatthew G. Knepley   PetscCall(DMDestroy(&subdm));
126836a3779SMatthew G. Knepley   PetscFunctionReturn(0);
127836a3779SMatthew G. Knepley }
128836a3779SMatthew G. Knepley 
1299371c9d4SSatish Balay int main(int argc, char **argv) {
130c4762a1bSJed Brown   DM        dm, subdm;
131*a44be8dcSMatthew G. Knepley   PetscReal height = -1.0;
132836a3779SMatthew G. Knepley   PetscBool volume = PETSC_TRUE, domain = PETSC_FALSE;
133c4762a1bSJed Brown 
134327415f7SBarry Smith   PetscFunctionBeginUser;
1359566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
136*a44be8dcSMatthew G. Knepley   PetscCall(PetscOptionsGetReal(NULL, NULL, "-height", &height, NULL));
137836a3779SMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-volume", &volume, NULL));
138d3ef4daaSMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-domain", &domain, NULL));
139836a3779SMatthew G. Knepley 
1409566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
141836a3779SMatthew G. Knepley   if (volume) {
142*a44be8dcSMatthew G. Knepley     PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_TRUE, height, &subdm));
1439566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(subdm));
1449566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
145*a44be8dcSMatthew G. Knepley     PetscCall(CreateVolumeSubmesh(dm, domain, PETSC_FALSE, height, &subdm));
1469566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(subdm));
1479566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&subdm));
148836a3779SMatthew G. Knepley   } else {
149836a3779SMatthew G. Knepley     PetscCall(TestBoundaryField(dm));
150836a3779SMatthew G. Knepley   }
1519566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
153b122ec5aSJacob Faibussowitsch   return 0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /*TEST
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   test:
159c4762a1bSJed Brown     suffix: 0
160c4762a1bSJed Brown     requires: triangle
161d3ef4daaSMatthew G. Knepley     args: -dm_coord_space 0 -sub_dm_plex_check_all \
162d3ef4daaSMatthew G. Knepley           -dm_view ascii::ascii_info_detail -sub_dm_view ascii::ascii_info_detail -map_view
163d3ef4daaSMatthew G. Knepley 
164d3ef4daaSMatthew G. Knepley   # These tests check that filtering is stable when boundary point ownership could change, so it needs 3 processes
165d3ef4daaSMatthew G. Knepley   testset:
166d3ef4daaSMatthew G. Knepley     nsize: 3
167d3ef4daaSMatthew G. Knepley     requires: parmetis
168d3ef4daaSMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 20,20 -petscpartitioner_type parmetis -dm_distribute_overlap 1 -sub_dm_distribute 0 \
169d3ef4daaSMatthew G. Knepley           -sub_dm_plex_check_all -dm_view -sub_dm_view
170d3ef4daaSMatthew G. Knepley 
171d3ef4daaSMatthew G. Knepley     test:
172d3ef4daaSMatthew G. Knepley       suffix: 1
173d3ef4daaSMatthew G. Knepley 
174d3ef4daaSMatthew G. Knepley     test:
175d3ef4daaSMatthew G. Knepley       suffix: 2
176d3ef4daaSMatthew G. Knepley       args: -domain
177c4762a1bSJed Brown 
178*a44be8dcSMatthew G. Knepley   # This set tests that global numberings can be made when some strata are missing on a process
179*a44be8dcSMatthew G. Knepley   testset:
180*a44be8dcSMatthew G. Knepley     nsize: 3
181*a44be8dcSMatthew G. Knepley     requires: hdf5
182*a44be8dcSMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_plex_box_faces 4,4 -petscpartitioner_type simple -sub_dm_distribute 0 \
183*a44be8dcSMatthew G. Knepley           -sub_dm_plex_check_all -sub_dm_view hdf5:subdm.h5
184*a44be8dcSMatthew G. Knepley 
185*a44be8dcSMatthew G. Knepley     test:
186*a44be8dcSMatthew G. Knepley       suffix: 3
187*a44be8dcSMatthew G. Knepley       args: -domain -height 0.625
188*a44be8dcSMatthew G. Knepley 
189836a3779SMatthew G. Knepley   # This test checks whether filter can extract a lower-dimensional manifold and output a field on it
190836a3779SMatthew G. Knepley   testset:
191836a3779SMatthew G. Knepley     args: -volume 0 -dm_plex_simplex 0 -sub_dm_view hdf5:subdm.h5 -vec_view hdf5:subdm.h5::append
192836a3779SMatthew G. Knepley     requires: hdf5
193836a3779SMatthew G. Knepley 
194836a3779SMatthew G. Knepley     test:
195836a3779SMatthew G. Knepley       suffix: surface_2d
196836a3779SMatthew G. Knepley       args: -dm_plex_box_faces 5,5
197836a3779SMatthew G. Knepley 
198836a3779SMatthew G. Knepley     test:
199836a3779SMatthew G. Knepley       suffix: surface_3d
200836a3779SMatthew G. Knepley       args: -dm_plex_box_faces 5,5,5
201836a3779SMatthew G. Knepley 
202c4762a1bSJed Brown TEST*/
203