xref: /petsc/src/dm/impls/moab/tests/ex3.cxx (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static char help[] = "Create a box mesh with DMMoab and test defining a tag on the mesh\n\n";
2 
3 #include <petscdmmoab.h>
4 
5 typedef struct {
6   DM            dm;                /* DM implementation using the MOAB interface */
7   PetscBool     debug;             /* The debugging level */
8   PetscLogEvent createMeshEvent;
9   /* Domain and mesh definition */
10   PetscInt      dim;                            /* The topological mesh dimension */
11   PetscInt      nele;                           /* Elements in each dimension */
12   PetscInt      degree;                         /* Degree of refinement */
13   PetscBool     simplex;                        /* Use simplex elements */
14   PetscInt      nlevels;                        /* Number of levels in mesh hierarchy */
15   PetscInt      nghost;                        /* Number of ghost layers in the mesh */
16   char          input_file[PETSC_MAX_PATH_LEN];   /* Import mesh from file */
17   char          output_file[PETSC_MAX_PATH_LEN];   /* Output mesh file name */
18   PetscBool     write_output;                        /* Write output mesh and data to file */
19 } AppCtx;
20 
21 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22 {
23   PetscFunctionBegin;
24   options->debug             = PETSC_FALSE;
25   options->nlevels           = 1;
26   options->nghost            = 1;
27   options->dim               = 2;
28   options->nele              = 5;
29   options->degree            = 2;
30   options->simplex           = PETSC_FALSE;
31   options->write_output      = PETSC_FALSE;
32   options->input_file[0]     = '\0';
33   PetscCall(PetscStrcpy(options->output_file,"ex3.h5m"));
34 
35   PetscOptionsBegin(comm, "", "Uniform Mesh Refinement Options", "DMMOAB");
36   PetscCall(PetscOptionsBool("-debug", "Enable debug messages", "ex2.cxx", options->debug, &options->debug, NULL));
37   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex3.cxx", options->dim, &options->dim, NULL,0,3));
38   PetscCall(PetscOptionsBoundedInt("-n", "The number of elements in each dimension", "ex3.cxx", options->nele, &options->nele, NULL,1));
39   PetscCall(PetscOptionsBoundedInt("-levels", "Number of levels in the hierarchy", "ex3.cxx", options->nlevels, &options->nlevels, NULL,0));
40   PetscCall(PetscOptionsBoundedInt("-degree", "Number of degrees at each level of refinement", "ex3.cxx", options->degree, &options->degree, NULL,0));
41   PetscCall(PetscOptionsBoundedInt("-ghost", "Number of ghost layers in the mesh", "ex3.cxx", options->nghost, &options->nghost, NULL,0));
42   PetscCall(PetscOptionsBool("-simplex", "Create simplices instead of tensor product elements", "ex3.cxx", options->simplex, &options->simplex, NULL));
43   PetscCall(PetscOptionsString("-input", "The input mesh file", "ex3.cxx", options->input_file, options->input_file, sizeof(options->input_file), NULL));
44   PetscCall(PetscOptionsString("-io", "Write out the mesh and solution that is defined on it (Default H5M format)", "ex3.cxx", options->output_file, options->output_file, sizeof(options->output_file), &options->write_output));
45   PetscOptionsEnd();
46 
47   PetscCall(PetscLogEventRegister("CreateMesh",          DM_CLASSID,   &options->createMeshEvent));
48   PetscFunctionReturn(0);
49 };
50 
51 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user)
52 {
53   size_t         len;
54   PetscMPIInt    rank;
55 
56   PetscFunctionBegin;
57   PetscCall(PetscLogEventBegin(user->createMeshEvent,0,0,0,0));
58   PetscCallMPI(MPI_Comm_rank(comm, &rank));
59   PetscCall(PetscStrlen(user->input_file, &len));
60   if (len) {
61     if (user->debug) PetscCall(PetscPrintf(comm, "Loading mesh from file: %s and creating the coarse level DM object.\n",user->input_file));
62     PetscCall(DMMoabLoadFromFile(comm, user->dim, user->nghost, user->input_file, "", &user->dm));
63   } else {
64     if (user->debug) PetscCall(PetscPrintf(comm, "Creating a %" PetscInt_FMT "-dimensional structured %s mesh of %" PetscInt_FMT "x%" PetscInt_FMT "x%" PetscInt_FMT " in memory and creating a DM object.\n",user->dim,(user->simplex?"simplex":"regular"),user->nele,user->nele,user->nele));
65     PetscCall(DMMoabCreateBoxMesh(comm, user->dim, user->simplex, NULL, user->nele, user->nghost, &user->dm));
66   }
67   PetscCall(PetscObjectSetName((PetscObject)user->dm, "Coarse Mesh"));
68   PetscCall(PetscLogEventEnd(user->createMeshEvent,0,0,0,0));
69   PetscFunctionReturn(0);
70 }
71 
72 int main(int argc, char **argv)
73 {
74   AppCtx         user;                 /* user-defined work context */
75   MPI_Comm       comm;
76   PetscInt       i;
77   Mat            R;
78   DM            *dmhierarchy;
79   PetscInt      *degrees;
80   PetscBool      createR;
81 
82   PetscFunctionBeginUser;
83   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
84   comm = PETSC_COMM_WORLD;
85   createR = PETSC_FALSE;
86 
87   PetscCall(ProcessOptions(comm, &user));
88   PetscCall(CreateMesh(comm, &user));
89   PetscCall(DMSetFromOptions(user.dm));
90 
91   /* SetUp the data structures for DMMOAB */
92   PetscCall(DMSetUp(user.dm));
93 
94   PetscCall(PetscMalloc(sizeof(DM)*(user.nlevels+1),&dmhierarchy));
95   for (i=0; i<=user.nlevels; i++) dmhierarchy[i] = NULL;
96 
97   // coarsest grid = 0
98   // finest grid = nlevels
99   dmhierarchy[0] = user.dm;
100   PetscObjectReference((PetscObject)user.dm);
101 
102   if (user.nlevels) {
103     PetscCall(PetscMalloc1(user.nlevels, &degrees));
104     for (i=0; i < user.nlevels; i++) degrees[i] = user.degree;
105     if (user.debug) PetscCall(PetscPrintf(comm, "Generate the MOAB mesh hierarchy with %" PetscInt_FMT " levels.\n", user.nlevels));
106     PetscCall(DMMoabGenerateHierarchy(user.dm,user.nlevels,degrees));
107 
108     PetscBool usehierarchy=PETSC_FALSE;
109     if (usehierarchy) PetscCall(DMRefineHierarchy(user.dm,user.nlevels,&dmhierarchy[1]));
110     else {
111       if (user.debug) {
112         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Level %" PetscInt_FMT "\n", 0));
113         PetscCall(DMView(user.dm, 0));
114       }
115       for (i=1; i<=user.nlevels; i++) {
116         if (user.debug) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Level %" PetscInt_FMT "\n", i));
117         PetscCall(DMRefine(dmhierarchy[i-1],MPI_COMM_NULL,&dmhierarchy[i]));
118         if (createR) PetscCall(DMCreateInterpolation(dmhierarchy[i-1],dmhierarchy[i],&R,NULL));
119         if (user.debug) {
120           PetscCall(DMView(dmhierarchy[i], 0));
121           if (createR) PetscCall(MatView(R,0));
122         }
123         /* Solvers could now set operator "R" to the multigrid PC object for level i
124             PCMGSetInterpolation(pc,i,R)
125         */
126         if (createR) PetscCall(MatDestroy(&R));
127       }
128     }
129   }
130 
131   if (user.write_output) {
132     if (user.debug) PetscCall(PetscPrintf(comm, "Output mesh hierarchy to file: %s.\n",user.output_file));
133     PetscCall(DMMoabOutput(dmhierarchy[user.nlevels],(const char*)user.output_file,""));
134   }
135 
136   for (i=0; i<=user.nlevels; i++) PetscCall(DMDestroy(&dmhierarchy[i]));
137   PetscCall(PetscFree(degrees));
138   PetscCall(PetscFree(dmhierarchy));
139   PetscCall(DMDestroy(&user.dm));
140   PetscCall(PetscFinalize());
141   return 0;
142 }
143 
144 /*TEST
145 
146      build:
147        requires: moab
148 
149      test:
150        args: -debug -n 2 -dim 2 -levels 2 -simplex
151        filter:  grep -v "DM_0x"
152 
153      test:
154        args: -debug -n 2 -dim 3 -levels 2
155        filter:  grep -v "DM_0x"
156        suffix: 1_2
157 
158      test:
159        args: -debug -n 2 -dim 3 -ghost 1 -levels 2
160        filter:  grep -v "DM_0x"
161        nsize: 2
162        suffix: 2_1
163 
164 TEST*/
165