xref: /petsc/src/snes/tutorials/ex27.c (revision 47bfdf3fa30783d319821566bb96996bc8fc6067)
1*47bfdf3fSMatthew G. Knepley static char help[] = "Poisson Problem in 2d and 3d with simplicial finite elements in both primal and mixed form.\n\
2*47bfdf3fSMatthew G. Knepley We solve the Poisson problem in a rectangular\n\
3*47bfdf3fSMatthew G. Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4*47bfdf3fSMatthew G. Knepley This example solves mixed form equation to get the flux field to calculate flux norm. We then use that for adaptive mesh refinement. \n\n\n";
5*47bfdf3fSMatthew G. Knepley 
6*47bfdf3fSMatthew G. Knepley /*
7*47bfdf3fSMatthew G. Knepley The primal (or original) Poisson problem, in the strong form, is given by,
8*47bfdf3fSMatthew G. Knepley 
9*47bfdf3fSMatthew G. Knepley \begin{align}
10*47bfdf3fSMatthew G. Knepley   - \nabla \cdot ( \nabla u ) = f
11*47bfdf3fSMatthew G. Knepley \end{align}
12*47bfdf3fSMatthew G. Knepley where $u$ is potential.
13*47bfdf3fSMatthew G. Knepley 
14*47bfdf3fSMatthew G. Knepley The weak form of this equation is
15*47bfdf3fSMatthew G. Knepley 
16*47bfdf3fSMatthew G. Knepley \begin{align}
17*47bfdf3fSMatthew G. Knepley   < \nabla v, \nabla u > - < v, \nabla u \cdot \hat{n} >_\Gamma - < v, f > = 0
18*47bfdf3fSMatthew G. Knepley \end{align}
19*47bfdf3fSMatthew G. Knepley 
20*47bfdf3fSMatthew G. Knepley The mixed Poisson problem, in the strong form, is given by,
21*47bfdf3fSMatthew G. Knepley 
22*47bfdf3fSMatthew G. Knepley \begin{align}
23*47bfdf3fSMatthew G. Knepley   q - \nabla u &= 0 \\
24*47bfdf3fSMatthew G. Knepley   - \nabla \cdot q &= f
25*47bfdf3fSMatthew G. Knepley \end{align}
26*47bfdf3fSMatthew G. Knepley where $u$ is the potential and $q$ is the flux.
27*47bfdf3fSMatthew G. Knepley 
28*47bfdf3fSMatthew G. Knepley The weak form of this equation is
29*47bfdf3fSMatthew G. Knepley 
30*47bfdf3fSMatthew G. Knepley \begin{align}
31*47bfdf3fSMatthew G. Knepley   < t, q > + < \nabla \cdot t, u > - < t \cdot \hat{n}, u >_\Gamma &= 0 \\
32*47bfdf3fSMatthew G. Knepley   <v, \nabla \cdot q> - < v, f > &= 0
33*47bfdf3fSMatthew G. Knepley \end{align}
34*47bfdf3fSMatthew G. Knepley 
35*47bfdf3fSMatthew G. Knepley We solve both primal and mixed problem and calculate the error in the flux norm, namely || e || = || q^m_h - \nabla u^p_h ||. Here superscript 'm' represents field from mixed form and 'p' represents field from the primal form.
36*47bfdf3fSMatthew G. Knepley 
37*47bfdf3fSMatthew G. Knepley The following boundary conditions are prescribed.
38*47bfdf3fSMatthew G. Knepley 
39*47bfdf3fSMatthew G. Knepley Primal problem:
40*47bfdf3fSMatthew G. Knepley \begin{align}
41*47bfdf3fSMatthew G. Knepley   u = u_0                    on \Gamma_D
42*47bfdf3fSMatthew G. Knepley   \nabla u \cdot \hat{n} = g on \Gamma_N
43*47bfdf3fSMatthew G. Knepley \end{align}
44*47bfdf3fSMatthew G. Knepley 
45*47bfdf3fSMatthew G. Knepley Mixed problem:
46*47bfdf3fSMatthew G. Knepley \begin{align}
47*47bfdf3fSMatthew G. Knepley   u = u_0             on \Gamma_D
48*47bfdf3fSMatthew G. Knepley   q \cdot \hat{n} = g on \Gamma_N
49*47bfdf3fSMatthew G. Knepley \end{align}
50*47bfdf3fSMatthew G. Knepley         __________\Gamma_D_____________
51*47bfdf3fSMatthew G. Knepley         |                              |
52*47bfdf3fSMatthew G. Knepley         |                              |
53*47bfdf3fSMatthew G. Knepley         |                              |
54*47bfdf3fSMatthew G. Knepley \Gamma_N                               \Gamma_N
55*47bfdf3fSMatthew G. Knepley         |                              |
56*47bfdf3fSMatthew G. Knepley         |                              |
57*47bfdf3fSMatthew G. Knepley         |                              |
58*47bfdf3fSMatthew G. Knepley         |_________\Gamma_D_____________|
59*47bfdf3fSMatthew G. Knepley 
60*47bfdf3fSMatthew G. Knepley To visualize the automated adaptation
61*47bfdf3fSMatthew G. Knepley 
62*47bfdf3fSMatthew G. Knepley   -dm_adapt_pre_view draw -dm_adapt_view draw -draw_pause -1 -geometry 0,0,1024,1024
63*47bfdf3fSMatthew G. Knepley 
64*47bfdf3fSMatthew G. Knepley and to compare with a naice gradient estimator use
65*47bfdf3fSMatthew G. Knepley 
66*47bfdf3fSMatthew G. Knepley   -adaptor_type gradient
67*47bfdf3fSMatthew G. Knepley 
68*47bfdf3fSMatthew G. Knepley To see a sequence of adaptations
69*47bfdf3fSMatthew G. Knepley 
70*47bfdf3fSMatthew G. Knepley   -snes_adapt_sequence 8 -adaptor_monitor_error draw::draw_lg
71*47bfdf3fSMatthew G. Knepley   -dm_adapt_pre_view draw -dm_adapt_iter_view draw -dm_adapt_view draw -draw_pause 1 -geometry 0,0,1024,1024
72*47bfdf3fSMatthew G. Knepley 
73*47bfdf3fSMatthew G. Knepley To get a better view of the by-hand process, use
74*47bfdf3fSMatthew G. Knepley 
75*47bfdf3fSMatthew G. Knepley   -dm_view hdf5:${PWD}/mesh.h5
76*47bfdf3fSMatthew G. Knepley   -primal_sol_vec_view hdf5:${PWD}/mesh.h5::append
77*47bfdf3fSMatthew G. Knepley   -flux_error_vec_view hdf5:${PWD}/mesh.h5::append
78*47bfdf3fSMatthew G. Knepley   -exact_error_vec_view hdf5:${PWD}/mesh.h5::append
79*47bfdf3fSMatthew G. Knepley   -refine_vec_view hdf5:${PWD}/mesh.h5::append
80*47bfdf3fSMatthew G. Knepley   -adapt_dm_view draw -draw_pause -1
81*47bfdf3fSMatthew G. Knepley 
82*47bfdf3fSMatthew G. Knepley This is also possible with the automated path
83*47bfdf3fSMatthew G. Knepley 
84*47bfdf3fSMatthew G. Knepley   -dm_view hdf5:${PWD}/mesh.h5
85*47bfdf3fSMatthew G. Knepley   -adapt_primal_sol_vec_view hdf5:${PWD}/mesh.h5::append
86*47bfdf3fSMatthew G. Knepley   -adapt_error_vec_view hdf5:${PWD}/mesh.h5::append
87*47bfdf3fSMatthew G. Knepley   -adapt_vec_view hdf5:${PWD}/mesh.h5::append
88*47bfdf3fSMatthew G. Knepley */
89*47bfdf3fSMatthew G. Knepley 
90*47bfdf3fSMatthew G. Knepley #include <petsc/private/petscfeimpl.h>
91*47bfdf3fSMatthew G. Knepley #include <petscdmplex.h>
92*47bfdf3fSMatthew G. Knepley #include <petscsnes.h>
93*47bfdf3fSMatthew G. Knepley #include <petscdmadaptor.h>
94*47bfdf3fSMatthew G. Knepley #include <petscds.h>
95*47bfdf3fSMatthew G. Knepley #include <petscviewerhdf5.h>
96*47bfdf3fSMatthew G. Knepley #include <petscbag.h>
97*47bfdf3fSMatthew G. Knepley 
98*47bfdf3fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SetupMixed(DMAdaptor, DM);
99*47bfdf3fSMatthew G. Knepley 
100*47bfdf3fSMatthew G. Knepley typedef enum {
101*47bfdf3fSMatthew G. Knepley   SOL_QUADRATIC,
102*47bfdf3fSMatthew G. Knepley   SOL_TRIG,
103*47bfdf3fSMatthew G. Knepley   SOL_SENSOR,
104*47bfdf3fSMatthew G. Knepley   SOL_UNKNOWN,
105*47bfdf3fSMatthew G. Knepley   NUM_SOL_TYPE
106*47bfdf3fSMatthew G. Knepley } SolType;
107*47bfdf3fSMatthew G. Knepley const char *SolTypeNames[NUM_SOL_TYPE + 4] = {"quadratic", "trig", "sensor", "unknown", "SolType", "SOL_", NULL};
108*47bfdf3fSMatthew G. Knepley 
109*47bfdf3fSMatthew G. Knepley typedef struct {
110*47bfdf3fSMatthew G. Knepley   PetscBag  param;
111*47bfdf3fSMatthew G. Knepley   SolType   solType;
112*47bfdf3fSMatthew G. Knepley   PetscBool byHand;
113*47bfdf3fSMatthew G. Knepley } AppCtx;
114*47bfdf3fSMatthew G. Knepley 
115*47bfdf3fSMatthew G. Knepley typedef struct {
116*47bfdf3fSMatthew G. Knepley   PetscReal k;
117*47bfdf3fSMatthew G. Knepley } Parameter;
118*47bfdf3fSMatthew G. Knepley 
119*47bfdf3fSMatthew G. Knepley /* Exact solution: u = x^2 + y^2 */
120*47bfdf3fSMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
121*47bfdf3fSMatthew G. Knepley {
122*47bfdf3fSMatthew G. Knepley   PetscInt d;
123*47bfdf3fSMatthew G. Knepley 
124*47bfdf3fSMatthew G. Knepley   u[0] = 0.0;
125*47bfdf3fSMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d] * x[d];
126*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
127*47bfdf3fSMatthew G. Knepley }
128*47bfdf3fSMatthew G. Knepley /* Exact solution: q = (2x, 2y) */
129*47bfdf3fSMatthew G. Knepley static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
130*47bfdf3fSMatthew G. Knepley {
131*47bfdf3fSMatthew G. Knepley   PetscInt c;
132*47bfdf3fSMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 2.0 * x[c];
133*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
134*47bfdf3fSMatthew G. Knepley }
135*47bfdf3fSMatthew G. Knepley 
136*47bfdf3fSMatthew G. Knepley /* Exact solution: u = sin( n \pi x ) * sin( n \pi y ) */
137*47bfdf3fSMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
138*47bfdf3fSMatthew G. Knepley {
139*47bfdf3fSMatthew G. Knepley   const PetscReal n = 5.5;
140*47bfdf3fSMatthew G. Knepley 
141*47bfdf3fSMatthew G. Knepley   u[0] = 1.0;
142*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) u[0] *= PetscSinReal(n * PETSC_PI * x[d]);
143*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
144*47bfdf3fSMatthew G. Knepley }
145*47bfdf3fSMatthew G. Knepley static PetscErrorCode trig_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
146*47bfdf3fSMatthew G. Knepley {
147*47bfdf3fSMatthew G. Knepley   const PetscReal n = 5.5;
148*47bfdf3fSMatthew G. Knepley 
149*47bfdf3fSMatthew G. Knepley   for (PetscInt c = 0; c < Nc; c++) u[c] = n * PETSC_PI * PetscCosReal(n * PETSC_PI * x[c]) * PetscSinReal(n * PETSC_PI * x[Nc - c - 1]);
150*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
151*47bfdf3fSMatthew G. Knepley }
152*47bfdf3fSMatthew G. Knepley 
153*47bfdf3fSMatthew G. Knepley /*
154*47bfdf3fSMatthew G. Knepley Classic hyperbolic sensor function for testing multi-scale anisotropic mesh adaptation:
155*47bfdf3fSMatthew G. Knepley 
156*47bfdf3fSMatthew G. Knepley   f:[-1, 1]x[-1, 1] \to R,
157*47bfdf3fSMatthew G. Knepley     f(x, y) = sin(50xy)/100 if |xy| > 2\pi/50 else sin(50xy)
158*47bfdf3fSMatthew G. Knepley 
159*47bfdf3fSMatthew G. Knepley (mapped to have domain [0,1] x [0,1] in this case).
160*47bfdf3fSMatthew G. Knepley */
161*47bfdf3fSMatthew G. Knepley static PetscErrorCode sensor_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
162*47bfdf3fSMatthew G. Knepley {
163*47bfdf3fSMatthew G. Knepley   const PetscReal xref = 2. * x[0] - 1.;
164*47bfdf3fSMatthew G. Knepley   const PetscReal yref = 2. * x[1] - 1.;
165*47bfdf3fSMatthew G. Knepley   const PetscReal xy   = xref * yref;
166*47bfdf3fSMatthew G. Knepley 
167*47bfdf3fSMatthew G. Knepley   u[0] = PetscSinReal(50. * xy);
168*47bfdf3fSMatthew G. Knepley   if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) u[0] *= 0.01;
169*47bfdf3fSMatthew G. Knepley 
170*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
171*47bfdf3fSMatthew G. Knepley }
172*47bfdf3fSMatthew G. Knepley 
173*47bfdf3fSMatthew G. Knepley /* Flux is (cos(50xy) * 50y/100, cos(50xy) * 50x/100) if |xy| > 2\pi/50 else (cos(50xy) * 50y, cos(50xy) * 50x) */
174*47bfdf3fSMatthew G. Knepley static PetscErrorCode sensor_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
175*47bfdf3fSMatthew G. Knepley {
176*47bfdf3fSMatthew G. Knepley   const PetscReal xref = 2. * x[0] - 1.;
177*47bfdf3fSMatthew G. Knepley   const PetscReal yref = 2. * x[1] - 1.;
178*47bfdf3fSMatthew G. Knepley   const PetscReal xy   = xref * yref;
179*47bfdf3fSMatthew G. Knepley 
180*47bfdf3fSMatthew G. Knepley   u[0] = 50. * yref * PetscCosReal(50. * xy) * 2.0;
181*47bfdf3fSMatthew G. Knepley   u[1] = 50. * xref * PetscCosReal(50. * xy) * 2.0;
182*47bfdf3fSMatthew G. Knepley   if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) {
183*47bfdf3fSMatthew G. Knepley     u[0] *= 0.01;
184*47bfdf3fSMatthew G. Knepley     u[1] *= 0.01;
185*47bfdf3fSMatthew G. Knepley   }
186*47bfdf3fSMatthew G. Knepley   return PETSC_SUCCESS;
187*47bfdf3fSMatthew G. Knepley }
188*47bfdf3fSMatthew G. Knepley 
189*47bfdf3fSMatthew G. Knepley /* We set up residuals and Jacobians for the primal problem. */
190*47bfdf3fSMatthew G. Knepley static void f0_quadratic_primal(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[])
191*47bfdf3fSMatthew G. Knepley {
192*47bfdf3fSMatthew G. Knepley   f0[0] = 4.0;
193*47bfdf3fSMatthew G. Knepley }
194*47bfdf3fSMatthew G. Knepley 
195*47bfdf3fSMatthew G. Knepley static void f0_trig_primal(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[])
196*47bfdf3fSMatthew G. Knepley {
197*47bfdf3fSMatthew G. Knepley   const PetscReal n = 5.5;
198*47bfdf3fSMatthew G. Knepley 
199*47bfdf3fSMatthew G. Knepley   f0[0] = -2.0 * PetscSqr(n * PETSC_PI) * PetscSinReal(n * PETSC_PI * x[0]) * PetscSinReal(n * PETSC_PI * x[1]);
200*47bfdf3fSMatthew G. Knepley }
201*47bfdf3fSMatthew G. Knepley 
202*47bfdf3fSMatthew G. Knepley static void f0_sensor_primal(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[])
203*47bfdf3fSMatthew G. Knepley {
204*47bfdf3fSMatthew G. Knepley   const PetscReal xref = 2. * x[0] - 1.;
205*47bfdf3fSMatthew G. Knepley   const PetscReal yref = 2. * x[1] - 1.;
206*47bfdf3fSMatthew G. Knepley   const PetscReal xy   = xref * yref;
207*47bfdf3fSMatthew G. Knepley 
208*47bfdf3fSMatthew G. Knepley   f0[0] = -2500.0 * PetscSinReal(50. * xy) * (xref * xref + yref * yref) * 4.0;
209*47bfdf3fSMatthew G. Knepley   if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) f0[0] *= 0.01;
210*47bfdf3fSMatthew G. Knepley }
211*47bfdf3fSMatthew G. Knepley 
212*47bfdf3fSMatthew G. Knepley static void f1_primal(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[])
213*47bfdf3fSMatthew G. Knepley {
214*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
215*47bfdf3fSMatthew G. Knepley 
216*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f1[d] = k * u_x[uOff_x[0] + d];
217*47bfdf3fSMatthew G. Knepley }
218*47bfdf3fSMatthew G. Knepley 
219*47bfdf3fSMatthew G. Knepley static void f0_quadratic_bd_primal(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[])
220*47bfdf3fSMatthew G. Knepley {
221*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
222*47bfdf3fSMatthew G. Knepley   PetscScalar     flux;
223*47bfdf3fSMatthew G. Knepley 
224*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, quadratic_q(dim, t, x, dim, &flux, NULL));
225*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux * n[d];
226*47bfdf3fSMatthew G. Knepley }
227*47bfdf3fSMatthew G. Knepley 
228*47bfdf3fSMatthew G. Knepley static void f0_trig_bd_primal(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[])
229*47bfdf3fSMatthew G. Knepley {
230*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
231*47bfdf3fSMatthew G. Knepley   PetscScalar     flux;
232*47bfdf3fSMatthew G. Knepley 
233*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, trig_q(dim, t, x, dim, &flux, NULL));
234*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux * n[d];
235*47bfdf3fSMatthew G. Knepley }
236*47bfdf3fSMatthew G. Knepley 
237*47bfdf3fSMatthew G. Knepley static void f0_sensor_bd_primal(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[])
238*47bfdf3fSMatthew G. Knepley {
239*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
240*47bfdf3fSMatthew G. Knepley   PetscScalar     flux;
241*47bfdf3fSMatthew G. Knepley 
242*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, sensor_q(dim, t, x, dim, &flux, NULL));
243*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux * n[d];
244*47bfdf3fSMatthew G. Knepley }
245*47bfdf3fSMatthew G. Knepley 
246*47bfdf3fSMatthew G. Knepley static void g3_primal_uu(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 g3[])
247*47bfdf3fSMatthew G. Knepley {
248*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
249*47bfdf3fSMatthew G. Knepley 
250*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = k;
251*47bfdf3fSMatthew G. Knepley }
252*47bfdf3fSMatthew G. Knepley 
253*47bfdf3fSMatthew G. Knepley /* Now we set up the residuals and Jacobians mixed problem. */
254*47bfdf3fSMatthew G. Knepley static void f0_mixed_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[])
255*47bfdf3fSMatthew G. Knepley {
256*47bfdf3fSMatthew G. Knepley   f0[0] = 4.0;
257*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d];
258*47bfdf3fSMatthew G. Knepley }
259*47bfdf3fSMatthew G. Knepley static void f0_mixed_trig_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[])
260*47bfdf3fSMatthew G. Knepley {
261*47bfdf3fSMatthew G. Knepley   PetscReal n = 5.5;
262*47bfdf3fSMatthew G. Knepley 
263*47bfdf3fSMatthew G. Knepley   f0[0] = -2.0 * PetscSqr(n * PETSC_PI) * PetscSinReal(n * PETSC_PI * x[0]) * PetscSinReal(n * PETSC_PI * x[1]);
264*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d];
265*47bfdf3fSMatthew G. Knepley }
266*47bfdf3fSMatthew G. Knepley static void f0_mixed_sensor_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[])
267*47bfdf3fSMatthew G. Knepley {
268*47bfdf3fSMatthew G. Knepley   const PetscReal xref = 2. * x[0] - 1.;
269*47bfdf3fSMatthew G. Knepley   const PetscReal yref = 2. * x[1] - 1.;
270*47bfdf3fSMatthew G. Knepley   const PetscReal xy   = xref * yref;
271*47bfdf3fSMatthew G. Knepley 
272*47bfdf3fSMatthew G. Knepley   f0[0] = -2500.0 * PetscSinReal(50. * xy) * (xref * xref + yref * yref) * 4.0;
273*47bfdf3fSMatthew G. Knepley   if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) f0[0] *= 0.01;
274*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d];
275*47bfdf3fSMatthew G. Knepley }
276*47bfdf3fSMatthew G. Knepley 
277*47bfdf3fSMatthew G. Knepley static void f0_mixed_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[])
278*47bfdf3fSMatthew G. Knepley {
279*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; d++) f0[d] = u[uOff[0] + d];
280*47bfdf3fSMatthew G. Knepley }
281*47bfdf3fSMatthew G. Knepley 
282*47bfdf3fSMatthew G. Knepley static void f1_mixed_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[])
283*47bfdf3fSMatthew G. Knepley {
284*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
285*47bfdf3fSMatthew G. Knepley 
286*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; d++) f1[d * dim + d] = k * u[uOff[1]];
287*47bfdf3fSMatthew G. Knepley }
288*47bfdf3fSMatthew G. Knepley 
289*47bfdf3fSMatthew G. Knepley static void f0_quadratic_bd_mixed_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[])
290*47bfdf3fSMatthew G. Knepley {
291*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
292*47bfdf3fSMatthew G. Knepley   PetscScalar     potential;
293*47bfdf3fSMatthew G. Knepley 
294*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, quadratic_u(dim, t, x, dim, &potential, NULL));
295*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * potential * n[d];
296*47bfdf3fSMatthew G. Knepley }
297*47bfdf3fSMatthew G. Knepley 
298*47bfdf3fSMatthew G. Knepley static void f0_trig_bd_mixed_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[])
299*47bfdf3fSMatthew G. Knepley {
300*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
301*47bfdf3fSMatthew G. Knepley   PetscScalar     potential;
302*47bfdf3fSMatthew G. Knepley 
303*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, trig_u(dim, t, x, dim, &potential, NULL));
304*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d * dim + d] = -k * potential * n[d];
305*47bfdf3fSMatthew G. Knepley }
306*47bfdf3fSMatthew G. Knepley 
307*47bfdf3fSMatthew G. Knepley static void f0_sensor_bd_mixed_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[])
308*47bfdf3fSMatthew G. Knepley {
309*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
310*47bfdf3fSMatthew G. Knepley   PetscScalar     potential;
311*47bfdf3fSMatthew G. Knepley 
312*47bfdf3fSMatthew G. Knepley   PetscCallAbort(PETSC_COMM_SELF, sensor_u(dim, t, x, dim, &potential, NULL));
313*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; ++d) f0[d * dim + d] = -k * potential * n[d];
314*47bfdf3fSMatthew G. Knepley }
315*47bfdf3fSMatthew G. Knepley 
316*47bfdf3fSMatthew G. Knepley /* <v, \nabla\cdot q> */
317*47bfdf3fSMatthew G. Knepley static void g1_mixed_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[])
318*47bfdf3fSMatthew G. Knepley {
319*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; d++) g1[d * dim + d] = -1.0;
320*47bfdf3fSMatthew G. Knepley }
321*47bfdf3fSMatthew G. Knepley 
322*47bfdf3fSMatthew G. Knepley /* < t, q> */
323*47bfdf3fSMatthew G. Knepley static void g0_mixed_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[])
324*47bfdf3fSMatthew G. Knepley {
325*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; d++) g0[d * dim + d] = 1.0;
326*47bfdf3fSMatthew G. Knepley }
327*47bfdf3fSMatthew G. Knepley 
328*47bfdf3fSMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */
329*47bfdf3fSMatthew G. Knepley static void g2_mixed_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[])
330*47bfdf3fSMatthew G. Knepley {
331*47bfdf3fSMatthew G. Knepley   const PetscReal k = PetscRealPart(constants[0]);
332*47bfdf3fSMatthew G. Knepley 
333*47bfdf3fSMatthew G. Knepley   for (PetscInt d = 0; d < dim; d++) g2[d * dim + d] = k;
334*47bfdf3fSMatthew G. Knepley }
335*47bfdf3fSMatthew G. Knepley 
336*47bfdf3fSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *user)
337*47bfdf3fSMatthew G. Knepley {
338*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
339*47bfdf3fSMatthew G. Knepley   PetscOptionsBegin(comm, "", "Flux norm error in primal poisson problem Options", "DMPLEX");
340*47bfdf3fSMatthew G. Knepley   user->byHand  = PETSC_TRUE;
341*47bfdf3fSMatthew G. Knepley   user->solType = SOL_QUADRATIC;
342*47bfdf3fSMatthew G. Knepley 
343*47bfdf3fSMatthew G. Knepley   PetscCall(PetscOptionsGetBool(NULL, NULL, "-by_hand", &user->byHand, NULL));
344*47bfdf3fSMatthew G. Knepley   PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex27.c", SolTypeNames, (PetscEnum)user->solType, (PetscEnum *)&user->solType, NULL));
345*47bfdf3fSMatthew G. Knepley   PetscOptionsEnd();
346*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
347*47bfdf3fSMatthew G. Knepley }
348*47bfdf3fSMatthew G. Knepley 
349*47bfdf3fSMatthew G. Knepley static PetscErrorCode SetupParameters(PetscBag bag, AppCtx *user)
350*47bfdf3fSMatthew G. Knepley {
351*47bfdf3fSMatthew G. Knepley   Parameter *param;
352*47bfdf3fSMatthew G. Knepley 
353*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
354*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagGetData(bag, (void **)&param));
355*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagSetName(bag, "par", "Poisson parameters"));
356*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagRegisterReal(bag, &param->k, 1.0, "k", "Thermal conductivity"));
357*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagSetFromOptions(bag));
358*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
359*47bfdf3fSMatthew G. Knepley }
360*47bfdf3fSMatthew G. Knepley 
361*47bfdf3fSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
362*47bfdf3fSMatthew G. Knepley {
363*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
364*47bfdf3fSMatthew G. Knepley   PetscCall(DMCreate(comm, dm));
365*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetType(*dm, DMPLEX));
366*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetFromOptions(*dm));
367*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetApplicationContext(*dm, &user));
368*47bfdf3fSMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
369*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
370*47bfdf3fSMatthew G. Knepley }
371*47bfdf3fSMatthew G. Knepley 
372*47bfdf3fSMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
373*47bfdf3fSMatthew G. Knepley {
374*47bfdf3fSMatthew G. Knepley   PetscDS       ds;
375*47bfdf3fSMatthew G. Knepley   DMLabel       label;
376*47bfdf3fSMatthew G. Knepley   PetscInt      id, bd;
377*47bfdf3fSMatthew G. Knepley   Parameter    *param;
378*47bfdf3fSMatthew G. Knepley   PetscWeakForm wf;
379*47bfdf3fSMatthew G. Knepley 
380*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
381*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
382*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "marker", &label));
383*47bfdf3fSMatthew G. Knepley 
384*47bfdf3fSMatthew G. Knepley   PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_primal_uu));
385*47bfdf3fSMatthew G. Knepley 
386*47bfdf3fSMatthew G. Knepley   switch (user->solType) {
387*47bfdf3fSMatthew G. Knepley   case SOL_QUADRATIC:
388*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_primal, f1_primal));
389*47bfdf3fSMatthew G. Knepley 
390*47bfdf3fSMatthew G. Knepley     id = 1;
391*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u, NULL, user, NULL));
392*47bfdf3fSMatthew G. Knepley 
393*47bfdf3fSMatthew G. Knepley     id = 2;
394*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
395*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
396*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_primal, 0, NULL));
397*47bfdf3fSMatthew G. Knepley 
398*47bfdf3fSMatthew G. Knepley     id = 3;
399*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u, NULL, user, NULL));
400*47bfdf3fSMatthew G. Knepley 
401*47bfdf3fSMatthew G. Knepley     id = 4;
402*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
403*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
404*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_primal, 0, NULL));
405*47bfdf3fSMatthew G. Knepley 
406*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user));
407*47bfdf3fSMatthew G. Knepley     break;
408*47bfdf3fSMatthew G. Knepley   case SOL_TRIG:
409*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_primal, f1_primal));
410*47bfdf3fSMatthew G. Knepley 
411*47bfdf3fSMatthew G. Knepley     id = 1;
412*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
413*47bfdf3fSMatthew G. Knepley 
414*47bfdf3fSMatthew G. Knepley     id = 2;
415*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
416*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
417*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_primal, 0, NULL));
418*47bfdf3fSMatthew G. Knepley 
419*47bfdf3fSMatthew G. Knepley     id = 3;
420*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL));
421*47bfdf3fSMatthew G. Knepley 
422*47bfdf3fSMatthew G. Knepley     id = 4;
423*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
424*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
425*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_primal, 0, NULL));
426*47bfdf3fSMatthew G. Knepley 
427*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user));
428*47bfdf3fSMatthew G. Knepley     break;
429*47bfdf3fSMatthew G. Knepley   case SOL_SENSOR:
430*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 0, f0_sensor_primal, f1_primal));
431*47bfdf3fSMatthew G. Knepley 
432*47bfdf3fSMatthew G. Knepley     id = 1;
433*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))sensor_u, NULL, user, NULL));
434*47bfdf3fSMatthew G. Knepley 
435*47bfdf3fSMatthew G. Knepley     id = 2;
436*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
437*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
438*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_primal, 0, NULL));
439*47bfdf3fSMatthew G. Knepley 
440*47bfdf3fSMatthew G. Knepley     id = 3;
441*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))sensor_u, NULL, user, NULL));
442*47bfdf3fSMatthew G. Knepley 
443*47bfdf3fSMatthew G. Knepley     id = 4;
444*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
445*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
446*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_primal, 0, NULL));
447*47bfdf3fSMatthew G. Knepley 
448*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, sensor_u, user));
449*47bfdf3fSMatthew G. Knepley     break;
450*47bfdf3fSMatthew G. Knepley   default:
451*47bfdf3fSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
452*47bfdf3fSMatthew G. Knepley   }
453*47bfdf3fSMatthew G. Knepley 
454*47bfdf3fSMatthew G. Knepley   /* Setup constants */
455*47bfdf3fSMatthew G. Knepley   {
456*47bfdf3fSMatthew G. Knepley     PetscCall(PetscBagGetData(user->param, (void **)&param));
457*47bfdf3fSMatthew G. Knepley     PetscScalar constants[1];
458*47bfdf3fSMatthew G. Knepley 
459*47bfdf3fSMatthew G. Knepley     constants[0] = param->k;
460*47bfdf3fSMatthew G. Knepley 
461*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetConstants(ds, 1, constants));
462*47bfdf3fSMatthew G. Knepley   }
463*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
464*47bfdf3fSMatthew G. Knepley }
465*47bfdf3fSMatthew G. Knepley 
466*47bfdf3fSMatthew G. Knepley static PetscErrorCode SetupPrimalDiscretization(DM dm, AppCtx *user)
467*47bfdf3fSMatthew G. Knepley {
468*47bfdf3fSMatthew G. Knepley   DM             cdm = dm;
469*47bfdf3fSMatthew G. Knepley   PetscFE        fe[1];
470*47bfdf3fSMatthew G. Knepley   DMPolytopeType ct;
471*47bfdf3fSMatthew G. Knepley   PetscInt       dim, cStart;
472*47bfdf3fSMatthew G. Knepley   MPI_Comm       comm;
473*47bfdf3fSMatthew G. Knepley 
474*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
475*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
476*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
477*47bfdf3fSMatthew G. Knepley 
478*47bfdf3fSMatthew G. Knepley   /* Create finite element */
479*47bfdf3fSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
480*47bfdf3fSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
481*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFECreateByCell(comm, dim, 1, ct, "primal_potential_", PETSC_DEFAULT, &fe[0]));
482*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe[0], "primal potential"));
483*47bfdf3fSMatthew G. Knepley 
484*47bfdf3fSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
485*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
486*47bfdf3fSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
487*47bfdf3fSMatthew G. Knepley   while (cdm) {
488*47bfdf3fSMatthew G. Knepley     PetscCall(DMCopyDisc(dm, cdm));
489*47bfdf3fSMatthew G. Knepley     PetscCall(DMGetCoarseDM(cdm, &cdm));
490*47bfdf3fSMatthew G. Knepley   }
491*47bfdf3fSMatthew G. Knepley 
492*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe[0]));
493*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
494*47bfdf3fSMatthew G. Knepley }
495*47bfdf3fSMatthew G. Knepley 
496*47bfdf3fSMatthew G. Knepley static PetscErrorCode SetupMixedProblem(DM dm, AppCtx *user)
497*47bfdf3fSMatthew G. Knepley {
498*47bfdf3fSMatthew G. Knepley   PetscDS       ds;
499*47bfdf3fSMatthew G. Knepley   DMLabel       label;
500*47bfdf3fSMatthew G. Knepley   PetscInt      id, bd;
501*47bfdf3fSMatthew G. Knepley   Parameter    *param;
502*47bfdf3fSMatthew G. Knepley   PetscWeakForm wf;
503*47bfdf3fSMatthew G. Knepley 
504*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
505*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetDS(dm, &ds));
506*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetLabel(dm, "marker", &label));
507*47bfdf3fSMatthew G. Knepley 
508*47bfdf3fSMatthew G. Knepley   /* Residual terms */
509*47bfdf3fSMatthew G. Knepley   PetscCall(PetscDSSetResidual(ds, 0, f0_mixed_q, f1_mixed_q));
510*47bfdf3fSMatthew G. Knepley 
511*47bfdf3fSMatthew G. Knepley   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mixed_qq, NULL, NULL, NULL));
512*47bfdf3fSMatthew G. Knepley   PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_mixed_qu, NULL));
513*47bfdf3fSMatthew G. Knepley 
514*47bfdf3fSMatthew G. Knepley   PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_mixed_uq, NULL, NULL));
515*47bfdf3fSMatthew G. Knepley 
516*47bfdf3fSMatthew G. Knepley   switch (user->solType) {
517*47bfdf3fSMatthew G. Knepley   case SOL_QUADRATIC:
518*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_quadratic_u, NULL));
519*47bfdf3fSMatthew G. Knepley 
520*47bfdf3fSMatthew G. Knepley     id = 1;
521*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
522*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
523*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_mixed_q, 0, NULL));
524*47bfdf3fSMatthew G. Knepley 
525*47bfdf3fSMatthew G. Knepley     id = 2;
526*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_q, NULL, user, NULL));
527*47bfdf3fSMatthew G. Knepley 
528*47bfdf3fSMatthew G. Knepley     id = 4;
529*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_q, NULL, user, NULL));
530*47bfdf3fSMatthew G. Knepley 
531*47bfdf3fSMatthew G. Knepley     id = 3;
532*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
533*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
534*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_mixed_q, 0, NULL));
535*47bfdf3fSMatthew G. Knepley 
536*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user));
537*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user));
538*47bfdf3fSMatthew G. Knepley     break;
539*47bfdf3fSMatthew G. Knepley   case SOL_TRIG:
540*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_trig_u, NULL));
541*47bfdf3fSMatthew G. Knepley 
542*47bfdf3fSMatthew G. Knepley     id = 1;
543*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
544*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
545*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_mixed_q, 0, NULL));
546*47bfdf3fSMatthew G. Knepley 
547*47bfdf3fSMatthew G. Knepley     id = 2;
548*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))trig_q, NULL, user, NULL));
549*47bfdf3fSMatthew G. Knepley 
550*47bfdf3fSMatthew G. Knepley     id = 4;
551*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))trig_q, NULL, user, NULL));
552*47bfdf3fSMatthew G. Knepley 
553*47bfdf3fSMatthew G. Knepley     id = 3;
554*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
555*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
556*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_mixed_q, 0, NULL));
557*47bfdf3fSMatthew G. Knepley 
558*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, trig_q, user));
559*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 1, trig_u, user));
560*47bfdf3fSMatthew G. Knepley     break;
561*47bfdf3fSMatthew G. Knepley   case SOL_SENSOR:
562*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_sensor_u, NULL));
563*47bfdf3fSMatthew G. Knepley 
564*47bfdf3fSMatthew G. Knepley     id = 1;
565*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
566*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
567*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_mixed_q, 0, NULL));
568*47bfdf3fSMatthew G. Knepley 
569*47bfdf3fSMatthew G. Knepley     id = 2;
570*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))sensor_q, NULL, user, NULL));
571*47bfdf3fSMatthew G. Knepley 
572*47bfdf3fSMatthew G. Knepley     id = 4;
573*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))sensor_q, NULL, user, NULL));
574*47bfdf3fSMatthew G. Knepley 
575*47bfdf3fSMatthew G. Knepley     id = 3;
576*47bfdf3fSMatthew G. Knepley     PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd));
577*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
578*47bfdf3fSMatthew G. Knepley     PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_mixed_q, 0, NULL));
579*47bfdf3fSMatthew G. Knepley 
580*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 0, sensor_q, user));
581*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetExactSolution(ds, 1, sensor_u, user));
582*47bfdf3fSMatthew G. Knepley     break;
583*47bfdf3fSMatthew G. Knepley   default:
584*47bfdf3fSMatthew G. Knepley     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
585*47bfdf3fSMatthew G. Knepley   }
586*47bfdf3fSMatthew G. Knepley 
587*47bfdf3fSMatthew G. Knepley   /* Setup constants */
588*47bfdf3fSMatthew G. Knepley   {
589*47bfdf3fSMatthew G. Knepley     PetscCall(PetscBagGetData(user->param, (void **)&param));
590*47bfdf3fSMatthew G. Knepley     PetscScalar constants[1];
591*47bfdf3fSMatthew G. Knepley 
592*47bfdf3fSMatthew G. Knepley     constants[0] = param->k;
593*47bfdf3fSMatthew G. Knepley 
594*47bfdf3fSMatthew G. Knepley     PetscCall(PetscDSSetConstants(ds, 1, constants));
595*47bfdf3fSMatthew G. Knepley   }
596*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
597*47bfdf3fSMatthew G. Knepley }
598*47bfdf3fSMatthew G. Knepley 
599*47bfdf3fSMatthew G. Knepley static PetscErrorCode SetupMixedDiscretization(DM dm, AppCtx *user)
600*47bfdf3fSMatthew G. Knepley {
601*47bfdf3fSMatthew G. Knepley   DM             cdm = dm;
602*47bfdf3fSMatthew G. Knepley   PetscFE        fe[2];
603*47bfdf3fSMatthew G. Knepley   DMPolytopeType ct;
604*47bfdf3fSMatthew G. Knepley   PetscInt       dim, cStart;
605*47bfdf3fSMatthew G. Knepley   MPI_Comm       comm;
606*47bfdf3fSMatthew G. Knepley 
607*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
608*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
609*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetDimension(dm, &dim));
610*47bfdf3fSMatthew G. Knepley 
611*47bfdf3fSMatthew G. Knepley   /* Create finite element */
612*47bfdf3fSMatthew G. Knepley   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
613*47bfdf3fSMatthew G. Knepley   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
614*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFECreateByCell(comm, dim, dim, ct, "mixed_flux_", PETSC_DEFAULT, &fe[0]));
615*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe[0], "mixed flux"));
616*47bfdf3fSMatthew G. Knepley   /* NOTE:  Set the same quadrature order as the primal problem here or use the command line option. */
617*47bfdf3fSMatthew G. Knepley 
618*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFECreateByCell(comm, dim, 1, ct, "mixed_potential_", PETSC_DEFAULT, &fe[1]));
619*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
620*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)fe[1], "mixed potential"));
621*47bfdf3fSMatthew G. Knepley 
622*47bfdf3fSMatthew G. Knepley   /* Set discretization and boundary conditions for each mesh */
623*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
624*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
625*47bfdf3fSMatthew G. Knepley   PetscCall(DMCreateDS(dm));
626*47bfdf3fSMatthew G. Knepley   while (cdm) {
627*47bfdf3fSMatthew G. Knepley     PetscCall(DMCopyDisc(dm, cdm));
628*47bfdf3fSMatthew G. Knepley     PetscCall(DMGetCoarseDM(cdm, &cdm));
629*47bfdf3fSMatthew G. Knepley   }
630*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe[0]));
631*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFEDestroy(&fe[1]));
632*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
633*47bfdf3fSMatthew G. Knepley }
634*47bfdf3fSMatthew G. Knepley 
635*47bfdf3fSMatthew G. Knepley PetscErrorCode SetupMixed(DMAdaptor adaptor, DM mdm)
636*47bfdf3fSMatthew G. Knepley {
637*47bfdf3fSMatthew G. Knepley   AppCtx *ctx;
638*47bfdf3fSMatthew G. Knepley 
639*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
640*47bfdf3fSMatthew G. Knepley   PetscCall(DMGetApplicationContext(mdm, (void **)&ctx));
641*47bfdf3fSMatthew G. Knepley   PetscCall(SetupMixedDiscretization(mdm, ctx));
642*47bfdf3fSMatthew G. Knepley   PetscCall(SetupMixedProblem(mdm, ctx));
643*47bfdf3fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
644*47bfdf3fSMatthew G. Knepley }
645*47bfdf3fSMatthew G. Knepley 
646*47bfdf3fSMatthew G. Knepley int main(int argc, char **argv)
647*47bfdf3fSMatthew G. Knepley {
648*47bfdf3fSMatthew G. Knepley   DM        dm, mdm;               /* problem specification */
649*47bfdf3fSMatthew G. Knepley   SNES      snes, msnes;           /* nonlinear solvers */
650*47bfdf3fSMatthew G. Knepley   Vec       u, mu;                 /* solution vectors */
651*47bfdf3fSMatthew G. Knepley   Vec       fluxError, exactError; /* Element wise error vector */
652*47bfdf3fSMatthew G. Knepley   PetscReal fluxNorm, exactNorm;   /* Flux error norm */
653*47bfdf3fSMatthew G. Knepley   AppCtx    user;                  /* user-defined work context */
654*47bfdf3fSMatthew G. Knepley 
655*47bfdf3fSMatthew G. Knepley   PetscFunctionBeginUser;
656*47bfdf3fSMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
657*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.param));
658*47bfdf3fSMatthew G. Knepley   PetscCall(SetupParameters(user.param, &user));
659*47bfdf3fSMatthew G. Knepley   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
660*47bfdf3fSMatthew G. Knepley 
661*47bfdf3fSMatthew G. Knepley   // Set up and solve primal problem
662*47bfdf3fSMatthew G. Knepley   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
663*47bfdf3fSMatthew G. Knepley   PetscCall(DMSetApplicationContext(dm, &user));
664*47bfdf3fSMatthew G. Knepley   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
665*47bfdf3fSMatthew G. Knepley   PetscCall(SNESSetDM(snes, dm));
666*47bfdf3fSMatthew G. Knepley 
667*47bfdf3fSMatthew G. Knepley   PetscCall(SetupPrimalDiscretization(dm, &user));
668*47bfdf3fSMatthew G. Knepley   PetscCall(SetupPrimalProblem(dm, &user));
669*47bfdf3fSMatthew G. Knepley   PetscCall(DMCreateGlobalVector(dm, &u));
670*47bfdf3fSMatthew G. Knepley   PetscCall(VecSet(u, 0.0));
671*47bfdf3fSMatthew G. Knepley   PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
672*47bfdf3fSMatthew G. Knepley   PetscCall(SNESSetFromOptions(snes));
673*47bfdf3fSMatthew G. Knepley 
674*47bfdf3fSMatthew G. Knepley   PetscCall(DMSNESCheckFromOptions(snes, u));
675*47bfdf3fSMatthew G. Knepley   PetscCall(SNESSolve(snes, NULL, u));
676*47bfdf3fSMatthew G. Knepley 
677*47bfdf3fSMatthew G. Knepley   // Needed if you allow SNES to refine
678*47bfdf3fSMatthew G. Knepley   PetscCall(SNESGetSolution(snes, &u));
679*47bfdf3fSMatthew G. Knepley   PetscCall(VecGetDM(u, &dm));
680*47bfdf3fSMatthew G. Knepley 
681*47bfdf3fSMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)u, "Primal Solution "));
682*47bfdf3fSMatthew G. Knepley   PetscCall(VecViewFromOptions(u, NULL, "-primal_sol_vec_view"));
683*47bfdf3fSMatthew G. Knepley 
684*47bfdf3fSMatthew G. Knepley   if (user.byHand) {
685*47bfdf3fSMatthew G. Knepley     // Set up and solve mixed problem
686*47bfdf3fSMatthew G. Knepley     PetscCall(DMClone(dm, &mdm));
687*47bfdf3fSMatthew G. Knepley     PetscCall(SNESCreate(PETSC_COMM_WORLD, &msnes));
688*47bfdf3fSMatthew G. Knepley     PetscCall(SNESSetOptionsPrefix(msnes, "mixed_"));
689*47bfdf3fSMatthew G. Knepley     PetscCall(SNESSetDM(msnes, mdm));
690*47bfdf3fSMatthew G. Knepley 
691*47bfdf3fSMatthew G. Knepley     PetscCall(SetupMixedDiscretization(mdm, &user));
692*47bfdf3fSMatthew G. Knepley     PetscCall(SetupMixedProblem(mdm, &user));
693*47bfdf3fSMatthew G. Knepley     PetscCall(DMCreateGlobalVector(mdm, &mu));
694*47bfdf3fSMatthew G. Knepley     PetscCall(VecSet(mu, 0.0));
695*47bfdf3fSMatthew G. Knepley     PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, &user));
696*47bfdf3fSMatthew G. Knepley     PetscCall(SNESSetFromOptions(msnes));
697*47bfdf3fSMatthew G. Knepley 
698*47bfdf3fSMatthew G. Knepley     PetscCall(DMSNESCheckFromOptions(msnes, mu));
699*47bfdf3fSMatthew G. Knepley     PetscCall(SNESSolve(msnes, NULL, mu));
700*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution "));
701*47bfdf3fSMatthew G. Knepley     PetscCall(VecViewFromOptions(mu, NULL, "-mixed_sol_vec_view"));
702*47bfdf3fSMatthew G. Knepley 
703*47bfdf3fSMatthew G. Knepley     // Create the error space of piecewise constants
704*47bfdf3fSMatthew G. Knepley     DM             dmErr;
705*47bfdf3fSMatthew G. Knepley     PetscFE        feErr;
706*47bfdf3fSMatthew G. Knepley     DMPolytopeType ct;
707*47bfdf3fSMatthew G. Knepley     PetscInt       dim, cStart;
708*47bfdf3fSMatthew G. Knepley 
709*47bfdf3fSMatthew G. Knepley     PetscCall(DMClone(dm, &dmErr));
710*47bfdf3fSMatthew G. Knepley     PetscCall(DMGetDimension(dmErr, &dim));
711*47bfdf3fSMatthew G. Knepley     PetscCall(DMPlexGetHeightStratum(dmErr, 0, &cStart, NULL));
712*47bfdf3fSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dmErr, cStart, &ct));
713*47bfdf3fSMatthew G. Knepley     PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &feErr));
714*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)feErr, "Error"));
715*47bfdf3fSMatthew G. Knepley     PetscCall(DMSetField(dmErr, 0, NULL, (PetscObject)feErr));
716*47bfdf3fSMatthew G. Knepley     PetscCall(PetscFEDestroy(&feErr));
717*47bfdf3fSMatthew G. Knepley     PetscCall(DMCreateDS(dmErr));
718*47bfdf3fSMatthew G. Knepley     PetscCall(DMViewFromOptions(dmErr, NULL, "-dmerr_view"));
719*47bfdf3fSMatthew G. Knepley 
720*47bfdf3fSMatthew G. Knepley     // Compute the flux norm
721*47bfdf3fSMatthew G. Knepley     PetscCall(DMGetGlobalVector(dmErr, &fluxError));
722*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)fluxError, "Flux Error"));
723*47bfdf3fSMatthew G. Knepley     PetscCall(DMGetGlobalVector(dmErr, &exactError));
724*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)exactError, "Analytical Error"));
725*47bfdf3fSMatthew G. Knepley     PetscCall(DMPlexComputeL2FluxDiffVec(u, 0, mu, 0, fluxError));
726*47bfdf3fSMatthew G. Knepley     {
727*47bfdf3fSMatthew G. Knepley       PetscDS             ds;
728*47bfdf3fSMatthew G. Knepley       PetscSimplePointFn *func[2] = {NULL, NULL};
729*47bfdf3fSMatthew G. Knepley       void               *ctx[2]  = {NULL, NULL};
730*47bfdf3fSMatthew G. Knepley 
731*47bfdf3fSMatthew G. Knepley       PetscCall(DMGetDS(mdm, &ds));
732*47bfdf3fSMatthew G. Knepley       PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctx[0]));
733*47bfdf3fSMatthew G. Knepley       PetscCall(DMPlexComputeL2DiffVec(mdm, 0.0, func, ctx, mu, exactError));
734*47bfdf3fSMatthew G. Knepley     }
735*47bfdf3fSMatthew G. Knepley     PetscCall(VecNorm(fluxError, NORM_2, &fluxNorm));
736*47bfdf3fSMatthew G. Knepley     PetscCall(VecNorm(exactError, NORM_2, &exactNorm));
737*47bfdf3fSMatthew G. Knepley     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Flux error norm = %g\t Exact flux error norm = %g\n", (double)fluxNorm, (double)exactNorm));
738*47bfdf3fSMatthew G. Knepley     PetscCall(VecViewFromOptions(fluxError, NULL, "-flux_error_vec_view"));
739*47bfdf3fSMatthew G. Knepley     PetscCall(VecViewFromOptions(exactError, NULL, "-exact_error_vec_view"));
740*47bfdf3fSMatthew G. Knepley 
741*47bfdf3fSMatthew G. Knepley     // Adaptive refinement based on calculated error
742*47bfdf3fSMatthew G. Knepley     DM        rdm;
743*47bfdf3fSMatthew G. Knepley     VecTagger refineTag;
744*47bfdf3fSMatthew G. Knepley     DMLabel   adaptLabel;
745*47bfdf3fSMatthew G. Knepley     IS        refineIS;
746*47bfdf3fSMatthew G. Knepley     Vec       ref;
747*47bfdf3fSMatthew G. Knepley 
748*47bfdf3fSMatthew G. Knepley     PetscCall(DMLabelCreate(PETSC_COMM_WORLD, "adapt", &adaptLabel));
749*47bfdf3fSMatthew G. Knepley     PetscCall(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN));
750*47bfdf3fSMatthew G. Knepley     PetscCall(VecTaggerCreate(PETSC_COMM_WORLD, &refineTag));
751*47bfdf3fSMatthew G. Knepley     PetscCall(VecTaggerSetFromOptions(refineTag));
752*47bfdf3fSMatthew G. Knepley     PetscCall(VecTaggerSetUp(refineTag));
753*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view"));
754*47bfdf3fSMatthew G. Knepley 
755*47bfdf3fSMatthew G. Knepley     PetscCall(VecTaggerComputeIS(refineTag, fluxError, &refineIS, NULL));
756*47bfdf3fSMatthew G. Knepley     PetscCall(VecTaggerDestroy(&refineTag));
757*47bfdf3fSMatthew G. Knepley     PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS));
758*47bfdf3fSMatthew G. Knepley     PetscCall(ISViewFromOptions(refineIS, NULL, "-refine_is_view"));
759*47bfdf3fSMatthew G. Knepley     PetscCall(ISDestroy(&refineIS));
760*47bfdf3fSMatthew G. Knepley 
761*47bfdf3fSMatthew G. Knepley     PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref));
762*47bfdf3fSMatthew G. Knepley     PetscCall(VecViewFromOptions(ref, NULL, "-refine_vec_view"));
763*47bfdf3fSMatthew G. Knepley     PetscCall(VecDestroy(&ref));
764*47bfdf3fSMatthew G. Knepley 
765*47bfdf3fSMatthew G. Knepley     // Mark adaptation phase with prefix ref_
766*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "adapt_"));
767*47bfdf3fSMatthew G. Knepley     PetscCall(DMAdaptLabel(dm, adaptLabel, &rdm));
768*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, NULL));
769*47bfdf3fSMatthew G. Knepley     PetscCall(PetscObjectSetName((PetscObject)rdm, "Adaptively Refined DM"));
770*47bfdf3fSMatthew G. Knepley     PetscCall(DMViewFromOptions(rdm, NULL, "-adapt_dm_view"));
771*47bfdf3fSMatthew G. Knepley     PetscCall(DMDestroy(&rdm));
772*47bfdf3fSMatthew G. Knepley     PetscCall(DMLabelDestroy(&adaptLabel));
773*47bfdf3fSMatthew G. Knepley 
774*47bfdf3fSMatthew G. Knepley     // Destroy the error structures
775*47bfdf3fSMatthew G. Knepley     PetscCall(DMRestoreGlobalVector(dmErr, &fluxError));
776*47bfdf3fSMatthew G. Knepley     PetscCall(DMRestoreGlobalVector(dmErr, &exactError));
777*47bfdf3fSMatthew G. Knepley     PetscCall(DMDestroy(&dmErr));
778*47bfdf3fSMatthew G. Knepley 
779*47bfdf3fSMatthew G. Knepley     // Destroy the mixed structures
780*47bfdf3fSMatthew G. Knepley     PetscCall(VecDestroy(&mu));
781*47bfdf3fSMatthew G. Knepley     PetscCall(DMDestroy(&mdm));
782*47bfdf3fSMatthew G. Knepley     PetscCall(SNESDestroy(&msnes));
783*47bfdf3fSMatthew G. Knepley   }
784*47bfdf3fSMatthew G. Knepley 
785*47bfdf3fSMatthew G. Knepley   // Destroy the primal structures
786*47bfdf3fSMatthew G. Knepley   PetscCall(VecDestroy(&u));
787*47bfdf3fSMatthew G. Knepley   PetscCall(DMDestroy(&dm));
788*47bfdf3fSMatthew G. Knepley   PetscCall(SNESDestroy(&snes));
789*47bfdf3fSMatthew G. Knepley   PetscCall(PetscBagDestroy(&user.param));
790*47bfdf3fSMatthew G. Knepley   PetscCall(PetscFinalize());
791*47bfdf3fSMatthew G. Knepley   return 0;
792*47bfdf3fSMatthew G. Knepley }
793*47bfdf3fSMatthew G. Knepley 
794*47bfdf3fSMatthew G. Knepley /*TEST
795*47bfdf3fSMatthew G. Knepley 
796*47bfdf3fSMatthew G. Knepley   # Tests using the explicit code above
797*47bfdf3fSMatthew G. Knepley   testset:
798*47bfdf3fSMatthew G. Knepley     suffix: 2d_p2_rt0p0_byhand
799*47bfdf3fSMatthew G. Knepley     requires: triangle
800*47bfdf3fSMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -dm_plex_box_lower 0,0 -dm_plex_box_upper 1,1 -dm_refine 3 \
801*47bfdf3fSMatthew G. Knepley           -primal_potential_petscspace_degree 2 \
802*47bfdf3fSMatthew G. Knepley           -mixed_potential_petscdualspace_lagrange_continuity 0 \
803*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_type ptrimmed \
804*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_components 2 \
805*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_ptrimmed_form_degree -1 \
806*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_order 1 \
807*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_form_degree -1 \
808*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_lagrange_trimmed true \
809*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscfe_default_quadrature_order 2 \
810*47bfdf3fSMatthew G. Knepley           -vec_tagger_type cdf -vec_tagger_box 0.9,1.0 \
811*47bfdf3fSMatthew G. Knepley             -tag_view \
812*47bfdf3fSMatthew G. Knepley             -adapt_dm_adaptor cellrefiner -adapt_dm_plex_transform_type refine_sbr \
813*47bfdf3fSMatthew G. Knepley           -dmsnes_check 0.001 -mixed_dmsnes_check 0.001 -pc_type jacobi -mixed_pc_type jacobi
814*47bfdf3fSMatthew G. Knepley     test:
815*47bfdf3fSMatthew G. Knepley       suffix: quadratic
816*47bfdf3fSMatthew G. Knepley       args: -sol_type quadratic
817*47bfdf3fSMatthew G. Knepley     test:
818*47bfdf3fSMatthew G. Knepley       suffix: trig
819*47bfdf3fSMatthew G. Knepley       args: -sol_type trig
820*47bfdf3fSMatthew G. Knepley     test:
821*47bfdf3fSMatthew G. Knepley       suffix: sensor
822*47bfdf3fSMatthew G. Knepley       args: -sol_type sensor
823*47bfdf3fSMatthew G. Knepley 
824*47bfdf3fSMatthew G. Knepley   # Tests using the embedded adaptor in SNES
825*47bfdf3fSMatthew G. Knepley   testset:
826*47bfdf3fSMatthew G. Knepley     suffix: 2d_p2_rt0p0
827*47bfdf3fSMatthew G. Knepley     requires: triangle defined(PETSC_HAVE_EXECUTABLE_EXPORT)
828*47bfdf3fSMatthew G. Knepley     args: -dm_plex_box_faces 1,1 -dm_plex_box_lower 0,0 -dm_plex_box_upper 1,1 -dm_refine 3 \
829*47bfdf3fSMatthew G. Knepley           -primal_potential_petscspace_degree 2 \
830*47bfdf3fSMatthew G. Knepley           -mixed_potential_petscdualspace_lagrange_continuity 0 \
831*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_type ptrimmed \
832*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_components 2 \
833*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscspace_ptrimmed_form_degree -1 \
834*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_order 1 \
835*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_form_degree -1 \
836*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscdualspace_lagrange_trimmed true \
837*47bfdf3fSMatthew G. Knepley           -mixed_flux_petscfe_default_quadrature_order 2 \
838*47bfdf3fSMatthew G. Knepley           -by_hand 0 \
839*47bfdf3fSMatthew G. Knepley           -refine_vec_tagger_type cdf -refine_vec_tagger_box 0.9,1.0 \
840*47bfdf3fSMatthew G. Knepley             -snes_adapt_view \
841*47bfdf3fSMatthew G. Knepley             -adapt_dm_adaptor cellrefiner -adapt_dm_plex_transform_type refine_sbr \
842*47bfdf3fSMatthew G. Knepley             -adaptor_criterion label -adaptor_type flux -adaptor_mixed_setup_function SetupMixed \
843*47bfdf3fSMatthew G. Knepley           -snes_adapt_sequence 1 -pc_type jacobi -mixed_pc_type jacobi
844*47bfdf3fSMatthew G. Knepley     test:
845*47bfdf3fSMatthew G. Knepley       suffix: quadratic
846*47bfdf3fSMatthew G. Knepley       args: -sol_type quadratic -adaptor_monitor_error
847*47bfdf3fSMatthew G. Knepley     test:
848*47bfdf3fSMatthew G. Knepley       suffix: trig
849*47bfdf3fSMatthew G. Knepley       args: -sol_type trig -adaptor_monitor_error
850*47bfdf3fSMatthew G. Knepley     test:
851*47bfdf3fSMatthew G. Knepley       suffix: sensor
852*47bfdf3fSMatthew G. Knepley       args: -sol_type sensor
853*47bfdf3fSMatthew G. Knepley 
854*47bfdf3fSMatthew G. Knepley TEST*/
855