Actual source code: ex3.cxx
petsc-3.10.5 2019-03-28
1: static char help[] = "Create a box mesh with DMMoab and test defining a tag on the mesh\n\n";
3: #include <petscdmmoab.h>
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;
21: PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
22: {
26: options->debug = PETSC_FALSE;
27: options->nlevels = 1;
28: options->nghost = 1;
29: options->dim = 2;
30: options->nele = 5;
31: options->degree = 2;
32: options->simplex = PETSC_FALSE;
33: options->write_output = PETSC_FALSE;
34: options->input_file[0] = '\0';
35: PetscStrcpy(options->output_file,"ex3.h5m");
37: PetscOptionsBegin(comm, "", "Uniform Mesh Refinement Options", "DMMOAB");
38: PetscOptionsBool("-debug", "Enable debug messages", "ex2.cxx", options->debug, &options->debug, NULL);
39: PetscOptionsInt("-dim", "The topological mesh dimension", "ex3.cxx", options->dim, &options->dim, NULL);
40: PetscOptionsInt("-n", "The number of elements in each dimension", "ex3.cxx", options->nele, &options->nele, NULL);
41: PetscOptionsInt("-levels", "Number of levels in the hierarchy", "ex3.cxx", options->nlevels, &options->nlevels, NULL);
42: PetscOptionsInt("-degree", "Number of degrees at each level of refinement", "ex3.cxx", options->degree, &options->degree, NULL);
43: PetscOptionsInt("-ghost", "Number of ghost layers in the mesh", "ex3.cxx", options->nghost, &options->nghost, NULL);
44: PetscOptionsBool("-simplex", "Create simplices instead of tensor product elements", "ex3.cxx", options->simplex, &options->simplex, NULL);
45: PetscOptionsString("-input", "The input mesh file", "ex3.cxx", options->input_file, options->input_file, PETSC_MAX_PATH_LEN, NULL);
46: PetscOptionsString("-io", "Write out the mesh and solution that is defined on it (Default H5M format)", "ex3.cxx", options->output_file, options->output_file, PETSC_MAX_PATH_LEN, &options->write_output);
47: PetscOptionsEnd();
49: PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent);
50: return(0);
51: };
53: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user)
54: {
55: size_t len;
56: PetscMPIInt rank;
60: PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
61: MPI_Comm_rank(comm, &rank);
62: PetscStrlen(user->input_file, &len);
63: if (len) {
64: if (user->debug) PetscPrintf(comm, "Loading mesh from file: %s and creating the coarse level DM object.\n",user->input_file);
65: DMMoabLoadFromFile(comm, user->dim, user->nghost, user->input_file, "", &user->dm);
66: }
67: else {
68: if (user->debug) {
69: PetscPrintf(comm, "Creating a %D-dimensional structured %s mesh of %Dx%Dx%D in memory and creating a DM object.\n",user->dim,(user->simplex?"simplex":"regular"),user->nele,user->nele,user->nele);
70: }
71: DMMoabCreateBoxMesh(comm, user->dim, user->simplex, NULL, user->nele, user->nghost, &user->dm);
72: }
74: PetscObjectSetName((PetscObject)user->dm, "Coarse Mesh");
75: PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
76: return(0);
77: }
79: int main(int argc, char **argv)
80: {
81: AppCtx user; /* user-defined work context */
82: MPI_Comm comm;
83: PetscInt i;
84: Mat R;
85: DM *dmhierarchy;
86: PetscInt *degrees;
87: PetscBool createR;
90: PetscInitialize(&argc, &argv, NULL, help);
92: /* Initialize input */
93: comm = PETSC_COMM_WORLD;
94: createR = PETSC_FALSE;
96: ProcessOptions(comm, &user);
97: CreateMesh(comm, &user);
98: DMSetFromOptions(user.dm);
100: /* SetUp the data structures for DMMOAB */
101: DMSetUp(user.dm);
103: PetscMalloc(sizeof(DM)*(user.nlevels+1),&dmhierarchy);
104: for (i=0; i<=user.nlevels; i++) dmhierarchy[i] = NULL;
106: // coarsest grid = 0
107: // finest grid = nlevels
108: dmhierarchy[0] = user.dm;
109: PetscObjectReference((PetscObject)user.dm);
111: if (user.nlevels) {
112: PetscMalloc1(user.nlevels, °rees);
113: for (i=0; i < user.nlevels; i++) degrees[i] = user.degree;
114: if (user.debug) PetscPrintf(comm, "Generate the MOAB mesh hierarchy with %D levels.\n", user.nlevels);
115: DMMoabGenerateHierarchy(user.dm,user.nlevels,degrees);
117: PetscBool usehierarchy=PETSC_FALSE;
118: if (usehierarchy) {
119: DMRefineHierarchy(user.dm,user.nlevels,&dmhierarchy[1]);
120: }
121: else {
122: if (user.debug) {
123: PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", 0);
124: DMView(user.dm, 0);
125: }
126: for (i=1; i<=user.nlevels; i++) {
127: if (user.debug) PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", i);
128: DMRefine(dmhierarchy[i-1],MPI_COMM_NULL,&dmhierarchy[i]);
129: if (createR) {
130: DMCreateInterpolation(dmhierarchy[i-1],dmhierarchy[i],&R,NULL);
131: }
132: if (user.debug) {
133: DMView(dmhierarchy[i], 0);
134: if (createR) {
135: MatView(R,0);
136: }
137: }
138: /* Solvers could now set operator "R" to the multigrid PC object for level i
139: PCMGSetInterpolation(pc,i,R)
140: */
141: if (createR) {
142: MatDestroy(&R);
143: }
144: }
145: }
146: }
148: if (user.write_output) {
149: if (user.debug) PetscPrintf(comm, "Output mesh hierarchy to file: %s.\n",user.output_file);
150: DMMoabOutput(dmhierarchy[user.nlevels],(const char*)user.output_file,"");
151: }
153: for (i=0; i<=user.nlevels; i++) {
154: DMDestroy(&dmhierarchy[i]);
155: }
156: PetscFree(degrees);
157: PetscFree(dmhierarchy);
158: DMDestroy(&user.dm);
159: PetscFinalize();
160: return 0;
161: }