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