xref: /petsc/src/dm/impls/plex/tests/ex72.c (revision 6e160c9fdcb0c66b2e49c5cf105016c31def5d46)
1*6e160c9fSMatthew G. Knepley static char help[] = "Tests for geometry models\n\n";
2*6e160c9fSMatthew G. Knepley 
3*6e160c9fSMatthew G. Knepley #include <petscdmplex.h>
4*6e160c9fSMatthew G. Knepley #include <petscds.h>
5*6e160c9fSMatthew G. Knepley 
6*6e160c9fSMatthew G. Knepley typedef struct {
7*6e160c9fSMatthew G. Knepley   PetscReal exactVol; // The exact volume of the shape
8*6e160c9fSMatthew G. Knepley   PetscReal volTol;   // The relative tolerance for checking the volume
9*6e160c9fSMatthew G. Knepley } AppCtx;
10*6e160c9fSMatthew G. Knepley 
11*6e160c9fSMatthew G. Knepley static void identity(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 f0[])
12*6e160c9fSMatthew G. Knepley {
13*6e160c9fSMatthew G. Knepley   f0[0] = 1.0;
14*6e160c9fSMatthew G. Knepley }
15*6e160c9fSMatthew G. Knepley 
16*6e160c9fSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
17*6e160c9fSMatthew G. Knepley {
18*6e160c9fSMatthew G. Knepley   PetscFunctionBegin;
19*6e160c9fSMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
20*6e160c9fSMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
21*6e160c9fSMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
22*6e160c9fSMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
23*6e160c9fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
24*6e160c9fSMatthew G. Knepley }
25*6e160c9fSMatthew G. Knepley 
26*6e160c9fSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
27*6e160c9fSMatthew G. Knepley {
28*6e160c9fSMatthew G. Knepley   options->exactVol = 4. * PETSC_PI / 3.;
29*6e160c9fSMatthew G. Knepley   options->volTol   = PETSC_SMALL;
30*6e160c9fSMatthew G. Knepley 
31*6e160c9fSMatthew G. Knepley   PetscFunctionBeginUser;
32*6e160c9fSMatthew G. Knepley   PetscOptionsBegin(comm, "", "Geometry Model Test Options", "DMPLEX");
33*6e160c9fSMatthew G. Knepley   PetscCall(PetscOptionsReal("-exact_vol", "Exact volume of the shape", __FILE__, options->exactVol, &options->exactVol, NULL));
34*6e160c9fSMatthew G. Knepley   PetscCall(PetscOptionsReal("-vol_tol", "Relative tolerance for checking the volume", __FILE__, options->volTol, &options->volTol, NULL));
35*6e160c9fSMatthew G. Knepley   PetscOptionsEnd();
36*6e160c9fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
37*6e160c9fSMatthew G. Knepley }
38*6e160c9fSMatthew G. Knepley 
39*6e160c9fSMatthew G. Knepley static PetscErrorCode CreateDiscretization(DM dm)
40*6e160c9fSMatthew G. Knepley {
41*6e160c9fSMatthew G. Knepley   PetscFE        fe;
42*6e160c9fSMatthew G. Knepley   PetscDS        ds;
43*6e160c9fSMatthew G. Knepley   DMPolytopeType ct;
44*6e160c9fSMatthew G. Knepley   PetscInt       dim, cStart;
45*6e160c9fSMatthew G. Knepley 
46*6e160c9fSMatthew G. Knepley   PetscFunctionBeginUser;
47*6e160c9fSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
48*6e160c9fSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
49*6e160c9fSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
50*6e160c9fSMatthew G. Knepley   PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, 1, ct, NULL, -1, &fe));
51*6e160c9fSMatthew G. Knepley   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
52*6e160c9fSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
53*6e160c9fSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
54*6e160c9fSMatthew G. Knepley   PetscCall(PetscDSSetObjective(ds, 0, identity));
55*6e160c9fSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe));
56*6e160c9fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
57*6e160c9fSMatthew G. Knepley }
58*6e160c9fSMatthew G. Knepley 
59*6e160c9fSMatthew G. Knepley int main(int argc, char **argv)
60*6e160c9fSMatthew G. Knepley {
61*6e160c9fSMatthew G. Knepley   DM          dm;
62*6e160c9fSMatthew G. Knepley   Vec         u;
63*6e160c9fSMatthew G. Knepley   PetscScalar volume;
64*6e160c9fSMatthew G. Knepley   AppCtx      user;
65*6e160c9fSMatthew G. Knepley 
66*6e160c9fSMatthew G. Knepley   PetscFunctionBeginUser;
67*6e160c9fSMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
68*6e160c9fSMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
69*6e160c9fSMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
70*6e160c9fSMatthew G. Knepley   PetscCall(CreateDiscretization(dm));
71*6e160c9fSMatthew G. Knepley   PetscCall(DMGetGlobalVector(dm, &u));
72*6e160c9fSMatthew G. Knepley   PetscCall(VecSet(u, 0.));
73*6e160c9fSMatthew G. Knepley   PetscCall(DMPlexComputeIntegralFEM(dm, u, &volume, NULL));
74*6e160c9fSMatthew G. Knepley   PetscCall(DMRestoreGlobalVector(dm, &u));
75*6e160c9fSMatthew G. Knepley   PetscCheck(PetscAbsScalar((volume - user.exactVol) / user.exactVol) < user.volTol, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Invalid volume %g != %g", (double)PetscRealPart(volume), (double)user.exactVol);
76*6e160c9fSMatthew G. Knepley   PetscCall(DMDestroy(&dm));
77*6e160c9fSMatthew G. Knepley   PetscCall(PetscFinalize());
78*6e160c9fSMatthew G. Knepley   return 0;
79*6e160c9fSMatthew G. Knepley }
80*6e160c9fSMatthew G. Knepley 
81*6e160c9fSMatthew G. Knepley /*TEST
82*6e160c9fSMatthew G. Knepley 
83*6e160c9fSMatthew G. Knepley   # -dm_refine 6 -vol_tol 6e-4 works
84*6e160c9fSMatthew G. Knepley   test:
85*6e160c9fSMatthew G. Knepley     suffix: ball_0
86*6e160c9fSMatthew G. Knepley     requires: ctetgen
87*6e160c9fSMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_shape ball -dm_refine 3 -dm_geom_model ball -vol_tol 5e-2
88*6e160c9fSMatthew G. Knepley 
89*6e160c9fSMatthew G. Knepley   # -dm_refine 4 -vol_tol 2e-3 works
90*6e160c9fSMatthew G. Knepley   test:
91*6e160c9fSMatthew G. Knepley     suffix: cylinder_0
92*6e160c9fSMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_shape cylinder -dm_plex_cylinder_num_refine 1 \
93*6e160c9fSMatthew G. Knepley           -dm_refine 1 -dm_geom_model cylinder -exact_vol 3.141592653589 -vol_tol 1e-1
94*6e160c9fSMatthew G. Knepley 
95*6e160c9fSMatthew G. Knepley TEST*/
96