Actual source code: ex3.cxx

  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: {

 25:   options->debug             = PETSC_FALSE;
 26:   options->nlevels           = 1;
 27:   options->nghost            = 1;
 28:   options->dim               = 2;
 29:   options->nele              = 5;
 30:   options->degree            = 2;
 31:   options->simplex           = PETSC_FALSE;
 32:   options->write_output      = PETSC_FALSE;
 33:   options->input_file[0]     = '\0';
 34:   PetscStrcpy(options->output_file,"ex3.h5m");

 36:   PetscOptionsBegin(comm, "", "Uniform Mesh Refinement Options", "DMMOAB");
 37:   PetscOptionsBool("-debug", "Enable debug messages", "ex2.cxx", options->debug, &options->debug, NULL);
 38:   PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex3.cxx", options->dim, &options->dim, NULL,0,3);
 39:   PetscOptionsBoundedInt("-n", "The number of elements in each dimension", "ex3.cxx", options->nele, &options->nele, NULL,1);
 40:   PetscOptionsBoundedInt("-levels", "Number of levels in the hierarchy", "ex3.cxx", options->nlevels, &options->nlevels, NULL,0);
 41:   PetscOptionsBoundedInt("-degree", "Number of degrees at each level of refinement", "ex3.cxx", options->degree, &options->degree, NULL,0);
 42:   PetscOptionsBoundedInt("-ghost", "Number of ghost layers in the mesh", "ex3.cxx", options->nghost, &options->nghost, NULL,0);
 43:   PetscOptionsBool("-simplex", "Create simplices instead of tensor product elements", "ex3.cxx", options->simplex, &options->simplex, NULL);
 44:   PetscOptionsString("-input", "The input mesh file", "ex3.cxx", options->input_file, options->input_file, sizeof(options->input_file), NULL);
 45:   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);
 46:   PetscOptionsEnd();

 48:   PetscLogEventRegister("CreateMesh",          DM_CLASSID,   &options->createMeshEvent);
 49:   return 0;
 50: };

 52: PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user)
 53: {
 54:   size_t         len;
 55:   PetscMPIInt    rank;

 57:   PetscLogEventBegin(user->createMeshEvent,0,0,0,0);
 58:   MPI_Comm_rank(comm, &rank);
 59:   PetscStrlen(user->input_file, &len);
 60:   if (len) {
 61:     if (user->debug) PetscPrintf(comm, "Loading mesh from file: %s and creating the coarse level DM object.\n",user->input_file);
 62:     DMMoabLoadFromFile(comm, user->dim, user->nghost, user->input_file, "", &user->dm);
 63:   } else {
 64:     if (user->debug) 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);
 65:     DMMoabCreateBoxMesh(comm, user->dim, user->simplex, NULL, user->nele, user->nghost, &user->dm);
 66:   }
 67:   PetscObjectSetName((PetscObject)user->dm, "Coarse Mesh");
 68:   PetscLogEventEnd(user->createMeshEvent,0,0,0,0);
 69:   return 0;
 70: }

 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;

 82:   PetscInitialize(&argc, &argv, NULL, help);
 83:   comm = PETSC_COMM_WORLD;
 84:   createR = PETSC_FALSE;

 86:   ProcessOptions(comm, &user);
 87:   CreateMesh(comm, &user);
 88:   DMSetFromOptions(user.dm);

 90:   /* SetUp the data structures for DMMOAB */
 91:   DMSetUp(user.dm);

 93:   PetscMalloc(sizeof(DM)*(user.nlevels+1),&dmhierarchy);
 94:   for (i=0; i<=user.nlevels; i++) dmhierarchy[i] = NULL;

 96:   // coarsest grid = 0
 97:   // finest grid = nlevels
 98:   dmhierarchy[0] = user.dm;
 99:   PetscObjectReference((PetscObject)user.dm);

101:   if (user.nlevels) {
102:     PetscMalloc1(user.nlevels, &degrees);
103:     for (i=0; i < user.nlevels; i++) degrees[i] = user.degree;
104:     if (user.debug) PetscPrintf(comm, "Generate the MOAB mesh hierarchy with %D levels.\n", user.nlevels);
105:     DMMoabGenerateHierarchy(user.dm,user.nlevels,degrees);

107:     PetscBool usehierarchy=PETSC_FALSE;
108:     if (usehierarchy) DMRefineHierarchy(user.dm,user.nlevels,&dmhierarchy[1]);
109:     else {
110:       if (user.debug) {
111:         PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", 0);
112:         DMView(user.dm, 0);
113:       }
114:       for (i=1; i<=user.nlevels; i++) {
115:         if (user.debug) PetscPrintf(PETSC_COMM_WORLD, "Level %D\n", i);
116:         DMRefine(dmhierarchy[i-1],MPI_COMM_NULL,&dmhierarchy[i]);
117:         if (createR) DMCreateInterpolation(dmhierarchy[i-1],dmhierarchy[i],&R,NULL);
118:         if (user.debug) {
119:           DMView(dmhierarchy[i], 0);
120:           if (createR) MatView(R,0);
121:         }
122:         /* Solvers could now set operator "R" to the multigrid PC object for level i
123:             PCMGSetInterpolation(pc,i,R)
124:         */
125:         if (createR) MatDestroy(&R);
126:       }
127:     }
128:   }

130:   if (user.write_output) {
131:     if (user.debug) PetscPrintf(comm, "Output mesh hierarchy to file: %s.\n",user.output_file);
132:     DMMoabOutput(dmhierarchy[user.nlevels],(const char*)user.output_file,"");
133:   }

135:   for (i=0; i<=user.nlevels; i++) DMDestroy(&dmhierarchy[i]);
136:   PetscFree(degrees);
137:   PetscFree(dmhierarchy);
138:   DMDestroy(&user.dm);
139:   PetscFinalize();
140:   return 0;
141: }

143: /*TEST

145:      build:
146:        requires: moab

148:      test:
149:        args: -debug -n 2 -dim 2 -levels 2 -simplex
150:        filter:  grep -v "DM_0x"

152:      test:
153:        args: -debug -n 2 -dim 3 -levels 2
154:        filter:  grep -v "DM_0x"
155:        suffix: 1_2

157:      test:
158:        args: -debug -n 2 -dim 3 -ghost 1 -levels 2
159:        filter:  grep -v "DM_0x"
160:        nsize: 2
161:        suffix: 2_1

163: TEST*/