1 static char help[] = "Tests DMPlex Gmsh reader.\n\n"; 2 3 #include <petscdmplex.h> 4 5 #if !defined(PETSC_GMSH_EXE) 6 #define PETSC_GMSH_EXE "gmsh" 7 #endif 8 9 #include <petscds.h> 10 11 static void one(PetscInt dim, PetscInt Nf, PetscInt NfAux, 12 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 13 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 14 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar value[]) 15 { 16 value[0] = (PetscReal)1; 17 } 18 19 static PetscErrorCode CreateFE(DM dm) 20 { 21 DM cdm; 22 PetscSpace P; 23 PetscDualSpace Q; 24 DM K; 25 PetscFE fe; 26 DMPolytopeType ptype; 27 28 PetscInt dim,k; 29 PetscBool isSimplex; 30 31 PetscDS ds; 32 33 PetscFunctionBeginUser; 34 PetscCall(DMGetCoordinateDM(dm, &cdm)); 35 PetscCall(DMGetField(cdm, 0, NULL, (PetscObject*) &fe)); 36 PetscCall(PetscFEGetBasisSpace(fe, &P)); 37 PetscCall(PetscFEGetDualSpace(fe, &Q)); 38 PetscCall(PetscDualSpaceGetDM(Q,&K)); 39 PetscCall(DMGetDimension(K,&dim)); 40 PetscCall(PetscSpaceGetDegree(P, &k, NULL)); 41 PetscCall(DMPlexGetCellType(K, 0, &ptype)); 42 switch (ptype) { 43 case DM_POLYTOPE_QUADRILATERAL: 44 case DM_POLYTOPE_HEXAHEDRON: 45 isSimplex = PETSC_FALSE; break; 46 default: 47 isSimplex = PETSC_TRUE; break; 48 } 49 50 PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe)); 51 PetscCall(PetscFESetName(fe, "scalar")); 52 PetscCall(DMAddField(dm, NULL, (PetscObject) fe)); 53 PetscCall(PetscFEDestroy(&fe)); 54 PetscCall(DMCreateDS(dm)); 55 56 PetscCall(DMGetDS(dm, &ds)); 57 PetscCall(PetscDSSetObjective(ds, 0, one)); 58 PetscFunctionReturn(0); 59 } 60 61 static PetscErrorCode CheckIntegral(DM dm, PetscReal integral, PetscReal tol) 62 { 63 Vec u; 64 PetscReal rval; 65 PetscScalar result; 66 PetscErrorCode ierr; 67 68 PetscFunctionBeginUser; 69 PetscCall(DMGetGlobalVector(dm, &u)); 70 PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, NULL)); 71 PetscCall(DMRestoreGlobalVector(dm, &u)); 72 rval = PetscRealPart(result); 73 if (integral > 0 && PetscAbsReal(integral - rval) > tol) { 74 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Calculated value %g != %g actual value (error %g > %g tol)\n", 75 (double) rval, (double) integral, (double) PetscAbsReal(integral - rval), (double) tol);PetscCall(ierr); 76 } 77 PetscFunctionReturn(0); 78 } 79 80 int main(int argc, char **argv) 81 { 82 DM dm; 83 const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex", 84 "vtx", "B2tri", "B2qua", "B3tet", "B3hex"}; 85 const char *const fmtlist[] = {"msh22", "msh40", "msh41"}; 86 PetscInt msh = 5; 87 PetscInt fmt = 2; 88 PetscBool bin = PETSC_TRUE; 89 PetscInt dim = 3; 90 PetscInt order = 2; 91 92 const char cmdtemplate[] = "%s -format %s %s -%d -order %d %s -o %s"; 93 char gmsh[PETSC_MAX_PATH_LEN] = PETSC_GMSH_EXE; 94 char tag[PETSC_MAX_PATH_LEN], path[PETSC_MAX_PATH_LEN]; 95 char geo[PETSC_MAX_PATH_LEN], geodir[PETSC_MAX_PATH_LEN] = "."; 96 char out[PETSC_MAX_PATH_LEN], outdir[PETSC_MAX_PATH_LEN] = "."; 97 char cmd[PETSC_MAX_PATH_LEN*4]; 98 PetscBool set,flg; 99 FILE *fp; 100 101 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 102 103 PetscCall(PetscStrncpy(geodir, "${PETSC_DIR}/share/petsc/datafiles/meshes", sizeof(geodir))); 104 PetscCall(PetscOptionsGetenv(PETSC_COMM_SELF, "GMSH", path, sizeof(path), &set)); 105 if (set) PetscCall(PetscStrncpy(gmsh, path, sizeof(gmsh))); 106 PetscCall(PetscOptionsGetString(NULL, NULL, "-gmsh", gmsh, sizeof(gmsh), NULL)); 107 PetscCall(PetscOptionsGetString(NULL, NULL, "-dir", geodir, sizeof(geodir), NULL)); 108 PetscCall(PetscOptionsGetString(NULL, NULL, "-out", outdir, sizeof(outdir), NULL)); 109 PetscCall(PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, PETSC_STATIC_ARRAY_LENGTH(mshlist), &msh, NULL)); 110 PetscCall(PetscOptionsGetEList(NULL, NULL, "-fmt", fmtlist, PETSC_STATIC_ARRAY_LENGTH(fmtlist), &fmt, NULL)); 111 PetscCall(PetscOptionsGetBool(NULL, NULL, "-bin", &bin, NULL)); 112 PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL)); 113 PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL)); 114 if (fmt == 1) bin = PETSC_FALSE; /* Recent Gmsh releases cannot generate msh40+binary format*/ 115 116 { /* This test requires Gmsh >= 4.2.0 */ 117 int inum = 0, major = 0, minor = 0, micro = 0; 118 PetscCall(PetscSNPrintf(cmd, sizeof(cmd), "%s -info", gmsh)); 119 PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 120 if (fp) {inum = fscanf(fp, "Version : %d.%d.%d", &major, &minor, µ);} 121 PetscCall(PetscPClose(PETSC_COMM_SELF, fp)); 122 if (inum != 3 || major < 4 || (major == 4 && minor < 2)) { 123 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Gmsh>=4.2.0 not available\n")); goto finish; 124 } 125 } 126 127 PetscCall(PetscSNPrintf(tag, sizeof(tag), "%s-%d-%d-%s%s", mshlist[msh], (int)dim, (int)order, fmtlist[fmt], bin?"-bin":"")); 128 PetscCall(PetscSNPrintf(geo, sizeof(geo), "%s/gmsh-%s.geo", geodir, mshlist[msh])); 129 PetscCall(PetscSNPrintf(out, sizeof(out), "%s/mesh-%s.msh", outdir, tag)); 130 PetscCall(PetscStrreplace(PETSC_COMM_SELF, geo, path, sizeof(path))); 131 PetscCall(PetscFixFilename(path, geo)); 132 PetscCall(PetscStrreplace(PETSC_COMM_SELF, out, path, sizeof(path))); 133 PetscCall(PetscFixFilename(path, out)); 134 PetscCall(PetscTestFile(geo, 'r', &flg)); 135 PetscCheck(flg,PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "File not found: %s", geo); 136 137 PetscCall(PetscSNPrintf(cmd, sizeof(cmd), cmdtemplate, gmsh, fmtlist[fmt], bin?"-bin":"", (int)dim, (int)order, geo, out)); 138 PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp)); 139 PetscCall(PetscPClose(PETSC_COMM_SELF, fp)); 140 141 PetscCall(DMPlexCreateFromFile(PETSC_COMM_SELF, out, "ex99_plex", PETSC_TRUE, &dm)); 142 PetscCall(PetscSNPrintf(tag, sizeof(tag), "mesh-%s", mshlist[msh])); 143 PetscCall(PetscObjectSetName((PetscObject)dm, tag)); 144 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 145 PetscCall(DMSetFromOptions(dm)); 146 { 147 PetscBool check; 148 PetscReal integral = 0, tol = (PetscReal)1.0e-4; 149 PetscCall(PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check)); 150 PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL)); 151 if (check) { 152 PetscCall(CreateFE(dm)); 153 PetscCall(CheckIntegral(dm, integral, tol)); 154 } 155 } 156 PetscCall(DMDestroy(&dm)); 157 158 finish: 159 PetscCall(PetscFinalize()); 160 return 0; 161 } 162 163 /*TEST 164 165 build: 166 requires: defined(PETSC_HAVE_POPEN) 167 168 test: 169 requires: defined(PETSC_GMSH_EXE) 170 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 171 args: -msh {{vtx}separate_output} 172 args: -order 1 173 args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 174 args: -dm_view ::ascii_info_detail 175 args: -dm_plex_check_all 176 args: -dm_plex_gmsh_highorder false 177 args: -dm_plex_gmsh_use_marker true 178 args: -dm_plex_gmsh_spacedim 3 179 180 test: 181 requires: defined(PETSC_GMSH_EXE) 182 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 183 args: -msh {{seg tri qua tet wed hex}separate_output} 184 args: -order {{1 2 3 7}} 185 args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}} 186 args: -dm_view ::ascii_info_detail 187 args: -dm_plex_check_all 188 args: -dm_plex_gmsh_highorder false 189 args: -dm_plex_gmsh_use_marker true 190 191 testset: 192 suffix: B2 # 2D ball 193 requires: defined(PETSC_GMSH_EXE) 194 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 195 args: -msh {{B2tri B2qua}} 196 args: -dim 2 -integral 3.141592653589793 # pi 197 args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 198 199 testset: 200 suffix: B2_bnd # 2D ball boundary 201 requires: defined(PETSC_GMSH_EXE) 202 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 203 args: -dm_plex_gmsh_spacedim 2 204 args: -msh {{B2tri B2qua}} 205 args: -dim 1 -integral 6.283185307179586 # 2*pi 206 args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05 207 208 testset: 209 suffix: B3 # 3D ball 210 requires: defined(PETSC_GMSH_EXE) 211 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 212 args: -msh {{B3tet B3hex}} 213 args: -dim 3 -integral 4.1887902047863905 # 4/3*pi 214 args: -order {{2 3 4 5}} -tol 0.20 215 216 testset: 217 suffix: B3_bnd # 3D ball boundary 218 requires: defined(PETSC_GMSH_EXE) 219 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 220 args: -dm_plex_gmsh_spacedim 3 221 args: -msh {{B3tet B3hex}} 222 args: -dim 2 -integral 12.566370614359172 # 4*pi 223 args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25 224 225 testset: 226 suffix: B_lin # linear discretizations 227 requires: defined(PETSC_GMSH_EXE) 228 args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes 229 args: -dm_plex_gmsh_highorder true 230 args: -dm_plex_gmsh_project true 231 args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output} 232 args: -dm_plex_gmsh_fe_view 233 args: -dm_plex_gmsh_project_fe_view 234 args: -order 1 -tol 1e-4 235 test: 236 suffix: dim-1 237 args: -dm_plex_gmsh_spacedim 2 238 args: -msh {{B2tri B2qua}separate_output} 239 args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2) 240 test: 241 suffix: dim-2 242 args: -dm_plex_gmsh_spacedim 2 243 args: -msh {{B2tri B2qua}separate_output} 244 args: -dim 2 -integral 2.000000000000000 # 2 245 test: 246 suffix: dim-2_msh-B3tet 247 args: -dm_plex_gmsh_spacedim 3 248 args: -msh B3tet -dim 2 -integral 9.914478 249 test: 250 suffix: dim-2_msh-B3hex 251 args: -dm_plex_gmsh_spacedim 3 252 args: -msh B3hex -dim 2 -integral 8.000000 253 test: 254 suffix: dim-3_msh-B3tet 255 args: -dm_plex_gmsh_spacedim 3 256 args: -msh B3tet -dim 3 -integral 2.666649 257 test: 258 suffix: dim-3_msh-B3hex 259 args: -dm_plex_gmsh_spacedim 3 260 args: -msh B3hex -dim 3 -integral 1.539600 261 262 TEST*/ 263