xref: /petsc/src/dm/impls/plex/tests/ex99.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
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   PetscErrorCode ierr;
3377215723SLisandro Dalcin 
3477215723SLisandro Dalcin   PetscFunctionBeginUser;
3577215723SLisandro Dalcin   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
361e1ea65dSPierre Jolivet   ierr = DMGetField(cdm, 0, NULL, (PetscObject*) &fe);CHKERRQ(ierr);
3777215723SLisandro Dalcin   ierr = PetscFEGetBasisSpace(fe, &P);CHKERRQ(ierr);
3877215723SLisandro Dalcin   ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr);
3977215723SLisandro Dalcin   ierr = PetscDualSpaceGetDM(Q,&K);CHKERRQ(ierr);
4077215723SLisandro Dalcin   ierr = DMGetDimension(K,&dim);CHKERRQ(ierr);
4177215723SLisandro Dalcin   ierr = PetscSpaceGetDegree(P, &k, NULL);CHKERRQ(ierr);
4277215723SLisandro Dalcin   ierr = DMPlexGetCellType(K, 0, &ptype);CHKERRQ(ierr);
4377215723SLisandro Dalcin   switch (ptype) {
4477215723SLisandro Dalcin   case DM_POLYTOPE_QUADRILATERAL:
4577215723SLisandro Dalcin   case DM_POLYTOPE_HEXAHEDRON:
4677215723SLisandro Dalcin     isSimplex = PETSC_FALSE; break;
4777215723SLisandro Dalcin   default:
4877215723SLisandro Dalcin     isSimplex = PETSC_TRUE; break;
4977215723SLisandro Dalcin   }
5077215723SLisandro Dalcin 
5177215723SLisandro Dalcin   ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, isSimplex, k, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
5277215723SLisandro Dalcin   ierr = PetscFESetName(fe, "scalar");CHKERRQ(ierr);
5377215723SLisandro Dalcin   ierr = DMAddField(dm, NULL, (PetscObject) fe);CHKERRQ(ierr);
5477215723SLisandro Dalcin   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
5577215723SLisandro Dalcin   ierr = DMCreateDS(dm);CHKERRQ(ierr);
5677215723SLisandro Dalcin 
5777215723SLisandro Dalcin   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
5877215723SLisandro Dalcin   ierr = PetscDSSetObjective(ds, 0, one);CHKERRQ(ierr);
5977215723SLisandro Dalcin   PetscFunctionReturn(0);
6077215723SLisandro Dalcin }
6177215723SLisandro Dalcin 
6277215723SLisandro Dalcin static PetscErrorCode CheckIntegral(DM dm, PetscReal integral, PetscReal tol)
6377215723SLisandro Dalcin {
6477215723SLisandro Dalcin   Vec            u;
6577215723SLisandro Dalcin   PetscReal      rval;
6677215723SLisandro Dalcin   PetscScalar    result;
6777215723SLisandro Dalcin   PetscErrorCode ierr;
6877215723SLisandro Dalcin 
6977215723SLisandro Dalcin   PetscFunctionBeginUser;
7077215723SLisandro Dalcin   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
7177215723SLisandro Dalcin   ierr = DMPlexComputeIntegralFEM(dm, u, &result, NULL);CHKERRQ(ierr);
7277215723SLisandro Dalcin   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
7377215723SLisandro Dalcin   rval = PetscRealPart(result);
7477215723SLisandro Dalcin   if (integral > 0 && PetscAbsReal(integral - rval) > tol) {
7577215723SLisandro Dalcin     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Calculated value %g != %g actual value (error %g > %g tol)\n",
7677215723SLisandro Dalcin                        (double) rval, (double) integral, (double) PetscAbsReal(integral - rval), (double) tol);CHKERRQ(ierr);
7777215723SLisandro Dalcin   }
7877215723SLisandro Dalcin   PetscFunctionReturn(0);
7977215723SLisandro Dalcin }
8077215723SLisandro Dalcin 
81402df9f0SLisandro Dalcin int main(int argc, char **argv)
82402df9f0SLisandro Dalcin {
83402df9f0SLisandro Dalcin   DM                dm;
8477215723SLisandro Dalcin   const char *const mshlist[] = {"seg", "tri", "qua", "tet", "wed", "hex",
85d4886c44SLisandro Dalcin                                  "vtx", "B2tri", "B2qua", "B3tet", "B3hex"};
86402df9f0SLisandro Dalcin   const char *const fmtlist[] = {"msh22", "msh40", "msh41"};
87402df9f0SLisandro Dalcin   PetscInt          msh = 5;
88402df9f0SLisandro Dalcin   PetscInt          fmt = 2;
8977215723SLisandro Dalcin   PetscBool         bin = PETSC_TRUE;
9077215723SLisandro Dalcin   PetscInt          dim = 3;
91402df9f0SLisandro Dalcin   PetscInt          order = 2;
92402df9f0SLisandro Dalcin 
9377215723SLisandro 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);
11177215723SLisandro Dalcin   ierr = PetscOptionsGetEList(NULL, NULL, "-msh", mshlist, (int)(sizeof(mshlist)/sizeof(mshlist[0])), &msh, NULL);CHKERRQ(ierr);
11277215723SLisandro 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);
11477215723SLisandro 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, &micro);}
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 
12977215723SLisandro 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 
13977215723SLisandro 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);
14877215723SLisandro Dalcin   {
14977215723SLisandro Dalcin     PetscBool check;
15077215723SLisandro Dalcin     PetscReal integral = 0, tol = (PetscReal)1.0e-4;
15177215723SLisandro Dalcin     ierr = PetscOptionsGetReal(NULL, NULL, "-integral", &integral, &check);CHKERRQ(ierr);
15277215723SLisandro Dalcin     ierr = PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL);CHKERRQ(ierr);
15377215723SLisandro Dalcin     if (check) {
15477215723SLisandro Dalcin       ierr = CreateFE(dm);CHKERRQ(ierr);
15577215723SLisandro Dalcin       ierr = CheckIntegral(dm, integral, tol);CHKERRQ(ierr);
15677215723SLisandro Dalcin     }
15777215723SLisandro 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:
168*dfd57a17SPierre Jolivet     requires: defined(PETSC_HAVE_POPEN)
169402df9f0SLisandro Dalcin 
170402df9f0SLisandro Dalcin   test:
171*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
172402df9f0SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
173d4886c44SLisandro Dalcin     args: -msh {{vtx}separate_output}
174d4886c44SLisandro Dalcin     args: -order 1
175d4886c44SLisandro Dalcin     args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}}
176d4886c44SLisandro Dalcin     args: -dm_view ::ascii_info_detail
177d4886c44SLisandro Dalcin     args: -dm_plex_check_all
178d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_highorder false
179d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_use_marker true
180d4886c44SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 3
181d4886c44SLisandro Dalcin 
182d4886c44SLisandro Dalcin   test:
183*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
184d4886c44SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
185402df9f0SLisandro Dalcin     args: -msh {{seg tri qua tet wed hex}separate_output}
186402df9f0SLisandro Dalcin     args: -order {{1 2 3 7}}
187402df9f0SLisandro Dalcin     args: -fmt {{msh22 msh40 msh41}} -bin {{0 1}}
188402df9f0SLisandro Dalcin     args: -dm_view ::ascii_info_detail
189402df9f0SLisandro Dalcin     args: -dm_plex_check_all
190066ea43fSLisandro Dalcin     args: -dm_plex_gmsh_highorder false
1918ec7b1ecSLisandro Dalcin     args: -dm_plex_gmsh_use_marker true
192402df9f0SLisandro Dalcin 
19377215723SLisandro Dalcin   testset:
19477215723SLisandro Dalcin     suffix: B2 # 2D ball
195*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
19677215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
19777215723SLisandro Dalcin     args: -msh {{B2tri B2qua}}
19877215723SLisandro Dalcin     args: -dim 2 -integral 3.141592653589793 # pi
19977215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05
20077215723SLisandro Dalcin 
20177215723SLisandro Dalcin   testset:
20277215723SLisandro Dalcin     suffix: B2_bnd # 2D ball boundary
203*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
20477215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
20577215723SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 2
20677215723SLisandro Dalcin     args: -msh {{B2tri B2qua}}
20777215723SLisandro Dalcin     args: -dim 1 -integral 6.283185307179586 # 2*pi
20877215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.05
20977215723SLisandro Dalcin 
21077215723SLisandro Dalcin   testset:
21177215723SLisandro Dalcin     suffix: B3 # 3D ball
212*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
21377215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
21477215723SLisandro Dalcin     args: -msh {{B3tet B3hex}}
21577215723SLisandro Dalcin     args: -dim 3 -integral 4.1887902047863905 # 4/3*pi
21677215723SLisandro Dalcin     args: -order {{2 3 4 5}} -tol 0.20
21777215723SLisandro Dalcin 
21877215723SLisandro Dalcin   testset:
21977215723SLisandro Dalcin     suffix: B3_bnd # 3D ball boundary
220*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
22177215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
22277215723SLisandro Dalcin     args: -dm_plex_gmsh_spacedim 3
22377215723SLisandro Dalcin     args: -msh {{B3tet B3hex}}
22477215723SLisandro Dalcin     args: -dim 2 -integral 12.566370614359172 # 4*pi
22577215723SLisandro Dalcin     args: -order {{2 3 4 5 6 7 8 9}} -tol 0.25
22677215723SLisandro Dalcin 
22777215723SLisandro Dalcin   testset:
22877215723SLisandro Dalcin     suffix: B_lin # linear discretizations
229*dfd57a17SPierre Jolivet     requires: defined(PETSC_GMSH_EXE)
23077215723SLisandro Dalcin     args: -dir ${wPETSC_DIR}/share/petsc/datafiles/meshes
23177215723SLisandro Dalcin     args: -dm_plex_gmsh_highorder true
23277215723SLisandro Dalcin     args: -dm_plex_gmsh_project true
23377215723SLisandro Dalcin     args: -dm_plex_gmsh_project_petscspace_degree {{1 2 3}separate_output}
23477215723SLisandro Dalcin     args: -dm_plex_gmsh_fe_view
23577215723SLisandro Dalcin     args: -dm_plex_gmsh_project_fe_view
23677215723SLisandro Dalcin     args: -order 1 -tol 1e-4
23777215723SLisandro Dalcin     test:
23877215723SLisandro Dalcin       suffix: dim-1
23977215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 2
24077215723SLisandro Dalcin       args: -msh {{B2tri B2qua}separate_output}
24177215723SLisandro Dalcin       args: -dim 1 -integral 5.656854249492381 # 4*sqrt(2)
24277215723SLisandro Dalcin     test:
24377215723SLisandro Dalcin       suffix: dim-2
24477215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 2
24577215723SLisandro Dalcin       args: -msh {{B2tri B2qua}separate_output}
24677215723SLisandro Dalcin       args: -dim 2 -integral 2.000000000000000 # 2
24777215723SLisandro Dalcin     test:
24877215723SLisandro Dalcin       suffix: dim-2_msh-B3tet
24977215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25077215723SLisandro Dalcin       args: -msh B3tet -dim 2 -integral 9.914478
25177215723SLisandro Dalcin     test:
25277215723SLisandro Dalcin       suffix: dim-2_msh-B3hex
25377215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25477215723SLisandro Dalcin       args: -msh B3hex -dim 2 -integral 8.000000
25577215723SLisandro Dalcin     test:
25677215723SLisandro Dalcin       suffix: dim-3_msh-B3tet
25777215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
25877215723SLisandro Dalcin       args: -msh B3tet -dim 3 -integral 2.666649
25977215723SLisandro Dalcin     test:
26077215723SLisandro Dalcin       suffix: dim-3_msh-B3hex
26177215723SLisandro Dalcin       args: -dm_plex_gmsh_spacedim 3
26277215723SLisandro Dalcin       args: -msh B3hex -dim 3 -integral 1.539600
26377215723SLisandro Dalcin 
264402df9f0SLisandro Dalcin TEST*/
265