1402df9f0SLisandro Dalcin static char help[] = "Tests DMPlex Gmsh reader.\n\n"; 2402df9f0SLisandro Dalcin 3402df9f0SLisandro Dalcin #include <petscdmplex.h> 4402df9f0SLisandro Dalcin 5402df9f0SLisandro Dalcin #if !defined(PETSC_GMSH_EXE) 6402df9f0SLisandro Dalcin #define PETSC_GMSH_EXE "gmsh" 7402df9f0SLisandro Dalcin #endif 8402df9f0SLisandro Dalcin 9*77215723SLisandro Dalcin #include <petscds.h> 10*77215723SLisandro Dalcin 11*77215723SLisandro Dalcin static void one(PetscInt dim, PetscInt Nf, PetscInt NfAux, 12*77215723SLisandro Dalcin const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 13*77215723SLisandro Dalcin const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 14*77215723SLisandro Dalcin PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar value[]) 15*77215723SLisandro Dalcin { 16*77215723SLisandro Dalcin value[0] = (PetscReal)1; 17*77215723SLisandro Dalcin } 18*77215723SLisandro Dalcin 19*77215723SLisandro Dalcin static PetscErrorCode CreateFE(DM dm) 20*77215723SLisandro Dalcin { 21*77215723SLisandro Dalcin DM cdm; 22*77215723SLisandro Dalcin PetscSpace P; 23*77215723SLisandro Dalcin PetscDualSpace Q; 24*77215723SLisandro Dalcin DM K; 25*77215723SLisandro Dalcin PetscFE fe; 26*77215723SLisandro Dalcin DMPolytopeType ptype; 27*77215723SLisandro Dalcin 28*77215723SLisandro Dalcin PetscInt dim,k; 29*77215723SLisandro Dalcin PetscBool isSimplex; 30*77215723SLisandro Dalcin 31*77215723SLisandro Dalcin PetscDS ds; 32*77215723SLisandro Dalcin PetscErrorCode ierr; 33*77215723SLisandro Dalcin 34*77215723SLisandro Dalcin PetscFunctionBeginUser; 35*77215723SLisandro Dalcin ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 36*77215723SLisandro Dalcin ierr = DMGetField(cdm, 0, NULL, (PetscObject*) &fe); 37*77215723SLisandro Dalcin ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr); 38*77215723SLisandro Dalcin ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 39*77215723SLisandro Dalcin ierr = PetscDualSpaceGetDM(Q,&K);CHKERRQ(ierr); 40*77215723SLisandro Dalcin ierr = DMGetDimension(K,&dim);CHKERRQ(ierr); 41*77215723SLisandro Dalcin ierr = PetscSpaceGetDegree(P, &k, NULL);CHKERRQ(ierr); 42*77215723SLisandro Dalcin ierr = DMPlexGetCellType(K, 0, &ptype);CHKERRQ(ierr); 43*77215723SLisandro Dalcin switch (ptype) { 44*77215723SLisandro Dalcin case DM_POLYTOPE_QUADRILATERAL: 45*77215723SLisandro Dalcin case DM_POLYTOPE_HEXAHEDRON: 46*77215723SLisandro Dalcin isSimplex = PETSC_FALSE; break; 47*77215723SLisandro Dalcin default: 48*77215723SLisandro Dalcin isSimplex = PETSC_TRUE; break; 49*77215723SLisandro Dalcin } 50*77215723SLisandro Dalcin 51*77215723SLisandro Dalcin ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe);CHKERRQ(ierr); 52*77215723SLisandro Dalcin ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr); 53*77215723SLisandro Dalcin ierr = DMAddField(dm, NULL, (PetscObject) fe);CHKERRQ(ierr); 54*77215723SLisandro Dalcin ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 55*77215723SLisandro Dalcin ierr = DMCreateDS(dm);CHKERRQ(ierr); 56*77215723SLisandro Dalcin 57*77215723SLisandro Dalcin ierr = DMGetDS(dm, &ds);CHKERRQ(ierr); 58*77215723SLisandro Dalcin ierr = PetscDSSetObjective(ds, 0, one);CHKERRQ(ierr); 59*77215723SLisandro Dalcin PetscFunctionReturn(0); 60*77215723SLisandro Dalcin } 61*77215723SLisandro Dalcin 62*77215723SLisandro Dalcin static PetscErrorCode CheckIntegral(DM dm, PetscReal integral, PetscReal tol) 63*77215723SLisandro Dalcin { 64*77215723SLisandro Dalcin Vec u; 65*77215723SLisandro Dalcin PetscReal rval; 66*77215723SLisandro Dalcin PetscScalar result; 67*77215723SLisandro Dalcin PetscErrorCode ierr; 68*77215723SLisandro Dalcin 69*77215723SLisandro Dalcin PetscFunctionBeginUser; 70*77215723SLisandro Dalcin ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 71*77215723SLisandro Dalcin ierr = DMPlexComputeIntegralFEM(dm, u, &result, NULL);CHKERRQ(ierr); 72*77215723SLisandro Dalcin ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 73*77215723SLisandro Dalcin rval = PetscRealPart(result); 74*77215723SLisandro Dalcin if (integral > 0 && PetscAbsReal(integral - rval) > tol) { 75*77215723SLisandro Dalcin ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Calculated value %g != %g actual value (error %g > %g tol)\n", 76*77215723SLisandro Dalcin (double) rval, (double) integral, (double) PetscAbsReal(integral - rval), (double) tol);CHKERRQ(ierr); 77*77215723SLisandro Dalcin } 78*77215723SLisandro Dalcin PetscFunctionReturn(0); 79*77215723SLisandro Dalcin } 80*77215723SLisandro Dalcin 81402df9f0SLisandro Dalcin int main(int argc, char **argv) 82402df9f0SLisandro Dalcin { 83402df9f0SLisandro Dalcin DM dm; 84*77215723SLisandro Dalcin const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex", 85*77215723SLisandro Dalcin "B2tri", "B2qua", "B3tet", "B3hex"}; 86402df9f0SLisandro Dalcin const char *const fmtlist[] = {"msh22", "msh40", "msh41"}; 87402df9f0SLisandro Dalcin PetscInt msh = 5; 88402df9f0SLisandro Dalcin PetscInt fmt = 2; 89*77215723SLisandro Dalcin PetscBool bin = PETSC_TRUE; 90*77215723SLisandro Dalcin PetscInt dim = 3; 91402df9f0SLisandro Dalcin PetscInt order = 2; 92402df9f0SLisandro Dalcin 93*77215723SLisandro Dalcin const char cmdtemplate[] = "%s -format %s %s -%d -order %d %s -o %s"; 94402df9f0SLisandro Dalcin char gmsh[PETSC_MAX_PATH_LEN] = PETSC_GMSH_EXE; 95402df9f0SLisandro Dalcin char tag[PETSC_MAX_PATH_LEN], path[PETSC_MAX_PATH_LEN]; 96402df9f0SLisandro Dalcin char geo[PETSC_MAX_PATH_LEN], geodir[PETSC_MAX_PATH_LEN] = "."; 97402df9f0SLisandro Dalcin char out[PETSC_MAX_PATH_LEN], outdir[PETSC_MAX_PATH_LEN] = "."; 98402df9f0SLisandro Dalcin char cmd[PETSC_MAX_PATH_LEN*4]; 99402df9f0SLisandro Dalcin PetscBool set,flg; 100402df9f0SLisandro Dalcin FILE *fp; 101402df9f0SLisandro Dalcin PetscErrorCode ierr; 102402df9f0SLisandro Dalcin 103402df9f0SLisandro Dalcin ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 104402df9f0SLisandro Dalcin 105402df9f0SLisandro Dalcin ierr = PetscStrncpy(geodir, "${PETSC_DIR}/share/petsc/datafiles/meshes", sizeof(geodir));CHKERRQ(ierr); 106402df9f0SLisandro Dalcin ierr = PetscOptionsGetenv(PETSC_COMM_SELF, "GMSH", path, sizeof(path), &set);CHKERRQ(ierr); 107402df9f0SLisandro Dalcin if (set) {ierr = PetscStrncpy(gmsh, path, sizeof(gmsh));CHKERRQ(ierr);} 108402df9f0SLisandro Dalcin ierr = PetscOptionsGetString(NULL, NULL, "-gmsh", gmsh, sizeof(gmsh), NULL);CHKERRQ(ierr); 109402df9f0SLisandro Dalcin ierr = PetscOptionsGetString(NULL, NULL, "-dir", geodir, sizeof(geodir), NULL);CHKERRQ(ierr); 110402df9f0SLisandro Dalcin ierr = PetscOptionsGetString(NULL, NULL, "-out", outdir, sizeof(outdir), NULL);CHKERRQ(ierr); 111*77215723SLisandro Dalcin ierr = PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, (int)(sizeof(mshlist)/sizeof(mshlist[0])), &msh, NULL);CHKERRQ(ierr); 112*77215723SLisandro Dalcin ierr = PetscOptionsGetEList(NULL, NULL, "-fmt", fmtlist, (int)(sizeof(fmtlist)/sizeof(fmtlist[0])), &fmt, NULL);CHKERRQ(ierr); 113402df9f0SLisandro Dalcin ierr = PetscOptionsGetBool(NULL, NULL, "-bin", &bin, NULL);CHKERRQ(ierr); 114*77215723SLisandro Dalcin ierr = PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL);CHKERRQ(ierr); 115402df9f0SLisandro Dalcin ierr = PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL);CHKERRQ(ierr); 116402df9f0SLisandro Dalcin if (fmt == 1) bin = PETSC_FALSE; /* Recent Gmsh releases cannot generate msh40+binary format*/ 117402df9f0SLisandro Dalcin 118402df9f0SLisandro Dalcin { /* This test requires Gmsh >= 4.2.0 */ 119402df9f0SLisandro Dalcin int inum = 0, major = 0, minor = 0, micro = 0; 120402df9f0SLisandro Dalcin ierr = PetscSNPrintf(cmd, sizeof(cmd), "%s -info", gmsh);CHKERRQ(ierr); 121402df9f0SLisandro Dalcin ierr = PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp);CHKERRQ(ierr); 122402df9f0SLisandro Dalcin if (fp) {inum = fscanf(fp, "Version : %d.%d.%d", &major, &minor, µ);} 123402df9f0SLisandro Dalcin ierr = PetscPClose(PETSC_COMM_SELF, fp);CHKERRQ(ierr); 124402df9f0SLisandro Dalcin if (inum != 3 || major < 4 || (major == 4 && minor < 2)) { 125402df9f0SLisandro Dalcin ierr = PetscPrintf(PETSC_COMM_SELF, "Gmsh>=4.2.0 not available\n");CHKERRQ(ierr); goto finish; 126402df9f0SLisandro Dalcin } 127402df9f0SLisandro Dalcin } 128402df9f0SLisandro Dalcin 129*77215723SLisandro Dalcin ierr = PetscSNPrintf(tag, sizeof(tag), "%s-%d-%d-%s%s", mshlist[msh], (int)dim, (int)order, fmtlist[fmt], bin?"-bin":"");CHKERRQ(ierr); 130402df9f0SLisandro Dalcin ierr = PetscSNPrintf(geo, sizeof(geo), "%s/gmsh-%s.geo", geodir, mshlist[msh]);CHKERRQ(ierr); 131402df9f0SLisandro Dalcin ierr = PetscSNPrintf(out, sizeof(out), "%s/mesh-%s.msh", outdir, tag);CHKERRQ(ierr); 132402df9f0SLisandro Dalcin ierr = PetscStrreplace(PETSC_COMM_SELF, geo, path, sizeof(path));CHKERRQ(ierr); 133402df9f0SLisandro Dalcin ierr = PetscFixFilename(path, geo);CHKERRQ(ierr); 134402df9f0SLisandro Dalcin ierr = PetscStrreplace(PETSC_COMM_SELF, out, path, sizeof(path));CHKERRQ(ierr); 135402df9f0SLisandro Dalcin ierr = PetscFixFilename(path, out);CHKERRQ(ierr); 136402df9f0SLisandro Dalcin ierr = PetscTestFile(geo, 'r', &flg);CHKERRQ(ierr); 137402df9f0SLisandro Dalcin if (!flg) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "File not found: %s", geo); 138402df9f0SLisandro Dalcin 139*77215723SLisandro Dalcin ierr = PetscSNPrintf(cmd, sizeof(cmd), cmdtemplate, gmsh, fmtlist[fmt], bin?"-bin":"", (int)dim, (int)order, geo, out);CHKERRQ(ierr); 140402df9f0SLisandro Dalcin ierr = PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp);CHKERRQ(ierr); 141402df9f0SLisandro Dalcin ierr = PetscPClose(PETSC_COMM_SELF, fp);CHKERRQ(ierr); 142402df9f0SLisandro Dalcin 143402df9f0SLisandro Dalcin ierr = DMPlexCreateFromFile(PETSC_COMM_SELF, out, PETSC_TRUE, &dm);CHKERRQ(ierr); 144402df9f0SLisandro Dalcin ierr = PetscSNPrintf(tag, sizeof(tag), "mesh-%s", mshlist[msh]);CHKERRQ(ierr); 145402df9f0SLisandro Dalcin ierr = PetscObjectSetName((PetscObject)dm, tag);CHKERRQ(ierr); 146402df9f0SLisandro Dalcin ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 147402df9f0SLisandro Dalcin ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 148*77215723SLisandro Dalcin { 149*77215723SLisandro Dalcin PetscBool check; 150*77215723SLisandro Dalcin PetscReal integral = 0, tol = (PetscReal)1.0e-4; 151*77215723SLisandro Dalcin ierr = PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check);CHKERRQ(ierr); 152*77215723SLisandro Dalcin ierr = PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL);CHKERRQ(ierr); 153*77215723SLisandro Dalcin if (check) { 154*77215723SLisandro Dalcin ierr = CreateFE(dm);CHKERRQ(ierr); 155*77215723SLisandro Dalcin ierr = CheckIntegral(dm, integral, tol);CHKERRQ(ierr); 156*77215723SLisandro Dalcin } 157*77215723SLisandro Dalcin } 158402df9f0SLisandro Dalcin ierr = DMDestroy(&dm);CHKERRQ(ierr); 159402df9f0SLisandro Dalcin 160402df9f0SLisandro Dalcin finish: 161402df9f0SLisandro Dalcin ierr = PetscFinalize(); 162402df9f0SLisandro Dalcin return ierr; 163402df9f0SLisandro Dalcin } 164402df9f0SLisandro Dalcin 165402df9f0SLisandro Dalcin /*TEST 166402df9f0SLisandro Dalcin 167402df9f0SLisandro Dalcin build: 168402df9f0SLisandro Dalcin requires: define(PETSC_HAVE_POPEN) 169402df9f0SLisandro Dalcin 170402df9f0SLisandro Dalcin test: 171402df9f0SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 172402df9f0SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 173402df9f0SLisandro Dalcin args: -msh {{seg tri qua tet wed hex}separate_output} 174402df9f0SLisandro Dalcin args: -order {{1 2 3 7}} 175402df9f0SLisandro Dalcin args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 176402df9f0SLisandro Dalcin args: -dm_view ::ascii_info_detail 177402df9f0SLisandro Dalcin args: -dm_plex_check_all 178066ea43fSLisandro Dalcin args: -dm_plex_gmsh_highorder false 179402df9f0SLisandro Dalcin 180*77215723SLisandro Dalcin 181*77215723SLisandro Dalcin testset: 182*77215723SLisandro Dalcin suffix: B2 # 2D ball 183*77215723SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 184*77215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 185*77215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 186*77215723SLisandro Dalcin args: -dim 2 -integral 3.141592653589793 # pi 187*77215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 188*77215723SLisandro Dalcin 189*77215723SLisandro Dalcin testset: 190*77215723SLisandro Dalcin suffix: B2_bnd # 2D ball boundary 191*77215723SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 192*77215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 193*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 194*77215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 195*77215723SLisandro Dalcin args: -dim 1 -integral 6.283185307179586 # 2*pi 196*77215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 197*77215723SLisandro Dalcin 198*77215723SLisandro Dalcin 199*77215723SLisandro Dalcin testset: 200*77215723SLisandro Dalcin suffix: B3 # 3D ball 201*77215723SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 202*77215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 203*77215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 204*77215723SLisandro Dalcin args: -dim 3 -integral 4.1887902047863905 # 4/3*pi 205*77215723SLisandro Dalcin args: -order {{2 3 4 5}} -tol 0.20 206*77215723SLisandro Dalcin 207*77215723SLisandro Dalcin testset: 208*77215723SLisandro Dalcin suffix: B3_bnd # 3D ball boundary 209*77215723SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 210*77215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 211*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 212*77215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 213*77215723SLisandro Dalcin args: -dim 2 -integral 12.566370614359172 # 4*pi 214*77215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25 215*77215723SLisandro Dalcin 216*77215723SLisandro Dalcin 217*77215723SLisandro Dalcin testset: 218*77215723SLisandro Dalcin suffix: B_lin # linear discretizations 219*77215723SLisandro Dalcin requires: define(PETSC_GMSH_EXE) 220*77215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 221*77215723SLisandro Dalcin args: -dm_plex_gmsh_highorder true 222*77215723SLisandro Dalcin args: -dm_plex_gmsh_project true 223*77215723SLisandro Dalcin args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output} 224*77215723SLisandro Dalcin args: -dm_plex_gmsh_fe_view 225*77215723SLisandro Dalcin args: -dm_plex_gmsh_project_fe_view 226*77215723SLisandro Dalcin args: -order 1 -tol 1e-4 227*77215723SLisandro Dalcin test: 228*77215723SLisandro Dalcin suffix: dim-1 229*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 230*77215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 231*77215723SLisandro Dalcin args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2) 232*77215723SLisandro Dalcin test: 233*77215723SLisandro Dalcin suffix: dim-2 234*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 235*77215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 236*77215723SLisandro Dalcin args: -dim 2 -integral 2.000000000000000 # 2 237*77215723SLisandro Dalcin test: 238*77215723SLisandro Dalcin suffix: dim-2_msh-B3tet 239*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 240*77215723SLisandro Dalcin args: -msh B3tet -dim 2 -integral 9.914478 241*77215723SLisandro Dalcin test: 242*77215723SLisandro Dalcin suffix: dim-2_msh-B3hex 243*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 244*77215723SLisandro Dalcin args: -msh B3hex -dim 2 -integral 8.000000 245*77215723SLisandro Dalcin test: 246*77215723SLisandro Dalcin suffix: dim-3_msh-B3tet 247*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 248*77215723SLisandro Dalcin args: -msh B3tet -dim 3 -integral 2.666649 249*77215723SLisandro Dalcin test: 250*77215723SLisandro Dalcin suffix: dim-3_msh-B3hex 251*77215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 252*77215723SLisandro Dalcin args: -msh B3hex -dim 3 -integral 1.539600 253*77215723SLisandro Dalcin 254402df9f0SLisandro Dalcin TEST*/ 255