xref: /petsc/src/dm/impls/plex/tests/ex99.c (revision d0609ced746bc51b019815ca91d747429db24893)
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;
349566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
359566063dSJacob Faibussowitsch   PetscCall(DMGetField(cdm, 0, NULL, (PetscObject*) &fe));
369566063dSJacob Faibussowitsch   PetscCall(PetscFEGetBasisSpace(fe, &P));
379566063dSJacob Faibussowitsch   PetscCall(PetscFEGetDualSpace(fe, &Q));
389566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(Q,&K));
399566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(K,&dim));
409566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(P, &k, NULL));
419566063dSJacob Faibussowitsch   PetscCall(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 
509566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe));
519566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(fe, "scalar"));
529566063dSJacob Faibussowitsch   PetscCall(DMAddField(dm, NULL, (PetscObject) fe));
539566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
549566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
5577215723SLisandro Dalcin 
569566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
579566063dSJacob Faibussowitsch   PetscCall(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 
6777215723SLisandro Dalcin   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &u));
699566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeIntegralFEM(dm, u, &result, NULL));
709566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &u));
7177215723SLisandro Dalcin   rval = PetscRealPart(result);
7277215723SLisandro Dalcin   if (integral > 0 && PetscAbsReal(integral - rval) > tol) {
73*d0609cedSBarry Smith     PetscCall(PetscPrintf(PetscObjectComm((PetscObject) dm), "Calculated value %g != %g actual value (error %g > %g tol)\n",
74*d0609cedSBarry Smith                           (double) rval, (double) integral, (double) PetscAbsReal(integral - rval), (double) tol));
7577215723SLisandro Dalcin   }
7677215723SLisandro Dalcin   PetscFunctionReturn(0);
7777215723SLisandro Dalcin }
7877215723SLisandro Dalcin 
79402df9f0SLisandro Dalcin int main(int argc, char **argv)
80402df9f0SLisandro Dalcin {
81402df9f0SLisandro Dalcin   DM                dm;
8277215723SLisandro Dalcin   const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex",
83d4886c44SLisandro Dalcin                                  "vtx", "B2tri", "B2qua", "B3tet", "B3hex"};
84402df9f0SLisandro Dalcin   const char *const fmtlist[] = {"msh22", "msh40", "msh41"};
85402df9f0SLisandro Dalcin   PetscInt          msh = 5;
86402df9f0SLisandro Dalcin   PetscInt          fmt = 2;
8777215723SLisandro Dalcin   PetscBool         bin = PETSC_TRUE;
8877215723SLisandro Dalcin   PetscInt          dim = 3;
89402df9f0SLisandro Dalcin   PetscInt          order = 2;
90402df9f0SLisandro Dalcin 
9177215723SLisandro Dalcin   const char        cmdtemplate[] = "%s -format %s %s -%d -order %d %s -o %s";
92402df9f0SLisandro Dalcin   char              gmsh[PETSC_MAX_PATH_LEN] = PETSC_GMSH_EXE;
93402df9f0SLisandro Dalcin   char              tag[PETSC_MAX_PATH_LEN], path[PETSC_MAX_PATH_LEN];
94402df9f0SLisandro Dalcin   char              geo[PETSC_MAX_PATH_LEN], geodir[PETSC_MAX_PATH_LEN] = ".";
95402df9f0SLisandro Dalcin   char              out[PETSC_MAX_PATH_LEN], outdir[PETSC_MAX_PATH_LEN] = ".";
96402df9f0SLisandro Dalcin   char              cmd[PETSC_MAX_PATH_LEN*4];
97402df9f0SLisandro Dalcin   PetscBool         set,flg;
98402df9f0SLisandro Dalcin   FILE              *fp;
99402df9f0SLisandro Dalcin 
1009566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
101402df9f0SLisandro Dalcin 
1029566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(geodir, "${PETSC_DIR}/share/petsc/datafiles/meshes", sizeof(geodir)));
1039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetenv(PETSC_COMM_SELF, "GMSH", path, sizeof(path), &set));
1049566063dSJacob Faibussowitsch   if (set) PetscCall(PetscStrncpy(gmsh, path, sizeof(gmsh)));
1059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-gmsh", gmsh, sizeof(gmsh), NULL));
1069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-dir", geodir, sizeof(geodir), NULL));
1079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-out", outdir, sizeof(outdir), NULL));
108dd39110bSPierre Jolivet   PetscCall(PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, PETSC_STATIC_ARRAY_LENGTH(mshlist), &msh, NULL));
109dd39110bSPierre Jolivet   PetscCall(PetscOptionsGetEList(NULL, NULL, "-fmt", fmtlist, PETSC_STATIC_ARRAY_LENGTH(fmtlist), &fmt, NULL));
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-bin", &bin, NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dim", &dim, NULL));
1129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-order", &order, NULL));
113402df9f0SLisandro Dalcin   if (fmt == 1) bin = PETSC_FALSE; /* Recent Gmsh releases cannot generate msh40+binary format*/
114402df9f0SLisandro Dalcin 
115402df9f0SLisandro Dalcin   { /* This test requires Gmsh >= 4.2.0 */
116402df9f0SLisandro Dalcin     int inum = 0, major = 0, minor = 0, micro = 0;
1179566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(cmd, sizeof(cmd), "%s -info", gmsh));
1189566063dSJacob Faibussowitsch     PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
119402df9f0SLisandro Dalcin     if (fp) {inum = fscanf(fp, "Version : %d.%d.%d", &major, &minor, &micro);}
1209566063dSJacob Faibussowitsch     PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
121402df9f0SLisandro Dalcin     if (inum != 3 || major < 4 || (major == 4 && minor < 2)) {
1229566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "Gmsh>=4.2.0 not available\n")); goto finish;
123402df9f0SLisandro Dalcin     }
124402df9f0SLisandro Dalcin   }
125402df9f0SLisandro Dalcin 
1269566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(tag, sizeof(tag), "%s-%d-%d-%s%s", mshlist[msh], (int)dim, (int)order, fmtlist[fmt], bin?"-bin":""));
1279566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(geo, sizeof(geo), "%s/gmsh-%s.geo", geodir, mshlist[msh]));
1289566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(out, sizeof(out), "%s/mesh-%s.msh", outdir, tag));
1299566063dSJacob Faibussowitsch   PetscCall(PetscStrreplace(PETSC_COMM_SELF, geo, path, sizeof(path)));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFixFilename(path, geo));
1319566063dSJacob Faibussowitsch   PetscCall(PetscStrreplace(PETSC_COMM_SELF, out, path, sizeof(path)));
1329566063dSJacob Faibussowitsch   PetscCall(PetscFixFilename(path, out));
1339566063dSJacob Faibussowitsch   PetscCall(PetscTestFile(geo, 'r', &flg));
13428b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_SELF, PETSC_ERR_USER_INPUT, "File not found: %s", geo);
135402df9f0SLisandro Dalcin 
1369566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(cmd, sizeof(cmd), cmdtemplate, gmsh, fmtlist[fmt], bin?"-bin":"", (int)dim, (int)order, geo, out));
1379566063dSJacob Faibussowitsch   PetscCall(PetscPOpen(PETSC_COMM_SELF, NULL, cmd, "r", &fp));
1389566063dSJacob Faibussowitsch   PetscCall(PetscPClose(PETSC_COMM_SELF, fp));
139402df9f0SLisandro Dalcin 
1409566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromFile(PETSC_COMM_SELF, out, "ex99_plex", PETSC_TRUE, &dm));
1419566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(tag, sizeof(tag), "mesh-%s", mshlist[msh]));
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)dm, tag));
1439566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
1449566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
14577215723SLisandro Dalcin   {
14677215723SLisandro Dalcin     PetscBool check;
14777215723SLisandro Dalcin     PetscReal integral = 0, tol = (PetscReal)1.0e-4;
1489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check));
1499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
15077215723SLisandro Dalcin     if (check) {
1519566063dSJacob Faibussowitsch       PetscCall(CreateFE(dm));
1529566063dSJacob Faibussowitsch       PetscCall(CheckIntegral(dm, integral, tol));
15377215723SLisandro Dalcin     }
15477215723SLisandro Dalcin   }
1559566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
156402df9f0SLisandro Dalcin 
157402df9f0SLisandro Dalcin finish:
1589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
159b122ec5aSJacob Faibussowitsch   return 0;
160402df9f0SLisandro Dalcin }
161402df9f0SLisandro Dalcin 
162402df9f0SLisandro Dalcin /*TEST
163402df9f0SLisandro Dalcin 
164402df9f0SLisandro Dalcin   build:
165dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_POPEN)
166402df9f0SLisandro Dalcin 
167402df9f0SLisandro Dalcin   test:
168dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
169402df9f0SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
170d4886c44SLisandro Dalcin     args: -msh {{vtx}separate_output}
171d4886c44SLisandro Dalcin     args: -order 1
172d4886c44SLisandro Dalcin     args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}}
173d4886c44SLisandro Dalcin     args: -dm_view ::ascii_info_detail
174d4886c44SLisandro Dalcin     args: -dm_plex_check_all
175d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_highorder false
176d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_use_marker true
177d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 3
178d4886c44SLisandro Dalcin 
179d4886c44SLisandro Dalcin   test:
180dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
181d4886c44SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
182402df9f0SLisandro Dalcin     args: -msh {{seg tri qua tet wed hex}separate_output}
183402df9f0SLisandro Dalcin     args: -order {{1 2 3 7}}
184402df9f0SLisandro Dalcin     args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}}
185402df9f0SLisandro Dalcin     args: -dm_view ::ascii_info_detail
186402df9f0SLisandro Dalcin     args: -dm_plex_check_all
187066ea43fSLisandro Dalcin     args: -dm_plex_gmsh_highorder false
1888ec7b1ecSLisandro Dalcin     args: -dm_plex_gmsh_use_marker true
189402df9f0SLisandro Dalcin 
19077215723SLisandro Dalcin   testset:
19177215723SLisandro Dalcin     suffix: B2 # 2D ball
192dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
19377215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
19477215723SLisandro Dalcin     args: -msh {{B2tri B2qua}}
19577215723SLisandro Dalcin     args: -dim 2 -integral 3.141592653589793 # pi
19677215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05
19777215723SLisandro Dalcin 
19877215723SLisandro Dalcin   testset:
19977215723SLisandro Dalcin     suffix: B2_bnd # 2D ball boundary
200dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
20177215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
20277215723SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 2
20377215723SLisandro Dalcin     args: -msh {{B2tri B2qua}}
20477215723SLisandro Dalcin     args: -dim 1 -integral 6.283185307179586 # 2*pi
20577215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05
20677215723SLisandro Dalcin 
20777215723SLisandro Dalcin   testset:
20877215723SLisandro Dalcin     suffix: B3 # 3D ball
209dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
21077215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
21177215723SLisandro Dalcin     args: -msh {{B3tet B3hex}}
21277215723SLisandro Dalcin     args: -dim 3 -integral 4.1887902047863905 # 4/3*pi
21377215723SLisandro Dalcin     args: -order {{2 3 4 5}} -tol 0.20
21477215723SLisandro Dalcin 
21577215723SLisandro Dalcin   testset:
21677215723SLisandro Dalcin     suffix: B3_bnd # 3D ball boundary
217dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
21877215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
21977215723SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 3
22077215723SLisandro Dalcin     args: -msh {{B3tet B3hex}}
22177215723SLisandro Dalcin     args: -dim 2 -integral 12.566370614359172 # 4*pi
22277215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25
22377215723SLisandro Dalcin 
22477215723SLisandro Dalcin   testset:
22577215723SLisandro Dalcin     suffix: B_lin # linear discretizations
226dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
22777215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
22877215723SLisandro Dalcin     args: -dm_plex_gmsh_highorder true
22977215723SLisandro Dalcin     args: -dm_plex_gmsh_project true
23077215723SLisandro Dalcin     args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output}
23177215723SLisandro Dalcin     args: -dm_plex_gmsh_fe_view
23277215723SLisandro Dalcin     args: -dm_plex_gmsh_project_fe_view
23377215723SLisandro Dalcin     args: -order 1 -tol 1e-4
23477215723SLisandro Dalcin     test:
23577215723SLisandro Dalcin       suffix: dim-1
23677215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 2
23777215723SLisandro Dalcin       args: -msh {{B2tri B2qua}separate_output}
23877215723SLisandro Dalcin       args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2)
23977215723SLisandro Dalcin     test:
24077215723SLisandro Dalcin       suffix: dim-2
24177215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 2
24277215723SLisandro Dalcin       args: -msh {{B2tri B2qua}separate_output}
24377215723SLisandro Dalcin       args: -dim 2 -integral 2.000000000000000 # 2
24477215723SLisandro Dalcin     test:
24577215723SLisandro Dalcin       suffix: dim-2_msh-B3tet
24677215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
24777215723SLisandro Dalcin       args: -msh B3tet -dim 2 -integral 9.914478
24877215723SLisandro Dalcin     test:
24977215723SLisandro Dalcin       suffix: dim-2_msh-B3hex
25077215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25177215723SLisandro Dalcin       args: -msh B3hex -dim 2 -integral 8.000000
25277215723SLisandro Dalcin     test:
25377215723SLisandro Dalcin       suffix: dim-3_msh-B3tet
25477215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25577215723SLisandro Dalcin       args: -msh B3tet -dim 3 -integral 2.666649
25677215723SLisandro Dalcin     test:
25777215723SLisandro Dalcin       suffix: dim-3_msh-B3hex
25877215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25977215723SLisandro Dalcin       args: -msh B3hex -dim 3 -integral 1.539600
26077215723SLisandro Dalcin 
261402df9f0SLisandro Dalcin TEST*/
262