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 977215723SLisandro Dalcin #include <petscds.h> 1077215723SLisandro Dalcin 11*9371c9d4SSatish Balay static void one(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar value[]) { 1277215723SLisandro Dalcin value[0] = (PetscReal)1; 1377215723SLisandro Dalcin } 1477215723SLisandro Dalcin 15*9371c9d4SSatish Balay static PetscErrorCode CreateFE(DM dm) { 1677215723SLisandro Dalcin DM cdm; 1777215723SLisandro Dalcin PetscSpace P; 1877215723SLisandro Dalcin PetscDualSpace Q; 1977215723SLisandro Dalcin DM K; 2077215723SLisandro Dalcin PetscFE fe; 2177215723SLisandro Dalcin DMPolytopeType ptype; 2277215723SLisandro Dalcin 2377215723SLisandro Dalcin PetscInt dim, k; 2477215723SLisandro Dalcin PetscBool isSimplex; 2577215723SLisandro Dalcin 2677215723SLisandro Dalcin PetscDS ds; 2777215723SLisandro Dalcin 2877215723SLisandro Dalcin PetscFunctionBeginUser; 299566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 309566063dSJacob Faibussowitsch PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe)); 319566063dSJacob Faibussowitsch PetscCall(PetscFEGetBasisSpace(fe, &P)); 329566063dSJacob Faibussowitsch PetscCall(PetscFEGetDualSpace(fe, &Q)); 339566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(Q, &K)); 349566063dSJacob Faibussowitsch PetscCall(DMGetDimension(K, &dim)); 359566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &k, NULL)); 369566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(K, 0, &ptype)); 3777215723SLisandro Dalcin switch (ptype) { 3877215723SLisandro Dalcin case DM_POLYTOPE_QUADRILATERAL: 39*9371c9d4SSatish Balay case DM_POLYTOPE_HEXAHEDRON: isSimplex = PETSC_FALSE; break; 40*9371c9d4SSatish Balay default: isSimplex = PETSC_TRUE; break; 4177215723SLisandro Dalcin } 4277215723SLisandro Dalcin 439566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe)); 449566063dSJacob Faibussowitsch PetscCall(PetscFESetName(fe, "scalar")); 459566063dSJacob Faibussowitsch PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 469566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 479566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 4877215723SLisandro Dalcin 499566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds)); 509566063dSJacob Faibussowitsch PetscCall(PetscDSSetObjective(ds, 0, one)); 5177215723SLisandro Dalcin PetscFunctionReturn(0); 5277215723SLisandro Dalcin } 5377215723SLisandro Dalcin 54*9371c9d4SSatish Balay static PetscErrorCode CheckIntegral(DM dm, PetscReal integral, PetscReal tol) { 5577215723SLisandro Dalcin Vec u; 5677215723SLisandro Dalcin PetscReal rval; 5777215723SLisandro Dalcin PetscScalar result; 5877215723SLisandro Dalcin 5977215723SLisandro Dalcin PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 619566063dSJacob Faibussowitsch PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, NULL)); 629566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 6377215723SLisandro Dalcin rval = PetscRealPart(result); 6477215723SLisandro Dalcin if (integral > 0 && PetscAbsReal(integral - rval) > tol) { 65*9371c9d4SSatish Balay PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm), "Calculated value %g != %g actual value (error %g > %g tol)\n", (double)rval, (double)integral, (double)PetscAbsReal(integral - rval), (double)tol)); 6677215723SLisandro Dalcin } 6777215723SLisandro Dalcin PetscFunctionReturn(0); 6877215723SLisandro Dalcin } 6977215723SLisandro Dalcin 70*9371c9d4SSatish Balay int main(int argc, char **argv) { 71402df9f0SLisandro Dalcin DM dm; 72*9371c9d4SSatish Balay const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex", "vtx", "B2tri", "B2qua", "B3tet", "B3hex"}; 73402df9f0SLisandro Dalcin const char *const fmtlist[] = {"msh22", "msh40", "msh41"}; 74402df9f0SLisandro Dalcin PetscInt msh = 5; 75402df9f0SLisandro Dalcin PetscInt fmt = 2; 7677215723SLisandro Dalcin PetscBool bin = PETSC_TRUE; 7777215723SLisandro Dalcin PetscInt dim = 3; 78402df9f0SLisandro Dalcin PetscInt order = 2; 79402df9f0SLisandro Dalcin 8077215723SLisandro Dalcin const char cmdtemplate[] = "%s -format %s %s -%d -order %d %s -o %s"; 81402df9f0SLisandro Dalcin char gmsh[PETSC_MAX_PATH_LEN] = PETSC_GMSH_EXE; 82402df9f0SLisandro Dalcin char tag[PETSC_MAX_PATH_LEN], path[PETSC_MAX_PATH_LEN]; 83402df9f0SLisandro Dalcin char geo[PETSC_MAX_PATH_LEN], geodir[PETSC_MAX_PATH_LEN] = "."; 84402df9f0SLisandro Dalcin char out[PETSC_MAX_PATH_LEN], outdir[PETSC_MAX_PATH_LEN] = "."; 85402df9f0SLisandro Dalcin char cmd[PETSC_MAX_PATH_LEN * 4]; 86402df9f0SLisandro Dalcin PetscBool set, flg; 87402df9f0SLisandro Dalcin FILE *fp; 88402df9f0SLisandro Dalcin 89327415f7SBarry Smith PetscFunctionBeginUser; 909566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 91402df9f0SLisandro Dalcin 929566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(geodir, "${PETSC_DIR}/share/petsc/datafiles/meshes", sizeof(geodir))); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetenv(PETSC_COMM_SELF, "GMSH", path, sizeof(path), &set)); 949566063dSJacob Faibussowitsch if (set) PetscCall(PetscStrncpy(gmsh, path, sizeof(gmsh))); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-gmsh", gmsh, sizeof(gmsh), NULL)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-dir", geodir, sizeof(geodir), NULL)); 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-out", outdir, sizeof(outdir), NULL)); 98dd39110bSPierre Jolivet PetscCall(PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, PETSC_STATIC_ARRAY_LENGTH(mshlist), &msh, NULL)); 99dd39110bSPierre Jolivet PetscCall(PetscOptionsGetEList(NULL, NULL, "-fmt", fmtlist, PETSC_STATIC_ARRAY_LENGTH(fmtlist), &fmt, NULL)); 1009566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-bin", &bin, NULL)); 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL)); 103402df9f0SLisandro Dalcin if (fmt == 1) bin = PETSC_FALSE; /* Recent Gmsh releases cannot generate msh40+binary format*/ 104402df9f0SLisandro Dalcin 105402df9f0SLisandro Dalcin { /* This test requires Gmsh >= 4.2.0 */ 106c1599abfSMatthew G. Knepley char space[PETSC_MAX_PATH_LEN]; 107402df9f0SLisandro Dalcin int inum = 0, major = 0, minor = 0, micro = 0; 1089566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(cmd, sizeof(cmd), "%s -info", gmsh)); 1099566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 110c1599abfSMatthew G. Knepley if (fp) { inum = fscanf(fp, "Version %s %d.%d.%d", space, &major, &minor, µ); } 1119566063dSJacob Faibussowitsch PetscCall(PetscPClose(PETSC_COMM_SELF, fp)); 112c1599abfSMatthew G. Knepley if (inum != 4 || major < 4 || (major == 4 && minor < 2)) { 113*9371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_SELF, "Gmsh>=4.2.0 not available\n")); 114*9371c9d4SSatish Balay goto finish; 115402df9f0SLisandro Dalcin } 116402df9f0SLisandro Dalcin } 117402df9f0SLisandro Dalcin 1189566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tag, sizeof(tag), "%s-%d-%d-%s%s", mshlist[msh], (int)dim, (int)order, fmtlist[fmt], bin ? "-bin" : "")); 1199566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(geo, sizeof(geo), "%s/gmsh-%s.geo", geodir, mshlist[msh])); 1209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(out, sizeof(out), "%s/mesh-%s.msh", outdir, tag)); 1219566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_SELF, geo, path, sizeof(path))); 1229566063dSJacob Faibussowitsch PetscCall(PetscFixFilename(path, geo)); 1239566063dSJacob Faibussowitsch PetscCall(PetscStrreplace(PETSC_COMM_SELF, out, path, sizeof(path))); 1249566063dSJacob Faibussowitsch PetscCall(PetscFixFilename(path, out)); 1259566063dSJacob Faibussowitsch PetscCall(PetscTestFile(geo, 'r', &flg)); 12628b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "File not found: %s", geo); 127402df9f0SLisandro Dalcin 1289566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(cmd, sizeof(cmd), cmdtemplate, gmsh, fmtlist[fmt], bin ? "-bin" : "", (int)dim, (int)order, geo, out)); 1299566063dSJacob Faibussowitsch PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 1309566063dSJacob Faibussowitsch PetscCall(PetscPClose(PETSC_COMM_SELF, fp)); 131402df9f0SLisandro Dalcin 1329566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromFile(PETSC_COMM_SELF, out, "ex99_plex", PETSC_TRUE, &dm)); 1339566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tag, sizeof(tag), "mesh-%s", mshlist[msh])); 1349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, tag)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); 136c1599abfSMatthew G. Knepley PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 13777215723SLisandro Dalcin { 13877215723SLisandro Dalcin PetscBool check; 13977215723SLisandro Dalcin PetscReal integral = 0, tol = (PetscReal)1.0e-4; 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL)); 14277215723SLisandro Dalcin if (check) { 1439566063dSJacob Faibussowitsch PetscCall(CreateFE(dm)); 1449566063dSJacob Faibussowitsch PetscCall(CheckIntegral(dm, integral, tol)); 14577215723SLisandro Dalcin } 14677215723SLisandro Dalcin } 1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 148402df9f0SLisandro Dalcin 149402df9f0SLisandro Dalcin finish: 1509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 151b122ec5aSJacob Faibussowitsch return 0; 152402df9f0SLisandro Dalcin } 153402df9f0SLisandro Dalcin 154402df9f0SLisandro Dalcin /*TEST 155402df9f0SLisandro Dalcin 156402df9f0SLisandro Dalcin build: 157dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_POPEN) 158402df9f0SLisandro Dalcin 159402df9f0SLisandro Dalcin test: 160dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 161402df9f0SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 162d4886c44SLisandro Dalcin args: -msh {{vtx}separate_output} 163d4886c44SLisandro Dalcin args: -order 1 164d4886c44SLisandro Dalcin args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 165d4886c44SLisandro Dalcin args: -dm_view ::ascii_info_detail 166d4886c44SLisandro Dalcin args: -dm_plex_check_all 167d4886c44SLisandro Dalcin args: -dm_plex_gmsh_highorder false 168c1599abfSMatthew G. Knepley args: -dm_plex_boundary_label marker 169d4886c44SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 170d4886c44SLisandro Dalcin 171d4886c44SLisandro Dalcin test: 172dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 173d4886c44SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 174402df9f0SLisandro Dalcin args: -msh {{seg tri qua tet wed hex}separate_output} 175402df9f0SLisandro Dalcin args: -order {{1 2 3 7}} 176402df9f0SLisandro Dalcin args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 177402df9f0SLisandro Dalcin args: -dm_view ::ascii_info_detail 178402df9f0SLisandro Dalcin args: -dm_plex_check_all 179066ea43fSLisandro Dalcin args: -dm_plex_gmsh_highorder false 180c1599abfSMatthew G. Knepley args: -dm_plex_boundary_label marker 181402df9f0SLisandro Dalcin 18277215723SLisandro Dalcin testset: 18377215723SLisandro Dalcin suffix: B2 # 2D ball 184dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 18577215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 18677215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 18777215723SLisandro Dalcin args: -dim 2 -integral 3.141592653589793 # pi 18877215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 18977215723SLisandro Dalcin 19077215723SLisandro Dalcin testset: 19177215723SLisandro Dalcin suffix: B2_bnd # 2D ball boundary 192dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 19377215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 19477215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 19577215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 19677215723SLisandro Dalcin args: -dim 1 -integral 6.283185307179586 # 2*pi 19777215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 19877215723SLisandro Dalcin 19977215723SLisandro Dalcin testset: 20077215723SLisandro Dalcin suffix: B3 # 3D ball 201dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 20277215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 20377215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 20477215723SLisandro Dalcin args: -dim 3 -integral 4.1887902047863905 # 4/3*pi 20577215723SLisandro Dalcin args: -order {{2 3 4 5}} -tol 0.20 20677215723SLisandro Dalcin 20777215723SLisandro Dalcin testset: 20877215723SLisandro Dalcin suffix: B3_bnd # 3D ball boundary 209dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 21077215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 21177215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 21277215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 21377215723SLisandro Dalcin args: -dim 2 -integral 12.566370614359172 # 4*pi 21477215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25 21577215723SLisandro Dalcin 21677215723SLisandro Dalcin testset: 21777215723SLisandro Dalcin suffix: B_lin # linear discretizations 218dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 21977215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 22077215723SLisandro Dalcin args: -dm_plex_gmsh_highorder true 22177215723SLisandro Dalcin args: -dm_plex_gmsh_project true 22277215723SLisandro Dalcin args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output} 22377215723SLisandro Dalcin args: -dm_plex_gmsh_fe_view 22477215723SLisandro Dalcin args: -dm_plex_gmsh_project_fe_view 22577215723SLisandro Dalcin args: -order 1 -tol 1e-4 22677215723SLisandro Dalcin test: 22777215723SLisandro Dalcin suffix: dim-1 22877215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 22977215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 23077215723SLisandro Dalcin args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2) 23177215723SLisandro Dalcin test: 23277215723SLisandro Dalcin suffix: dim-2 23377215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 23477215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 23577215723SLisandro Dalcin args: -dim 2 -integral 2.000000000000000 # 2 23677215723SLisandro Dalcin test: 23777215723SLisandro Dalcin suffix: dim-2_msh-B3tet 23877215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 23977215723SLisandro Dalcin args: -msh B3tet -dim 2 -integral 9.914478 24077215723SLisandro Dalcin test: 24177215723SLisandro Dalcin suffix: dim-2_msh-B3hex 24277215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 24377215723SLisandro Dalcin args: -msh B3hex -dim 2 -integral 8.000000 24477215723SLisandro Dalcin test: 24577215723SLisandro Dalcin suffix: dim-3_msh-B3tet 24677215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 24777215723SLisandro Dalcin args: -msh B3tet -dim 3 -integral 2.666649 24877215723SLisandro Dalcin test: 24977215723SLisandro Dalcin suffix: dim-3_msh-B3hex 25077215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 25177215723SLisandro Dalcin args: -msh B3hex -dim 3 -integral 1.539600 25277215723SLisandro Dalcin 253402df9f0SLisandro Dalcin TEST*/ 254