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 1177215723SLisandro Dalcin static void one(PetscInt dim, PetscInt Nf, PetscInt NfAux, 1277215723SLisandro Dalcin const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 1377215723SLisandro Dalcin const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 1477215723SLisandro Dalcin PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar value[]) 1577215723SLisandro Dalcin { 1677215723SLisandro Dalcin value[0] = (PetscReal)1; 1777215723SLisandro Dalcin } 1877215723SLisandro Dalcin 1977215723SLisandro Dalcin static PetscErrorCode CreateFE(DM dm) 2077215723SLisandro Dalcin { 2177215723SLisandro Dalcin DM cdm; 2277215723SLisandro Dalcin PetscSpace P; 2377215723SLisandro Dalcin PetscDualSpace Q; 2477215723SLisandro Dalcin DM K; 2577215723SLisandro Dalcin PetscFE fe; 2677215723SLisandro Dalcin DMPolytopeType ptype; 2777215723SLisandro Dalcin 2877215723SLisandro Dalcin PetscInt dim,k; 2977215723SLisandro Dalcin PetscBool isSimplex; 3077215723SLisandro Dalcin 3177215723SLisandro Dalcin PetscDS ds; 3277215723SLisandro Dalcin 3377215723SLisandro Dalcin PetscFunctionBeginUser; 345f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetField(cdm, 0, NULL, (PetscObject*) &fe)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetBasisSpace(fe, &P)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEGetDualSpace(fe, &Q)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDualSpaceGetDM(Q,&K)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDimension(K,&dim)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSpaceGetDegree(P, &k, NULL)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCellType(K, 0, &ptype)); 4277215723SLisandro Dalcin switch (ptype) { 4377215723SLisandro Dalcin case DM_POLYTOPE_QUADRILATERAL: 4477215723SLisandro Dalcin case DM_POLYTOPE_HEXAHEDRON: 4577215723SLisandro Dalcin isSimplex = PETSC_FALSE; break; 4677215723SLisandro Dalcin default: 4777215723SLisandro Dalcin isSimplex = PETSC_TRUE; break; 4877215723SLisandro Dalcin } 4977215723SLisandro Dalcin 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFESetName(fe, "scalar")); 525f80ce2aSJacob Faibussowitsch CHKERRQ(DMAddField(dm, NULL, (PetscObject) fe)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFEDestroy(&fe)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDS(dm)); 5577215723SLisandro Dalcin 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetDS(dm, &ds)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDSSetObjective(ds, 0, one)); 5877215723SLisandro Dalcin PetscFunctionReturn(0); 5977215723SLisandro Dalcin } 6077215723SLisandro Dalcin 6177215723SLisandro Dalcin static PetscErrorCode CheckIntegral(DM dm, PetscReal integral, PetscReal tol) 6277215723SLisandro Dalcin { 6377215723SLisandro Dalcin Vec u; 6477215723SLisandro Dalcin PetscReal rval; 6577215723SLisandro Dalcin PetscScalar result; 6677215723SLisandro Dalcin PetscErrorCode ierr; 6777215723SLisandro Dalcin 6877215723SLisandro Dalcin PetscFunctionBeginUser; 695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetGlobalVector(dm, &u)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexComputeIntegralFEM(dm, u, &result, NULL)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreGlobalVector(dm, &u)); 7277215723SLisandro Dalcin rval = PetscRealPart(result); 7377215723SLisandro Dalcin if (integral > 0 && PetscAbsReal(integral - rval) > tol) { 7477215723SLisandro Dalcin ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Calculated value %g != %g actual value (error %g > %g tol)\n", 7577215723SLisandro Dalcin (double) rval, (double) integral, (double) PetscAbsReal(integral - rval), (double) tol);CHKERRQ(ierr); 7677215723SLisandro Dalcin } 7777215723SLisandro Dalcin PetscFunctionReturn(0); 7877215723SLisandro Dalcin } 7977215723SLisandro Dalcin 80402df9f0SLisandro Dalcin int main(int argc, char **argv) 81402df9f0SLisandro Dalcin { 82402df9f0SLisandro Dalcin DM dm; 8377215723SLisandro Dalcin const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex", 84d4886c44SLisandro Dalcin "vtx", "B2tri", "B2qua", "B3tet", "B3hex"}; 85402df9f0SLisandro Dalcin const char *const fmtlist[] = {"msh22", "msh40", "msh41"}; 86402df9f0SLisandro Dalcin PetscInt msh = 5; 87402df9f0SLisandro Dalcin PetscInt fmt = 2; 8877215723SLisandro Dalcin PetscBool bin = PETSC_TRUE; 8977215723SLisandro Dalcin PetscInt dim = 3; 90402df9f0SLisandro Dalcin PetscInt order = 2; 91402df9f0SLisandro Dalcin 9277215723SLisandro Dalcin const char cmdtemplate[] = "%s -format %s %s -%d -order %d %s -o %s"; 93402df9f0SLisandro Dalcin char gmsh[PETSC_MAX_PATH_LEN] = PETSC_GMSH_EXE; 94402df9f0SLisandro Dalcin char tag[PETSC_MAX_PATH_LEN], path[PETSC_MAX_PATH_LEN]; 95402df9f0SLisandro Dalcin char geo[PETSC_MAX_PATH_LEN], geodir[PETSC_MAX_PATH_LEN] = "."; 96402df9f0SLisandro Dalcin char out[PETSC_MAX_PATH_LEN], outdir[PETSC_MAX_PATH_LEN] = "."; 97402df9f0SLisandro Dalcin char cmd[PETSC_MAX_PATH_LEN*4]; 98402df9f0SLisandro Dalcin PetscBool set,flg; 99402df9f0SLisandro Dalcin FILE *fp; 100402df9f0SLisandro Dalcin PetscErrorCode ierr; 101402df9f0SLisandro Dalcin 102402df9f0SLisandro Dalcin ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 103402df9f0SLisandro Dalcin 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(geodir, "${PETSC_DIR}/share/petsc/datafiles/meshes", sizeof(geodir))); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetenv(PETSC_COMM_SELF, "GMSH", path, sizeof(path), &set)); 1065f80ce2aSJacob Faibussowitsch if (set) CHKERRQ(PetscStrncpy(gmsh, path, sizeof(gmsh))); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-gmsh", gmsh, sizeof(gmsh), NULL)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-dir", geodir, sizeof(geodir), NULL)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-out", outdir, sizeof(outdir), NULL)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, (int)(sizeof(mshlist)/sizeof(mshlist[0])), &msh, NULL)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEList(NULL, NULL, "-fmt", fmtlist, (int)(sizeof(fmtlist)/sizeof(fmtlist[0])), &fmt, NULL)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-bin", &bin, NULL)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL)); 115402df9f0SLisandro Dalcin if (fmt == 1) bin = PETSC_FALSE; /* Recent Gmsh releases cannot generate msh40+binary format*/ 116402df9f0SLisandro Dalcin 117402df9f0SLisandro Dalcin { /* This test requires Gmsh >= 4.2.0 */ 118402df9f0SLisandro Dalcin int inum = 0, major = 0, minor = 0, micro = 0; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(cmd, sizeof(cmd), "%s -info", gmsh)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 121402df9f0SLisandro Dalcin if (fp) {inum = fscanf(fp, "Version : %d.%d.%d", &major, &minor, µ);} 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPClose(PETSC_COMM_SELF, fp)); 123402df9f0SLisandro Dalcin if (inum != 3 || major < 4 || (major == 4 && minor < 2)) { 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "Gmsh>=4.2.0 not available\n")); goto finish; 125402df9f0SLisandro Dalcin } 126402df9f0SLisandro Dalcin } 127402df9f0SLisandro Dalcin 1285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(tag, sizeof(tag), "%s-%d-%d-%s%s", mshlist[msh], (int)dim, (int)order, fmtlist[fmt], bin?"-bin":"")); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(geo, sizeof(geo), "%s/gmsh-%s.geo", geodir, mshlist[msh])); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(out, sizeof(out), "%s/mesh-%s.msh", outdir, tag)); 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrreplace(PETSC_COMM_SELF, geo, path, sizeof(path))); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFixFilename(path, geo)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrreplace(PETSC_COMM_SELF, out, path, sizeof(path))); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFixFilename(path, out)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTestFile(geo, 'r', &flg)); 136*28b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "File not found: %s", geo); 137402df9f0SLisandro Dalcin 1385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(cmd, sizeof(cmd), cmdtemplate, gmsh, fmtlist[fmt], bin?"-bin":"", (int)dim, (int)order, geo, out)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPClose(PETSC_COMM_SELF, fp)); 141402df9f0SLisandro Dalcin 1425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromFile(PETSC_COMM_SELF, out, "ex99_plex", PETSC_TRUE, &dm)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(tag, sizeof(tag), "mesh-%s", mshlist[msh])); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)dm, tag)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view")); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(dm)); 14777215723SLisandro Dalcin { 14877215723SLisandro Dalcin PetscBool check; 14977215723SLisandro Dalcin PetscReal integral = 0, tol = (PetscReal)1.0e-4; 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL)); 15277215723SLisandro Dalcin if (check) { 1535f80ce2aSJacob Faibussowitsch CHKERRQ(CreateFE(dm)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(CheckIntegral(dm, integral, tol)); 15577215723SLisandro Dalcin } 15677215723SLisandro Dalcin } 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 158402df9f0SLisandro Dalcin 159402df9f0SLisandro Dalcin finish: 160402df9f0SLisandro Dalcin ierr = PetscFinalize(); 161402df9f0SLisandro Dalcin return ierr; 162402df9f0SLisandro Dalcin } 163402df9f0SLisandro Dalcin 164402df9f0SLisandro Dalcin /*TEST 165402df9f0SLisandro Dalcin 166402df9f0SLisandro Dalcin build: 167dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_POPEN) 168402df9f0SLisandro Dalcin 169402df9f0SLisandro Dalcin test: 170dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 171402df9f0SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 172d4886c44SLisandro Dalcin args: -msh {{vtx}separate_output} 173d4886c44SLisandro Dalcin args: -order 1 174d4886c44SLisandro Dalcin args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 175d4886c44SLisandro Dalcin args: -dm_view ::ascii_info_detail 176d4886c44SLisandro Dalcin args: -dm_plex_check_all 177d4886c44SLisandro Dalcin args: -dm_plex_gmsh_highorder false 178d4886c44SLisandro Dalcin args: -dm_plex_gmsh_use_marker true 179d4886c44SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 180d4886c44SLisandro Dalcin 181d4886c44SLisandro Dalcin test: 182dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 183d4886c44SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 184402df9f0SLisandro Dalcin args: -msh {{seg tri qua tet wed hex}separate_output} 185402df9f0SLisandro Dalcin args: -order {{1 2 3 7}} 186402df9f0SLisandro Dalcin args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 187402df9f0SLisandro Dalcin args: -dm_view ::ascii_info_detail 188402df9f0SLisandro Dalcin args: -dm_plex_check_all 189066ea43fSLisandro Dalcin args: -dm_plex_gmsh_highorder false 1908ec7b1ecSLisandro Dalcin args: -dm_plex_gmsh_use_marker true 191402df9f0SLisandro Dalcin 19277215723SLisandro Dalcin testset: 19377215723SLisandro Dalcin suffix: B2 # 2D ball 194dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 19577215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 19677215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 19777215723SLisandro Dalcin args: -dim 2 -integral 3.141592653589793 # pi 19877215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 19977215723SLisandro Dalcin 20077215723SLisandro Dalcin testset: 20177215723SLisandro Dalcin suffix: B2_bnd # 2D ball boundary 202dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 20377215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 20477215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 20577215723SLisandro Dalcin args: -msh {{B2tri B2qua}} 20677215723SLisandro Dalcin args: -dim 1 -integral 6.283185307179586 # 2*pi 20777215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 20877215723SLisandro Dalcin 20977215723SLisandro Dalcin testset: 21077215723SLisandro Dalcin suffix: B3 # 3D ball 211dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 21277215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 21377215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 21477215723SLisandro Dalcin args: -dim 3 -integral 4.1887902047863905 # 4/3*pi 21577215723SLisandro Dalcin args: -order {{2 3 4 5}} -tol 0.20 21677215723SLisandro Dalcin 21777215723SLisandro Dalcin testset: 21877215723SLisandro Dalcin suffix: B3_bnd # 3D ball boundary 219dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 22077215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 22177215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 22277215723SLisandro Dalcin args: -msh {{B3tet B3hex}} 22377215723SLisandro Dalcin args: -dim 2 -integral 12.566370614359172 # 4*pi 22477215723SLisandro Dalcin args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25 22577215723SLisandro Dalcin 22677215723SLisandro Dalcin testset: 22777215723SLisandro Dalcin suffix: B_lin # linear discretizations 228dfd57a17SPierre Jolivet requires: defined(PETSC_GMSH_EXE) 22977215723SLisandro Dalcin args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 23077215723SLisandro Dalcin args: -dm_plex_gmsh_highorder true 23177215723SLisandro Dalcin args: -dm_plex_gmsh_project true 23277215723SLisandro Dalcin args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output} 23377215723SLisandro Dalcin args: -dm_plex_gmsh_fe_view 23477215723SLisandro Dalcin args: -dm_plex_gmsh_project_fe_view 23577215723SLisandro Dalcin args: -order 1 -tol 1e-4 23677215723SLisandro Dalcin test: 23777215723SLisandro Dalcin suffix: dim-1 23877215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 23977215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 24077215723SLisandro Dalcin args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2) 24177215723SLisandro Dalcin test: 24277215723SLisandro Dalcin suffix: dim-2 24377215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 2 24477215723SLisandro Dalcin args: -msh {{B2tri B2qua}separate_output} 24577215723SLisandro Dalcin args: -dim 2 -integral 2.000000000000000 # 2 24677215723SLisandro Dalcin test: 24777215723SLisandro Dalcin suffix: dim-2_msh-B3tet 24877215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 24977215723SLisandro Dalcin args: -msh B3tet -dim 2 -integral 9.914478 25077215723SLisandro Dalcin test: 25177215723SLisandro Dalcin suffix: dim-2_msh-B3hex 25277215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 25377215723SLisandro Dalcin args: -msh B3hex -dim 2 -integral 8.000000 25477215723SLisandro Dalcin test: 25577215723SLisandro Dalcin suffix: dim-3_msh-B3tet 25677215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 25777215723SLisandro Dalcin args: -msh B3tet -dim 3 -integral 2.666649 25877215723SLisandro Dalcin test: 25977215723SLisandro Dalcin suffix: dim-3_msh-B3hex 26077215723SLisandro Dalcin args: -dm_plex_gmsh_spacedim 3 26177215723SLisandro Dalcin args: -msh B3hex -dim 3 -integral 1.539600 26277215723SLisandro Dalcin 263402df9f0SLisandro Dalcin TEST*/ 264