xref: /petsc/src/snes/tutorials/ex24.c (revision 910b42cbfb1b6524ef4594118dd5e013aee7da5d)
1c4762a1bSJed Brown static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the Poisson problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown This example supports automatic convergence estimation\n\
5c4762a1bSJed Brown and Hdiv elements.\n\n\n";
6c4762a1bSJed Brown 
769b95281SMatthew G. Knepley /*
869b95281SMatthew G. Knepley 
969b95281SMatthew G. Knepley The mixed form of Poisson's equation, e.g. https://firedrakeproject.org/demos/poisson_mixed.py.html, is given
1069b95281SMatthew G. Knepley in the strong form by
1169b95281SMatthew G. Knepley \begin{align}
1269b95281SMatthew G. Knepley   \vb{q} - \nabla u   &= 0 \\
1369b95281SMatthew G. Knepley   \nabla \cdot \vb{q} &= f
1469b95281SMatthew G. Knepley \end{align}
1569b95281SMatthew G. Knepley where $u$ is the potential, as in the original problem, but we also discretize the gradient of potential $\vb{q}$,
1669b95281SMatthew G. Knepley or flux, directly. The weak form is then
1769b95281SMatthew G. Knepley \begin{align}
1869b95281SMatthew G. Knepley   <t, \vb{q}> + <\nabla \cdot t, u> - <t_n, u>_\Gamma &= 0 \\
1969b95281SMatthew G. Knepley   <v, \nabla \cdot \vb{q}> &= <v, f>
2069b95281SMatthew G. Knepley \end{align}
2169b95281SMatthew G. Knepley For the original Poisson problem, the Dirichlet boundary forces an essential boundary condition on the potential space,
2269b95281SMatthew G. Knepley and the Neumann boundary gives a natural boundary condition in the weak form. In the mixed formulation, the Neumann
2369b95281SMatthew G. Knepley boundary gives an essential boundary condition on the flux space, $\vb{q} \cdot \vu{n} = h$, and the Dirichlet condition
2469b95281SMatthew G. Knepley becomes a natural condition in the weak form, <t_n, g>_\Gamma.
2569b95281SMatthew G. Knepley */
2669b95281SMatthew G. Knepley 
27c4762a1bSJed Brown #include <petscdmplex.h>
28c4762a1bSJed Brown #include <petscsnes.h>
29c4762a1bSJed Brown #include <petscds.h>
30c4762a1bSJed Brown #include <petscconvest.h>
31c4762a1bSJed Brown 
329371c9d4SSatish Balay typedef enum {
339371c9d4SSatish Balay   SOL_LINEAR,
349371c9d4SSatish Balay   SOL_QUADRATIC,
359371c9d4SSatish Balay   SOL_QUARTIC,
369371c9d4SSatish Balay   SOL_QUARTIC_NEUMANN,
379371c9d4SSatish Balay   SOL_UNKNOWN,
389371c9d4SSatish Balay   NUM_SOLTYPE
399371c9d4SSatish Balay } SolType;
4069b95281SMatthew G. Knepley const char *SolTypeNames[NUM_SOLTYPE + 3] = {"linear", "quadratic", "quartic", "quartic_neumann", "unknown", "SolType", "SOL_", NULL};
41c4762a1bSJed Brown 
42c4762a1bSJed Brown typedef struct {
43c4762a1bSJed Brown   SolType solType; /* The type of exact solution */
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46d71ae5a4SJacob Faibussowitsch static void f0_u(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[])
47d71ae5a4SJacob Faibussowitsch {
4869b95281SMatthew G. Knepley   PetscInt d;
4969b95281SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += u_x[uOff_x[0] + d * dim + d];
5069b95281SMatthew G. Knepley }
5169b95281SMatthew G. Knepley 
52c4762a1bSJed Brown /* 2D Dirichlet potential example
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   u = x
55c4762a1bSJed Brown   q = <1, 0>
56c4762a1bSJed Brown   f = 0
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   We will need a boundary integral of u over \Gamma.
59c4762a1bSJed Brown */
60d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown   u[0] = x[0];
633ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
64c4762a1bSJed Brown }
65d71ae5a4SJacob Faibussowitsch static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
66d71ae5a4SJacob Faibussowitsch {
67c4762a1bSJed Brown   PetscInt c;
68c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0;
693ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72d71ae5a4SJacob Faibussowitsch static void f0_linear_u(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[])
73d71ae5a4SJacob Faibussowitsch {
7469b95281SMatthew G. Knepley   f0[0] = 0.0;
7569b95281SMatthew G. Knepley   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
7669b95281SMatthew G. Knepley }
77d71ae5a4SJacob Faibussowitsch static void f0_bd_linear_q(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
78d71ae5a4SJacob Faibussowitsch {
7969b95281SMatthew G. Knepley   PetscScalar potential;
8069b95281SMatthew G. Knepley   PetscInt    d;
8169b95281SMatthew G. Knepley 
823ba16761SJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, linear_u(dim, t, x, dim, &potential, NULL));
8369b95281SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[d] = -potential * n[d];
8469b95281SMatthew G. Knepley }
8569b95281SMatthew G. Knepley 
86c4762a1bSJed Brown /* 2D Dirichlet potential example
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   u = x^2 + y^2
89c4762a1bSJed Brown   q = <2x, 2y>
90c4762a1bSJed Brown   f = 4
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   We will need a boundary integral of u over \Gamma.
93c4762a1bSJed Brown */
94d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
95d71ae5a4SJacob Faibussowitsch {
96c4762a1bSJed Brown   PetscInt d;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   u[0] = 0.0;
99c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[0] += x[d] * x[d];
1003ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
101c4762a1bSJed Brown }
102d71ae5a4SJacob Faibussowitsch static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
103d71ae5a4SJacob Faibussowitsch {
104c4762a1bSJed Brown   PetscInt c;
105c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) u[c] = 2.0 * x[c];
1063ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
107c4762a1bSJed Brown }
108c4762a1bSJed Brown 
109d71ae5a4SJacob Faibussowitsch static void f0_quadratic_u(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[])
110d71ae5a4SJacob Faibussowitsch {
11169b95281SMatthew G. Knepley   f0[0] = -4.0;
11269b95281SMatthew G. Knepley   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
11369b95281SMatthew G. Knepley }
114d71ae5a4SJacob Faibussowitsch static void f0_bd_quadratic_q(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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
115d71ae5a4SJacob Faibussowitsch {
11669b95281SMatthew G. Knepley   PetscScalar potential;
11769b95281SMatthew G. Knepley   PetscInt    d;
11869b95281SMatthew G. Knepley 
1193ba16761SJacob Faibussowitsch   PetscCallAbort(PETSC_COMM_SELF, quadratic_u(dim, t, x, dim, &potential, NULL));
12069b95281SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[d] = -potential * n[d];
12169b95281SMatthew G. Knepley }
12269b95281SMatthew G. Knepley 
123c4762a1bSJed Brown /* 2D Dirichlet potential example
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   u = x (1-x) y (1-y)
126c4762a1bSJed Brown   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
127c4762a1bSJed Brown   f = -y (1-y) - x (1-x)
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   u|_\Gamma = 0 so that the boundary integral vanishes
130c4762a1bSJed Brown */
131d71ae5a4SJacob Faibussowitsch static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown   PetscInt d;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   u[0] = 1.0;
136c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[0] *= x[d] * (1.0 - x[d]);
1373ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
138c4762a1bSJed Brown }
139d71ae5a4SJacob Faibussowitsch static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
140d71ae5a4SJacob Faibussowitsch {
141c4762a1bSJed Brown   PetscInt c, d;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
144c4762a1bSJed Brown     u[c] = 1.0;
145c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
146c4762a1bSJed Brown       if (c == d) u[c] *= 1 - 2.0 * x[d];
147c4762a1bSJed Brown       else u[c] *= x[d] * (1.0 - x[d]);
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown   }
1503ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153d71ae5a4SJacob Faibussowitsch static void f0_quartic_u(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[])
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   PetscInt d;
156c4762a1bSJed Brown   f0[0] = 0.0;
15769b95281SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += 2.0 * x[d] * (1.0 - x[d]);
15869b95281SMatthew G. Knepley   f0_u(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, f0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
16169b95281SMatthew G. Knepley /* 2D Dirichlet potential example
16269b95281SMatthew G. Knepley 
16369b95281SMatthew G. Knepley   u = x (1-x) y (1-y)
16469b95281SMatthew G. Knepley   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
16569b95281SMatthew G. Knepley   f = -y (1-y) - x (1-x)
16669b95281SMatthew G. Knepley 
16769b95281SMatthew G. Knepley   du/dn_\Gamma = {(1-2x) y (1-y), x (1-x) (1-2y)} . n produces an essential condition on q
16869b95281SMatthew G. Knepley */
16969b95281SMatthew G. Knepley 
170d71ae5a4SJacob Faibussowitsch static void f0_q(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[])
171d71ae5a4SJacob Faibussowitsch {
172c4762a1bSJed Brown   PetscInt c;
17369b95281SMatthew G. Knepley   for (c = 0; c < dim; ++c) f0[c] = u[uOff[0] + c];
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown /* <\nabla\cdot w, u> = <\nabla w, Iu> */
177d71ae5a4SJacob Faibussowitsch static void f1_q(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 f1[])
178d71ae5a4SJacob Faibussowitsch {
17969b95281SMatthew G. Knepley   PetscInt c;
18069b95281SMatthew G. Knepley   for (c = 0; c < dim; ++c) f1[c * dim + c] = u[uOff[1]];
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
18369b95281SMatthew G. Knepley /* <t, q> */
184d71ae5a4SJacob Faibussowitsch static void g0_qq(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
185d71ae5a4SJacob Faibussowitsch {
186c4762a1bSJed Brown   PetscInt c;
187c4762a1bSJed Brown   for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0;
188c4762a1bSJed Brown }
189c4762a1bSJed Brown 
19069b95281SMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */
191d71ae5a4SJacob Faibussowitsch static void g2_qu(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
192d71ae5a4SJacob Faibussowitsch {
193c4762a1bSJed Brown   PetscInt d;
194c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d * dim + d] = 1.0;
195c4762a1bSJed Brown }
196c4762a1bSJed Brown 
19769b95281SMatthew G. Knepley /* <v, \nabla\cdot q> */
198d71ae5a4SJacob Faibussowitsch static void g1_uq(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, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
199d71ae5a4SJacob Faibussowitsch {
200c4762a1bSJed Brown   PetscInt d;
20169b95281SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
205d71ae5a4SJacob Faibussowitsch {
206c4762a1bSJed Brown   PetscFunctionBeginUser;
207c4762a1bSJed Brown   options->solType = SOL_LINEAR;
208d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");
2099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum)options->solType, (PetscEnum *)&options->solType, NULL));
210d0609cedSBarry Smith   PetscOptionsEnd();
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
215d71ae5a4SJacob Faibussowitsch {
216c4762a1bSJed Brown   PetscFunctionBeginUser;
2179566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2189566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Example Mesh"));
2209566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*dm, user));
2219566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
2229566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224c4762a1bSJed Brown }
225c4762a1bSJed Brown 
226d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
227d71ae5a4SJacob Faibussowitsch {
22845480ffeSMatthew G. Knepley   PetscDS        ds;
22945480ffeSMatthew G. Knepley   DMLabel        label;
23045480ffeSMatthew G. Knepley   PetscWeakForm  wf;
231c4762a1bSJed Brown   const PetscInt id = 1;
23245480ffeSMatthew G. Knepley   PetscInt       bd;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   PetscFunctionBeginUser;
2359566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
2369566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2379566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(ds, 0, f0_q, f1_q));
2389566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_qq, NULL, NULL, NULL));
2399566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_qu, NULL));
2409566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_uq, NULL, NULL));
2419371c9d4SSatish Balay   switch (user->solType) {
242c4762a1bSJed Brown   case SOL_LINEAR:
2439566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 1, f0_linear_u, NULL));
2449566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
2459566063dSJacob Faibussowitsch     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2469566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_linear_q, 0, NULL));
2479566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, linear_q, user));
2489566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 1, linear_u, user));
249c4762a1bSJed Brown     break;
250c4762a1bSJed Brown   case SOL_QUADRATIC:
2519566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 1, f0_quadratic_u, NULL));
2529566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
2539566063dSJacob Faibussowitsch     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2549566063dSJacob Faibussowitsch     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_quadratic_q, 0, NULL));
2559566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user));
2569566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user));
257c4762a1bSJed Brown     break;
258c4762a1bSJed Brown   case SOL_QUARTIC:
2599566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL));
2609566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user));
2619566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user));
262c4762a1bSJed Brown     break;
26369b95281SMatthew G. Knepley   case SOL_QUARTIC_NEUMANN:
2649566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 1, f0_quartic_u, NULL));
2659566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, quartic_q, user));
2669566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 1, quartic_u, user));
2679566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "Flux condition", label, 1, &id, 0, 0, NULL, (void (*)(void))quartic_q, NULL, user, NULL));
26869b95281SMatthew G. Knepley     break;
269d71ae5a4SJacob Faibussowitsch   default:
270d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
271c4762a1bSJed Brown   }
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
276d71ae5a4SJacob Faibussowitsch {
277c4762a1bSJed Brown   DM             cdm = dm;
278c4762a1bSJed Brown   PetscFE        feq, feu;
27969b95281SMatthew G. Knepley   DMPolytopeType ct;
28069b95281SMatthew G. Knepley   PetscBool      simplex;
28169b95281SMatthew G. Knepley   PetscInt       dim, cStart;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   PetscFunctionBeginUser;
2849566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2859566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
28769b95281SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
288c4762a1bSJed Brown   /* Create finite element */
2899566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "field_", -1, &feq));
2909566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)feq, "field"));
2919566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "potential_", -1, &feu));
2929566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)feu, "potential"));
2939566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(feq, feu));
294c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
2959566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)feq));
2969566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)feu));
2979566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2989566063dSJacob Faibussowitsch   PetscCall((*setup)(dm, user));
299c4762a1bSJed Brown   while (cdm) {
3009566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
3019566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
302c4762a1bSJed Brown   }
3039566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&feq));
3049566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&feu));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306c4762a1bSJed Brown }
307c4762a1bSJed Brown 
308d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
309d71ae5a4SJacob Faibussowitsch {
310c4762a1bSJed Brown   DM     dm;   /* Problem specification */
311c4762a1bSJed Brown   SNES   snes; /* Nonlinear solver */
312c4762a1bSJed Brown   Vec    u;    /* Solutions */
313c4762a1bSJed Brown   AppCtx user; /* User-defined work context */
314c4762a1bSJed Brown 
315327415f7SBarry Smith   PetscFunctionBeginUser;
3169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3179566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
318c4762a1bSJed Brown   /* Primal system */
3199566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
3209566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
3219566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
3229566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, SetupPrimalProblem, &user));
3239566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
3249566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
3259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "potential"));
3266493148fSStefano Zampini   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
3279566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
3289566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
3299566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
3309566063dSJacob Faibussowitsch   PetscCall(SNESGetSolution(snes, &u));
3319566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-potential_view"));
332c4762a1bSJed Brown   /* Cleanup */
3339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3349566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
3359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
337b122ec5aSJacob Faibussowitsch   return 0;
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown /*TEST
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   test:
34369b95281SMatthew G. Knepley     suffix: 2d_bdm1_p0
344c4762a1bSJed Brown     requires: triangle
345c4762a1bSJed Brown     args: -sol_type linear \
34669b95281SMatthew G. Knepley           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \
347c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
348c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
349c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
350c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
351c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
352c4762a1bSJed Brown   test:
35307a75e27SStefano Zampini     nsize: 4
35407a75e27SStefano Zampini     suffix: 2d_bdm1_p0_bddc
355*910b42cbSStefano Zampini     requires: triangle !single
35607a75e27SStefano Zampini     args: -sol_type linear \
35707a75e27SStefano Zampini           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 1 \
35807a75e27SStefano Zampini           -dmsnes_check .001 -snes_error_if_not_converged \
35907a75e27SStefano Zampini           -ksp_error_if_not_converged -ksp_type cg \
36007a75e27SStefano Zampini           -petscpartitioner_type simple -dm_mat_type is \
36107a75e27SStefano Zampini           -pc_type bddc -pc_bddc_use_local_mat_graph 0 \
36207a75e27SStefano Zampini           -pc_bddc_benign_trick -pc_bddc_nonetflux -pc_bddc_detect_disconnected -pc_bddc_use_qr_single \
36307a75e27SStefano Zampini           -pc_bddc_coarse_redundant_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_dirichlet_pc_type svd
36407a75e27SStefano Zampini 
36507a75e27SStefano Zampini   test:
36607a75e27SStefano Zampini     nsize: 9
367*910b42cbSStefano Zampini     requires: !single
36807a75e27SStefano Zampini     suffix: 2d_rt1_p0_bddc
36907a75e27SStefano Zampini     args: -sol_type quadratic \
37007a75e27SStefano Zampini           -potential_petscspace_degree 0 \
37107a75e27SStefano Zampini           -potential_petscdualspace_lagrange_use_moments \
37207a75e27SStefano Zampini           -potential_petscdualspace_lagrange_moment_order 2 \
37307a75e27SStefano Zampini           -field_petscfe_default_quadrature_order 1 \
37407a75e27SStefano Zampini           -field_petscspace_degree 1 \
37507a75e27SStefano Zampini           -field_petscspace_type sum \
37607a75e27SStefano Zampini           -field_petscspace_variables 2 \
37707a75e27SStefano Zampini           -field_petscspace_components 2 \
37807a75e27SStefano Zampini           -field_petscspace_sum_spaces 2 \
37907a75e27SStefano Zampini           -field_petscspace_sum_concatenate true \
38007a75e27SStefano Zampini           -field_sumcomp_0_petscspace_variables 2 \
38107a75e27SStefano Zampini           -field_sumcomp_0_petscspace_type tensor \
38207a75e27SStefano Zampini           -field_sumcomp_0_petscspace_tensor_spaces 2 \
38307a75e27SStefano Zampini           -field_sumcomp_0_petscspace_tensor_uniform false \
38407a75e27SStefano Zampini           -field_sumcomp_0_tensorcomp_0_petscspace_degree 1 \
38507a75e27SStefano Zampini           -field_sumcomp_0_tensorcomp_1_petscspace_degree 0 \
38607a75e27SStefano Zampini           -field_sumcomp_1_petscspace_variables 2 \
38707a75e27SStefano Zampini           -field_sumcomp_1_petscspace_type tensor \
38807a75e27SStefano Zampini           -field_sumcomp_1_petscspace_tensor_spaces 2 \
38907a75e27SStefano Zampini           -field_sumcomp_1_petscspace_tensor_uniform false \
39007a75e27SStefano Zampini           -field_sumcomp_1_tensorcomp_0_petscspace_degree 0 \
39107a75e27SStefano Zampini           -field_sumcomp_1_tensorcomp_1_petscspace_degree 1 \
39207a75e27SStefano Zampini           -field_petscdualspace_form_degree -1 \
39307a75e27SStefano Zampini           -field_petscdualspace_order 1 \
39407a75e27SStefano Zampini           -field_petscdualspace_lagrange_trimmed true \
39507a75e27SStefano Zampini           -dm_plex_box_faces 3,3 dm_refine 1 -dm_plex_simplex 0 \
39607a75e27SStefano Zampini           -dmsnes_check .001 -snes_error_if_not_converged \
39707a75e27SStefano Zampini           -ksp_error_if_not_converged -ksp_type cg \
39807a75e27SStefano Zampini           -petscpartitioner_type simple -dm_mat_type is \
39907a75e27SStefano Zampini           -pc_type bddc -pc_bddc_use_local_mat_graph 0 \
40007a75e27SStefano Zampini           -pc_bddc_benign_trick -pc_bddc_nonetflux -pc_bddc_detect_disconnected -pc_bddc_use_qr_single \
40107a75e27SStefano Zampini           -pc_bddc_coarse_redundant_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_dirichlet_pc_type svd
40207a75e27SStefano Zampini 
40307a75e27SStefano Zampini   test:
40469b95281SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.0, 1.0]
40569b95281SMatthew G. Knepley     # Using -sol_type quadratic -dm_refine 2 -convest_num_refine 3 we get L_2 convergence rate: [2.9, 1.0]
40669b95281SMatthew G. Knepley     suffix: 2d_bdm1_p0_conv
407c4762a1bSJed Brown     requires: triangle
408c4762a1bSJed Brown     args: -sol_type quartic \
409c4762a1bSJed Brown           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
410c4762a1bSJed Brown           -snes_error_if_not_converged \
411c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
412c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
413c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
414c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
415c4762a1bSJed Brown 
416c4762a1bSJed Brown   test:
41769b95281SMatthew G. Knepley     # HDF5 output: -dm_view hdf5:${PETSC_DIR}/sol.h5 -potential_view hdf5:${PETSC_DIR}/sol.h5::append -exact_vec_view hdf5:${PETSC_DIR}/sol.h5::append
41869b95281SMatthew G. Knepley     # VTK output: -potential_view vtk: -exact_vec_view vtk:
41969b95281SMatthew G. Knepley     # VTU output: -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu
42069b95281SMatthew G. Knepley     suffix: 2d_q2_p0
421c4762a1bSJed Brown     requires: triangle
42230602db0SMatthew G. Knepley     args: -sol_type linear -dm_plex_simplex 0 \
42369b95281SMatthew G. Knepley           -field_petscspace_degree 2 -dm_refine 0 \
424c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
425c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
426c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
427c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
428c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
429c4762a1bSJed Brown   test:
43069b95281SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [0.0057, 1.0]
43169b95281SMatthew G. Knepley     suffix: 2d_q2_p0_conv
432c4762a1bSJed Brown     requires: triangle
43330602db0SMatthew G. Knepley     args: -sol_type linear -dm_plex_simplex 0 \
434c4762a1bSJed Brown           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
43569b95281SMatthew G. Knepley           -snes_error_if_not_converged \
436c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
437c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
438c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
439c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
44069b95281SMatthew G. Knepley   test:
44169b95281SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 3 we get L_2 convergence rate: [-0.014, 0.0066]
44269b95281SMatthew G. Knepley     suffix: 2d_q2_p0_neumann_conv
44369b95281SMatthew G. Knepley     requires: triangle
44430602db0SMatthew G. Knepley     args: -sol_type quartic_neumann -dm_plex_simplex 0 \
44569b95281SMatthew G. Knepley           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
44669b95281SMatthew G. Knepley           -snes_error_if_not_converged \
44769b95281SMatthew G. Knepley           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
44869b95281SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
44969b95281SMatthew G. Knepley             -fieldsplit_field_pc_type lu \
45069b95281SMatthew G. Knepley             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type svd
45169b95281SMatthew G. Knepley 
452c4762a1bSJed Brown TEST*/
453c4762a1bSJed Brown 
454c4762a1bSJed Brown /* These tests will be active once tensor P^- is working
455c4762a1bSJed Brown 
456c4762a1bSJed Brown   test:
457c4762a1bSJed Brown     suffix: 2d_bdmq1_p0_0
458c4762a1bSJed Brown     requires: triangle
45930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -sol_type linear \
460c4762a1bSJed Brown           -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
461c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
462c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
463c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
464c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
465c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
466c4762a1bSJed Brown   test:
467c4762a1bSJed Brown     suffix: 2d_bdmq1_p0_2
468c4762a1bSJed Brown     requires: triangle
46930602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -sol_type quartic \
470c4762a1bSJed Brown           -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
471c4762a1bSJed Brown           -dmsnes_check .001 -snes_error_if_not_converged \
472c4762a1bSJed Brown           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
473c4762a1bSJed Brown           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
474c4762a1bSJed Brown             -fieldsplit_field_pc_type lu \
475c4762a1bSJed Brown             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
476c4762a1bSJed Brown 
477c4762a1bSJed Brown */
478