xref: /petsc/src/dm/impls/plex/tutorials/ex13.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a) !
1 static char help[] = "Create a Plex Schwarz P surface with quads\n\n";
2 
3 #include <petscdmplex.h>
4 
5 int main(int argc, char **argv) {
6   DM             dm;
7   PetscInt       extent[3] = {1, 1, 1}, refine = 0, layers = 0, three;
8   PetscReal      thickness   = 0.;
9   PetscBool      distribute  = PETSC_TRUE;
10   DMBoundaryType periodic[3] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE};
11   DMPlexTPSType  tps_type    = DMPLEX_TPS_SCHWARZ_P;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Schwarz P Example", NULL);
16   PetscCall(PetscOptionsIntArray("-extent", "Number of replicas for each of three dimensions", NULL, extent, (three = 3, &three), NULL));
17   PetscCall(PetscOptionsInt("-refine", "Number of refinements", NULL, refine, &refine, NULL));
18   PetscCall(PetscOptionsEnumArray("-periodic", "Periodicity in each of three dimensions", NULL, DMBoundaryTypes, (PetscEnum *)periodic, (three = 3, &three), NULL));
19   PetscCall(PetscOptionsBool("-distribute", "Distribute TPS manifold prior to refinement and extrusion", NULL, distribute, &distribute, NULL));
20   PetscCall(PetscOptionsInt("-layers", "Number of layers in volumetric extrusion (or zero to not extrude)", NULL, layers, &layers, NULL));
21   PetscCall(PetscOptionsReal("-thickness", "Thickness of volumetric extrusion", NULL, thickness, &thickness, NULL));
22   PetscCall(PetscOptionsEnum("-tps_type", "Type of triply-periodic surface", NULL, DMPlexTPSTypes, (PetscEnum)tps_type, (PetscEnum *)&tps_type, NULL));
23   PetscOptionsEnd();
24   PetscCall(DMPlexCreateTPSMesh(PETSC_COMM_WORLD, tps_type, extent, periodic, distribute, refine, layers, thickness, &dm));
25   PetscCall(PetscObjectSetName((PetscObject)dm, "TPS"));
26   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
27   PetscCall(DMDestroy(&dm));
28   PetscCall(PetscFinalize());
29   return 0;
30 }
31 
32 /*TEST
33 
34   test:
35     suffix: 0
36     args: -extent 1,2,3 -dm_view -refine 0
37   test:
38     suffix: 1
39     args: -extent 2,3,1 -dm_view -refine 1
40 
41   test:
42     suffix: gyroid_0
43     args: -extent 1,2,3 -dm_view -refine 0 -tps_type gyroid
44   test:
45     suffix: gyroid_1
46     args: -extent 2,3,1 -dm_view -refine 1 -tps_type gyroid
47   test:
48     suffix: extrude_0
49     args: -extent 2,3,1 -dm_view -refine 0 -layers 3 -thickness .2
50   test:
51     suffix: extrude_1_dist
52     nsize: 2
53     args: -extent 2,1,1 -dm_view -refine 1 -layers 3 -thickness .2
54 
55 TEST*/
56