xref: /petsc/src/dm/dt/fe/tests/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1d21efd2eSMatthew G. Knepley static const char help[] = "Tests for determining whether a new finite element works";
2d21efd2eSMatthew G. Knepley 
3d21efd2eSMatthew G. Knepley /*
4d21efd2eSMatthew G. Knepley   Use -interpolation_view and -l2_projection_view to look at the interpolants.
5d21efd2eSMatthew G. Knepley */
6d21efd2eSMatthew G. Knepley 
7d21efd2eSMatthew G. Knepley #include <petscdmplex.h>
8d21efd2eSMatthew G. Knepley #include <petscfe.h>
9d21efd2eSMatthew G. Knepley #include <petscds.h>
10d21efd2eSMatthew G. Knepley #include <petscsnes.h>
11d21efd2eSMatthew G. Knepley 
12d21efd2eSMatthew G. Knepley static void constant(PetscInt dim, PetscInt Nf, PetscInt NfAux,
13d21efd2eSMatthew G. Knepley                      const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
14d21efd2eSMatthew G. Knepley                      const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
15d21efd2eSMatthew G. Knepley                      PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
16d21efd2eSMatthew G. Knepley {
17d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
18d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.;
19d21efd2eSMatthew G. Knepley }
20d21efd2eSMatthew G. Knepley 
21d21efd2eSMatthew G. Knepley static void linear(PetscInt dim, PetscInt Nf, PetscInt NfAux,
22d21efd2eSMatthew G. Knepley                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
23d21efd2eSMatthew G. Knepley                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
24d21efd2eSMatthew G. Knepley                    PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
25d21efd2eSMatthew G. Knepley {
26d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
27d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c];
28d21efd2eSMatthew G. Knepley }
29d21efd2eSMatthew G. Knepley 
30d21efd2eSMatthew G. Knepley static void quadratic(PetscInt dim, PetscInt Nf, PetscInt NfAux,
31d21efd2eSMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
32d21efd2eSMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
33d21efd2eSMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
34d21efd2eSMatthew G. Knepley {
35d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
36d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.*x[c]*x[c];
37d21efd2eSMatthew G. Knepley }
38d21efd2eSMatthew G. Knepley 
39d21efd2eSMatthew G. Knepley static void trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
40d21efd2eSMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
41d21efd2eSMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
42d21efd2eSMatthew G. Knepley                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
43d21efd2eSMatthew G. Knepley {
44d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
45d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += PetscCosReal(2.*PETSC_PI*x[c]);
46d21efd2eSMatthew G. Knepley }
47d21efd2eSMatthew G. Knepley 
48e239af90SMatthew G. Knepley /*
49e239af90SMatthew G. Knepley  The prime basis for the Wheeler-Yotov-Xue prism.
50e239af90SMatthew G. Knepley  */
51e239af90SMatthew G. Knepley static void prime(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52e239af90SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53e239af90SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54e239af90SMatthew G. Knepley                   PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
55e239af90SMatthew G. Knepley {
56e239af90SMatthew G. Knepley   PetscReal x = X[0], y = X[1], z = X[2], b = 1 + x + y + z;
57e239af90SMatthew G. Knepley   f0[0] += b + 2.0*x*z + 2.0*y*z + x*y + x*x;
58e239af90SMatthew G. Knepley   f0[1] += b + 2.0*x*z + 2.0*y*z + x*y + y*y;
59e239af90SMatthew G. Knepley   f0[2] += b - 3.0*x*z - 3.0*y*z -   2.0*z*z;
60e239af90SMatthew G. Knepley }
61e239af90SMatthew G. Knepley 
62e239af90SMatthew G. Knepley static const char    *names[]     = {"constant", "linear", "quadratic", "trig", "prime"};
63e239af90SMatthew G. Knepley static PetscPointFunc functions[] = { constant,   linear,   quadratic,   trig,   prime };
64d21efd2eSMatthew G. Knepley 
65d21efd2eSMatthew G. Knepley typedef struct {
66d21efd2eSMatthew G. Knepley   PetscPointFunc exactSol;
67e239af90SMatthew G. Knepley   PetscReal shear,flatten;
68d21efd2eSMatthew G. Knepley } AppCtx;
69d21efd2eSMatthew G. Knepley 
70d21efd2eSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
71d21efd2eSMatthew G. Knepley {
72d21efd2eSMatthew G. Knepley   char           name[PETSC_MAX_PATH_LEN] = "constant";
73d21efd2eSMatthew G. Knepley   PetscInt       Nfunc = sizeof(names)/sizeof(char *), i;
74d21efd2eSMatthew G. Knepley   PetscErrorCode ierr;
75d21efd2eSMatthew G. Knepley 
76d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
77d21efd2eSMatthew G. Knepley   options->exactSol = NULL;
78e239af90SMatthew G. Knepley   options->shear    = 0.;
79e239af90SMatthew G. Knepley   options->flatten  = 1.;
80d21efd2eSMatthew G. Knepley 
81d21efd2eSMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "FE Test Options", "PETSCFE");CHKERRQ(ierr);
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsString("-func", "Function to project into space", "", name, name, PETSC_MAX_PATH_LEN, NULL));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-shear", "Factor by which to shear along the x-direction", "", options->shear, &(options->shear), NULL));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-flatten", "Factor by which to flatten", "", options->flatten, &(options->flatten), NULL));
85d21efd2eSMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
86d21efd2eSMatthew G. Knepley 
87d21efd2eSMatthew G. Knepley   for (i = 0; i < Nfunc; ++i) {
88d21efd2eSMatthew G. Knepley     PetscBool flg;
89d21efd2eSMatthew G. Knepley 
905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscStrcmp(name, names[i], &flg));
91d21efd2eSMatthew G. Knepley     if (flg) {options->exactSol = functions[i]; break;}
92d21efd2eSMatthew G. Knepley   }
93e239af90SMatthew G. Knepley   PetscCheck(options->exactSol, comm, PETSC_ERR_ARG_WRONG, "Invalid test function %s", name);
94d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
95d21efd2eSMatthew G. Knepley }
96d21efd2eSMatthew G. Knepley 
97d21efd2eSMatthew G. Knepley /* The exact solution is the negative of the f0 contribution */
98d21efd2eSMatthew G. Knepley static PetscErrorCode exactSolution(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
99d21efd2eSMatthew G. Knepley {
100d21efd2eSMatthew G. Knepley   AppCtx  *user    = (AppCtx *) ctx;
101e239af90SMatthew G. Knepley   PetscInt uOff[2] = {0, Nc};
102d21efd2eSMatthew G. Knepley 
103d21efd2eSMatthew G. Knepley   user->exactSol(dim, 1, 0, uOff, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, time, x, 0, NULL, u);
104d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) u[c] *= -1.;
105d21efd2eSMatthew G. Knepley   return 0;
106d21efd2eSMatthew G. Knepley }
107d21efd2eSMatthew G. Knepley 
108d21efd2eSMatthew G. Knepley static void f0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
109d21efd2eSMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
110d21efd2eSMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111d21efd2eSMatthew G. Knepley                PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
112d21efd2eSMatthew G. Knepley {
113d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
114d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) f0[c] += u[c];
115d21efd2eSMatthew G. Knepley }
116d21efd2eSMatthew G. Knepley 
117d21efd2eSMatthew G. Knepley static void g0(PetscInt dim, PetscInt Nf, PetscInt NfAux,
118d21efd2eSMatthew G. Knepley                const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
119d21efd2eSMatthew G. Knepley                const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
120d21efd2eSMatthew G. Knepley                PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
121d21efd2eSMatthew G. Knepley {
122d21efd2eSMatthew G. Knepley   const PetscInt Nc = uOff[1] - uOff[0];
123d21efd2eSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; ++c) g0[c*Nc+c] = 1.0;
124d21efd2eSMatthew G. Knepley }
125d21efd2eSMatthew G. Knepley 
126d21efd2eSMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
127d21efd2eSMatthew G. Knepley {
128d21efd2eSMatthew G. Knepley   PetscDS        ds;
129d21efd2eSMatthew G. Knepley   PetscWeakForm  wf;
130d21efd2eSMatthew G. Knepley 
131d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetWeakForm(ds, &wf));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0, 0, NULL));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, user->exactSol, 0, NULL));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0, 0, NULL, 0, NULL, 0, NULL));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetExactSolution(ds, 0, exactSolution, user));
138d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
139d21efd2eSMatthew G. Knepley }
140d21efd2eSMatthew G. Knepley 
141d21efd2eSMatthew G. Knepley static PetscErrorCode SetupDiscretization(DM dm, const char name[], AppCtx *user)
142d21efd2eSMatthew G. Knepley {
143d21efd2eSMatthew G. Knepley   DM             cdm = dm;
144d21efd2eSMatthew G. Knepley   PetscFE        fe;
145d21efd2eSMatthew G. Knepley   char           prefix[PETSC_MAX_PATH_LEN];
146d21efd2eSMatthew G. Knepley 
147d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateFEDefault(dm, 1, name ? prefix : NULL, -1, &fe));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, name ? name : "Solution"));
151d21efd2eSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupProblem(dm, user));
155d21efd2eSMatthew G. Knepley   while (cdm) {
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, cdm));
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
158d21efd2eSMatthew G. Knepley   }
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
160d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
161d21efd2eSMatthew G. Knepley }
162d21efd2eSMatthew G. Knepley 
163d21efd2eSMatthew G. Knepley /* This test tells us whether the given function is contained in the approximation space */
164d21efd2eSMatthew G. Knepley static PetscErrorCode CheckInterpolation(DM dm, AppCtx *user)
165d21efd2eSMatthew G. Knepley {
166d21efd2eSMatthew G. Knepley   PetscSimplePointFunc exactSol[1];
167d21efd2eSMatthew G. Knepley   void                *exactCtx[1];
168d21efd2eSMatthew G. Knepley   PetscDS              ds;
169d21efd2eSMatthew G. Knepley   Vec                  u;
170d21efd2eSMatthew G. Knepley   PetscReal            error, tol = PETSC_SMALL;
171d21efd2eSMatthew G. Knepley   MPI_Comm             comm;
172d21efd2eSMatthew G. Knepley 
173d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &u));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMProjectFunction(dm, 0.0, exactSol, exactCtx, INSERT_ALL_VALUES, u));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-interpolation_view"));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &u));
1825f80ce2aSJacob Faibussowitsch   if (error > tol) CHKERRQ(PetscPrintf(comm, "Interpolation tests FAIL at tolerance %g error %g\n", (double) tol, (double) error));
1835f80ce2aSJacob Faibussowitsch   else             CHKERRQ(PetscPrintf(comm, "Interpolation tests pass at tolerance %g\n", (double) tol));
184d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
185d21efd2eSMatthew G. Knepley }
186d21efd2eSMatthew G. Knepley 
187d21efd2eSMatthew G. Knepley /* This test tells us whether the element is unisolvent (the mass matrix has full rank), and what rate of convergence we achieve */
188d21efd2eSMatthew G. Knepley static PetscErrorCode CheckL2Projection(DM dm, AppCtx *user)
189d21efd2eSMatthew G. Knepley {
190d21efd2eSMatthew G. Knepley   PetscSimplePointFunc exactSol[1];
191d21efd2eSMatthew G. Knepley   void                *exactCtx[1];
192d21efd2eSMatthew G. Knepley   SNES                 snes;
193d21efd2eSMatthew G. Knepley   PetscDS              ds;
194d21efd2eSMatthew G. Knepley   Vec                  u;
195d21efd2eSMatthew G. Knepley   PetscReal            error, tol = PETSC_SMALL;
196d21efd2eSMatthew G. Knepley   MPI_Comm             comm;
197d21efd2eSMatthew G. Knepley 
198d21efd2eSMatthew G. Knepley   PetscFunctionBeginUser;
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalVector(dm, &u));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm, &snes));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes, dm));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(u, 0.0));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "solution"));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexSetSNESLocalFEM(dm, user, user, user));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSNESCheckFromOptions(snes, u));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes, NULL, u));
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-l2_projection_view"));
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreGlobalVector(dm, &u));
2155f80ce2aSJacob Faibussowitsch   if (error > tol) CHKERRQ(PetscPrintf(comm, "L2 projection tests FAIL at tolerance %g error %g\n", (double) tol, (double) error));
2165f80ce2aSJacob Faibussowitsch   else             CHKERRQ(PetscPrintf(comm, "L2 projection tests pass at tolerance %g\n", (double) tol));
217d21efd2eSMatthew G. Knepley   PetscFunctionReturn(0);
218d21efd2eSMatthew G. Knepley }
219d21efd2eSMatthew G. Knepley 
220e239af90SMatthew G. Knepley /* Distorts the mesh by shearing in the x-direction and flattening, factors provided in the options. */
221e239af90SMatthew G. Knepley static PetscErrorCode DistortMesh(DM dm, AppCtx *user)
222e239af90SMatthew G. Knepley {
223e239af90SMatthew G. Knepley   Vec            coordinates;
224e239af90SMatthew G. Knepley   PetscScalar   *ca;
225e239af90SMatthew G. Knepley   PetscInt       dE, n, i;
226e239af90SMatthew G. Knepley 
227e239af90SMatthew G. Knepley   PetscFunctionBeginUser;
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDim(dm, &dE));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(dm, &coordinates));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(coordinates, &n));
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(coordinates, &ca));
232e239af90SMatthew G. Knepley   for (i = 0; i < (n/dE); ++i) {
233e239af90SMatthew G. Knepley     ca[i*dE+0] += user->shear*ca[i*dE+0];
234e239af90SMatthew G. Knepley     ca[i*dE+1] *= user->flatten;
235e239af90SMatthew G. Knepley   }
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(coordinates, &ca));
237e239af90SMatthew G. Knepley   PetscFunctionReturn(0);
238e239af90SMatthew G. Knepley }
239e239af90SMatthew G. Knepley 
240d21efd2eSMatthew G. Knepley int main(int argc, char **argv)
241d21efd2eSMatthew G. Knepley {
242d21efd2eSMatthew G. Knepley   DM             dm;
243d21efd2eSMatthew G. Knepley   AppCtx         user;
244d21efd2eSMatthew G. Knepley   PetscMPIInt    size;
245d21efd2eSMatthew G. Knepley 
246*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
2475f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
248e239af90SMatthew G. Knepley   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only.");
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(DistortMesh(dm,&user));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-dm_view"));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, NULL, &user));
256d21efd2eSMatthew G. Knepley 
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckInterpolation(dm, &user));
2585f80ce2aSJacob Faibussowitsch   CHKERRQ(CheckL2Projection(dm, &user));
259d21efd2eSMatthew G. Knepley 
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
261*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
262*b122ec5aSJacob Faibussowitsch   return 0;
263d21efd2eSMatthew G. Knepley }
264d21efd2eSMatthew G. Knepley 
265d21efd2eSMatthew G. Knepley /*TEST
266d21efd2eSMatthew G. Knepley 
267d21efd2eSMatthew G. Knepley   testset:
268d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -petscspace_degree 1\
269d21efd2eSMatthew G. Knepley           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
270d21efd2eSMatthew G. Knepley 
271d21efd2eSMatthew G. Knepley     test:
272d21efd2eSMatthew G. Knepley       suffix: p1_0
273d21efd2eSMatthew G. Knepley       args: -func {{constant linear}}
274d21efd2eSMatthew G. Knepley 
275d21efd2eSMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 4 gives convergence rate 2.0
276d21efd2eSMatthew G. Knepley     test:
277d21efd2eSMatthew G. Knepley       suffix: p1_1
278d21efd2eSMatthew G. Knepley       args: -func {{quadratic trig}} \
279d21efd2eSMatthew G. Knepley             -snes_convergence_estimate -convest_num_refine 2
280d21efd2eSMatthew G. Knepley 
281d21efd2eSMatthew G. Knepley   testset:
282e239af90SMatthew G. Knepley     requires: !complex double
283d21efd2eSMatthew G. Knepley     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
284d21efd2eSMatthew G. Knepley             -petscspace_type sum \
285d21efd2eSMatthew G. Knepley             -petscspace_variables 3 \
286d21efd2eSMatthew G. Knepley             -petscspace_components 3 \
287d21efd2eSMatthew G. Knepley             -petscspace_sum_spaces 2 \
288d21efd2eSMatthew G. Knepley             -petscspace_sum_concatenate false \
289d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_variables 3 \
290d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_components 3 \
291d21efd2eSMatthew G. Knepley               -sumcomp_0_petscspace_degree 1 \
292d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_variables 3 \
293d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_components 3 \
294d21efd2eSMatthew G. Knepley               -sumcomp_1_petscspace_type wxy \
295d21efd2eSMatthew G. Knepley             -petscdualspace_form_degree 0 \
296d21efd2eSMatthew G. Knepley             -petscdualspace_order 1 \
297d21efd2eSMatthew G. Knepley             -petscdualspace_components 3 \
298d21efd2eSMatthew G. Knepley           -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
299d21efd2eSMatthew G. Knepley 
300d21efd2eSMatthew G. Knepley     test:
301d21efd2eSMatthew G. Knepley       suffix: wxy_0
302d21efd2eSMatthew G. Knepley       args: -func constant
303d21efd2eSMatthew G. Knepley 
304d21efd2eSMatthew G. Knepley     test:
305d21efd2eSMatthew G. Knepley       suffix: wxy_1
306d21efd2eSMatthew G. Knepley       args: -func linear
307d21efd2eSMatthew G. Knepley 
308e239af90SMatthew G. Knepley     test:
309e239af90SMatthew G. Knepley       suffix: wxy_2
310e239af90SMatthew G. Knepley       args: -func prime
311e239af90SMatthew G. Knepley 
312e239af90SMatthew G. Knepley     test:
313e239af90SMatthew G. Knepley       suffix: wxy_3
314e239af90SMatthew G. Knepley       args: -func linear -shear 1 -flatten 1e-5
315e239af90SMatthew G. Knepley 
316d21efd2eSMatthew G. Knepley TEST*/
317