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 PetscErrorCode ierr; 24 25 PetscFunctionBegin; 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 PetscCall(PetscStrcpy(options->output_file,"ex3.h5m")); 36 37 ierr = PetscOptionsBegin(comm, "", "Uniform Mesh Refinement Options", "DMMOAB");PetscCall(ierr); 38 PetscCall(PetscOptionsBool("-debug", "Enable debug messages", "ex2.cxx", options->debug, &options->debug, NULL)); 39 PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex3.cxx", options->dim, &options->dim, NULL,0,3)); 40 PetscCall(PetscOptionsBoundedInt("-n", "The number of elements in each dimension", "ex3.cxx", options->nele, &options->nele, NULL,1)); 41 PetscCall(PetscOptionsBoundedInt("-levels", "Number of levels in the hierarchy", "ex3.cxx", options->nlevels, &options->nlevels, NULL,0)); 42 PetscCall(PetscOptionsBoundedInt("-degree", "Number of degrees at each level of refinement", "ex3.cxx", options->degree, &options->degree, NULL,0)); 43 PetscCall(PetscOptionsBoundedInt("-ghost", "Number of ghost layers in the mesh", "ex3.cxx", options->nghost, &options->nghost, NULL,0)); 44 PetscCall(PetscOptionsBool("-simplex", "Create simplices instead of tensor product elements", "ex3.cxx", options->simplex, &options->simplex, NULL)); 45 PetscCall(PetscOptionsString("-input", "The input mesh file", "ex3.cxx", options->input_file, options->input_file, sizeof(options->input_file), NULL)); 46 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)); 47 ierr = PetscOptionsEnd();PetscCall(ierr); 48 49 PetscCall(PetscLogEventRegister("CreateMesh", DM_CLASSID, &options->createMeshEvent)); 50 PetscFunctionReturn(0); 51 }; 52 53 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user) 54 { 55 size_t len; 56 PetscMPIInt rank; 57 58 PetscFunctionBegin; 59 PetscCall(PetscLogEventBegin(user->createMeshEvent,0,0,0,0)); 60 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 61 PetscCall(PetscStrlen(user->input_file, &len)); 62 if (len) { 63 if (user->debug) PetscCall(PetscPrintf(comm, "Loading mesh from file: %s and creating the coarse level DM object.\n",user->input_file)); 64 PetscCall(DMMoabLoadFromFile(comm, user->dim, user->nghost, user->input_file, "", &user->dm)); 65 } else { 66 if (user->debug) PetscCall(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)); 67 PetscCall(DMMoabCreateBoxMesh(comm, user->dim, user->simplex, NULL, user->nele, user->nghost, &user->dm)); 68 } 69 PetscCall(PetscObjectSetName((PetscObject)user->dm, "Coarse Mesh")); 70 PetscCall(PetscLogEventEnd(user->createMeshEvent,0,0,0,0)); 71 PetscFunctionReturn(0); 72 } 73 74 int main(int argc, char **argv) 75 { 76 AppCtx user; /* user-defined work context */ 77 MPI_Comm comm; 78 PetscInt i; 79 Mat R; 80 DM *dmhierarchy; 81 PetscInt *degrees; 82 PetscBool createR; 83 84 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 85 comm = PETSC_COMM_WORLD; 86 createR = PETSC_FALSE; 87 88 PetscCall(ProcessOptions(comm, &user)); 89 PetscCall(CreateMesh(comm, &user)); 90 PetscCall(DMSetFromOptions(user.dm)); 91 92 /* SetUp the data structures for DMMOAB */ 93 PetscCall(DMSetUp(user.dm)); 94 95 PetscCall(PetscMalloc(sizeof(DM)*(user.nlevels+1),&dmhierarchy)); 96 for (i=0; i<=user.nlevels; i++) dmhierarchy[i] = NULL; 97 98 // coarsest grid = 0 99 // finest grid = nlevels 100 dmhierarchy[0] = user.dm; 101 PetscObjectReference((PetscObject)user.dm); 102 103 if (user.nlevels) { 104 PetscCall(PetscMalloc1(user.nlevels, °rees)); 105 for (i=0; i < user.nlevels; i++) degrees[i] = user.degree; 106 if (user.debug) PetscCall(PetscPrintf(comm, "Generate the MOAB mesh hierarchy with %D levels.\n", user.nlevels)); 107 PetscCall(DMMoabGenerateHierarchy(user.dm,user.nlevels,degrees)); 108 109 PetscBool usehierarchy=PETSC_FALSE; 110 if (usehierarchy) PetscCall(DMRefineHierarchy(user.dm,user.nlevels,&dmhierarchy[1])); 111 else { 112 if (user.debug) { 113 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", 0)); 114 PetscCall(DMView(user.dm, 0)); 115 } 116 for (i=1; i<=user.nlevels; i++) { 117 if (user.debug) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", i)); 118 PetscCall(DMRefine(dmhierarchy[i-1],MPI_COMM_NULL,&dmhierarchy[i])); 119 if (createR) PetscCall(DMCreateInterpolation(dmhierarchy[i-1],dmhierarchy[i],&R,NULL)); 120 if (user.debug) { 121 PetscCall(DMView(dmhierarchy[i], 0)); 122 if (createR) PetscCall(MatView(R,0)); 123 } 124 /* Solvers could now set operator "R" to the multigrid PC object for level i 125 PCMGSetInterpolation(pc,i,R) 126 */ 127 if (createR) PetscCall(MatDestroy(&R)); 128 } 129 } 130 } 131 132 if (user.write_output) { 133 if (user.debug) PetscCall(PetscPrintf(comm, "Output mesh hierarchy to file: %s.\n",user.output_file)); 134 PetscCall(DMMoabOutput(dmhierarchy[user.nlevels],(const char*)user.output_file,"")); 135 } 136 137 for (i=0; i<=user.nlevels; i++) PetscCall(DMDestroy(&dmhierarchy[i])); 138 PetscCall(PetscFree(degrees)); 139 PetscCall(PetscFree(dmhierarchy)); 140 PetscCall(DMDestroy(&user.dm)); 141 PetscCall(PetscFinalize()); 142 return 0; 143 } 144 145 /*TEST 146 147 build: 148 requires: moab 149 150 test: 151 args: -debug -n 2 -dim 2 -levels 2 -simplex 152 filter: grep -v "DM_0x" 153 154 test: 155 args: -debug -n 2 -dim 3 -levels 2 156 filter: grep -v "DM_0x" 157 suffix: 1_2 158 159 test: 160 args: -debug -n 2 -dim 3 -ghost 1 -levels 2 161 filter: grep -v "DM_0x" 162 nsize: 2 163 suffix: 2_1 164 165 TEST*/ 166