xref: /petsc/src/snes/tutorials/ex76.c (revision 30602db00db74b7e41a0c75e517aefe6711423f0)
1649ef022SMatthew Knepley static char help[] = "Low Mach Flow in 2d and 3d channels with finite elements.\n\
2649ef022SMatthew Knepley We solve the Low Mach flow problem in a rectangular\n\
3649ef022SMatthew Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4649ef022SMatthew Knepley 
5649ef022SMatthew Knepley /*F
6649ef022SMatthew Knepley This Low Mach flow is a steady-state isoviscous Navier-Stokes flow. We discretize using the
7649ef022SMatthew Knepley finite element method on an unstructured mesh. The weak form equations are
8649ef022SMatthew Knepley 
9649ef022SMatthew Knepley \begin{align*}
10649ef022SMatthew Knepley     < q, \nabla\cdot u >                                                                                     = 0
11649ef022SMatthew Knepley     <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
12649ef022SMatthew Knepley     < w, u \cdot \nabla T > - < \nabla w, \alpha \nabla T > - < w, Q >                                       = 0
13649ef022SMatthew Knepley \end{align*}
14649ef022SMatthew Knepley 
15649ef022SMatthew Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
16649ef022SMatthew Knepley 
17649ef022SMatthew Knepley For visualization, use
18649ef022SMatthew Knepley 
19649ef022SMatthew Knepley   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
20649ef022SMatthew Knepley F*/
21649ef022SMatthew Knepley 
22649ef022SMatthew Knepley #include <petscdmplex.h>
23649ef022SMatthew Knepley #include <petscsnes.h>
24649ef022SMatthew Knepley #include <petscds.h>
25649ef022SMatthew Knepley #include <petscbag.h>
26649ef022SMatthew Knepley 
27649ef022SMatthew Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType;
28649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic",  "unknown"};
29649ef022SMatthew Knepley 
30649ef022SMatthew Knepley typedef struct {
31649ef022SMatthew Knepley   PetscReal nu;      /* Kinematic viscosity */
32649ef022SMatthew Knepley   PetscReal theta;   /* Angle of pipe wall to x-axis */
33649ef022SMatthew Knepley   PetscReal alpha;   /* Thermal diffusivity */
34649ef022SMatthew Knepley   PetscReal T_in;    /* Inlet temperature*/
35649ef022SMatthew Knepley } Parameter;
36649ef022SMatthew Knepley 
37649ef022SMatthew Knepley typedef struct {
38*30602db0SMatthew G. Knepley   PetscBool showError;
39*30602db0SMatthew G. Knepley   PetscBag  bag;
40649ef022SMatthew Knepley   SolType   solType;
41649ef022SMatthew Knepley } AppCtx;
42649ef022SMatthew Knepley 
43649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
44649ef022SMatthew Knepley {
45649ef022SMatthew Knepley   PetscInt d;
46649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
47649ef022SMatthew Knepley   return 0;
48649ef022SMatthew Knepley }
49649ef022SMatthew Knepley 
50649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51649ef022SMatthew Knepley {
52649ef022SMatthew Knepley   PetscInt d;
53649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
54649ef022SMatthew Knepley   return 0;
55649ef022SMatthew Knepley }
56649ef022SMatthew Knepley 
57649ef022SMatthew Knepley /*
58649ef022SMatthew Knepley   CASE: quadratic
59649ef022SMatthew Knepley   In 2D we use exact solution:
60649ef022SMatthew Knepley 
61649ef022SMatthew Knepley     u = x^2 + y^2
62649ef022SMatthew Knepley     v = 2x^2 - 2xy
63649ef022SMatthew Knepley     p = x + y - 1
64649ef022SMatthew Knepley     T = x + y
65649ef022SMatthew Knepley     f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -4\nu + 1>
66649ef022SMatthew Knepley     Q = 3x^2 + y^2 - 2xy
67649ef022SMatthew Knepley 
68649ef022SMatthew Knepley   so that
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley (1)  \nabla \cdot u  = 2x - 2x = 0
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley (2)  u \cdot \nabla u - \nu \Delta u + \nabla p - f
73649ef022SMatthew Knepley      = <2x^3 + 4x^2y -2xy^2, 4xy^2 + 2x^2y - 2y^3> -\nu <4, 4> + <1, 1> - <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -         4\nu + 1>  = 0
74649ef022SMatthew Knepley 
75649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0
76649ef022SMatthew Knepley */
77649ef022SMatthew Knepley 
78649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
79649ef022SMatthew Knepley {
80649ef022SMatthew Knepley   u[0] = X[0]*X[0] + X[1]*X[1];
81649ef022SMatthew Knepley   u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
82649ef022SMatthew Knepley   return 0;
83649ef022SMatthew Knepley }
84649ef022SMatthew Knepley 
85649ef022SMatthew Knepley static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
86649ef022SMatthew Knepley {
87649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
88649ef022SMatthew Knepley   return 0;
89649ef022SMatthew Knepley }
90649ef022SMatthew Knepley 
91649ef022SMatthew Knepley static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
92649ef022SMatthew Knepley {
93649ef022SMatthew Knepley   T[0] = X[0] + X[1];
94649ef022SMatthew Knepley   return 0;
95649ef022SMatthew Knepley }
96649ef022SMatthew Knepley 
97649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
98649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
99649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
101649ef022SMatthew Knepley {
102649ef022SMatthew Knepley   PetscInt                   c, d;
103649ef022SMatthew Knepley   PetscInt                   Nc = dim;
104649ef022SMatthew Knepley   const PetscReal    nu = PetscRealPart(constants[0]);
105649ef022SMatthew Knepley 
106649ef022SMatthew Knepley   for (c=0; c<Nc; ++c) {
107649ef022SMatthew Knepley     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
108649ef022SMatthew Knepley   }
109649ef022SMatthew Knepley   f0[0] -= (2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 1);
110649ef022SMatthew Knepley   f0[1] -= (4*X[0]*X[1]*X[1] + 2*X[0]*X[0]*X[1] - 2*X[1]*X[1]*X[1] - 4.0*nu + 1);
111649ef022SMatthew Knepley }
112649ef022SMatthew Knepley 
113649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
117649ef022SMatthew Knepley {
118649ef022SMatthew Knepley   PetscInt d;
119649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
120649ef022SMatthew Knepley   f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]);
121649ef022SMatthew Knepley }
122649ef022SMatthew Knepley 
123649ef022SMatthew Knepley 
124649ef022SMatthew Knepley /*
125649ef022SMatthew Knepley   CASE: cubic
126649ef022SMatthew Knepley   In 2D we use exact solution:
127649ef022SMatthew Knepley 
128649ef022SMatthew Knepley     u = x^3 + y^3
129649ef022SMatthew Knepley     v = 2x^3 - 3x^2y
130649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
131649ef022SMatthew Knepley     T = 1/2 x^2 + 1/2 y^2
132649ef022SMatthew Knepley     f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y>
133649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2
134649ef022SMatthew Knepley 
135649ef022SMatthew Knepley   so that
136649ef022SMatthew Knepley 
137649ef022SMatthew Knepley   \nabla \cdot u = 3x^2 - 3x^2 = 0
138649ef022SMatthew Knepley 
139649ef022SMatthew Knepley   u \cdot \nabla u - \nu \Delta u + \nabla p - f
140649ef022SMatthew Knepley   = <3x^5 + 6x^3y^2 - 6x^2y^3, 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> = 0
141649ef022SMatthew Knepley 
142649ef022SMatthew Knepley   u \cdot \nabla T - \alpha\Delta T - Q = (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2)   = 0
143649ef022SMatthew Knepley */
144649ef022SMatthew Knepley 
145649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
146649ef022SMatthew Knepley {
147649ef022SMatthew Knepley   u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
148649ef022SMatthew Knepley   u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
149649ef022SMatthew Knepley   return 0;
150649ef022SMatthew Knepley }
151649ef022SMatthew Knepley 
152649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
153649ef022SMatthew Knepley {
154649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
155649ef022SMatthew Knepley   return 0;
156649ef022SMatthew Knepley }
157649ef022SMatthew Knepley 
158649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
159649ef022SMatthew Knepley {
160649ef022SMatthew Knepley   T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
161649ef022SMatthew Knepley   return 0;
162649ef022SMatthew Knepley }
163649ef022SMatthew Knepley 
164649ef022SMatthew Knepley 
165649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
166649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
167649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
169649ef022SMatthew Knepley {
170649ef022SMatthew Knepley   PetscInt                   c, d;
171649ef022SMatthew Knepley   PetscInt                   Nc = dim;
172649ef022SMatthew Knepley   const PetscReal    nu = PetscRealPart(constants[0]);
173649ef022SMatthew Knepley 
174649ef022SMatthew Knepley   for (c=0; c<Nc; ++c) {
175649ef022SMatthew Knepley     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
176649ef022SMatthew Knepley   }
177649ef022SMatthew Knepley   f0[0] -= (3*X[0]*X[0]*X[0]*X[0]*X[0] + 6*X[0]*X[0]*X[0]*X[1]*X[1] - 6*X[0]*X[0]*X[1]*X[1]*X[1] - (6*X[0]+6*X[1])*nu + 3*X[0]);
178649ef022SMatthew Knepley   f0[1] -= (6*X[0]*X[0]*X[1]*X[1]*X[1] + 3*X[0]*X[0]*X[0]*X[0]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1]);
179649ef022SMatthew Knepley }
180649ef022SMatthew Knepley 
181649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
182649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
183649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
184649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
185649ef022SMatthew Knepley {
186649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
187649ef022SMatthew Knepley   PetscInt        d;
188649ef022SMatthew Knepley 
189649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
190649ef022SMatthew Knepley   f0[0] -= (X[0]*X[0]*X[0]*X[0] + X[0]*X[1]*X[1]*X[1] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] - 2.0*alpha);
191649ef022SMatthew Knepley }
192649ef022SMatthew Knepley 
193649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
194649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
195649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
197649ef022SMatthew Knepley {
198649ef022SMatthew Knepley   PetscInt d;
199649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
200649ef022SMatthew Knepley }
201649ef022SMatthew Knepley 
202649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
203649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
204649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
205649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
206649ef022SMatthew Knepley {
207649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
208649ef022SMatthew Knepley   const PetscInt  Nc = dim;
209649ef022SMatthew Knepley   PetscInt        c, d;
210649ef022SMatthew Knepley 
211649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
212649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
213649ef022SMatthew Knepley       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
214649ef022SMatthew Knepley       //f1[c*dim+d] = nu*u_x[c*dim+d];
215649ef022SMatthew Knepley     }
216649ef022SMatthew Knepley     f1[c*dim+c] -= u[uOff[1]];
217649ef022SMatthew Knepley   }
218649ef022SMatthew Knepley }
219649ef022SMatthew Knepley 
220649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
224649ef022SMatthew Knepley {
225649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
226649ef022SMatthew Knepley   PetscInt d;
227649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
228649ef022SMatthew Knepley }
229649ef022SMatthew Knepley 
230649ef022SMatthew Knepley /*Jacobians*/
231649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
232649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
233649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
234649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
235649ef022SMatthew Knepley {
236649ef022SMatthew Knepley   PetscInt d;
237649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
238649ef022SMatthew Knepley }
239649ef022SMatthew Knepley 
240649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
241649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
242649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
243649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
244649ef022SMatthew Knepley {
245649ef022SMatthew Knepley   const PetscInt  Nc = dim;
246649ef022SMatthew Knepley    PetscInt            c, d;
247649ef022SMatthew Knepley 
248649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
249649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
250649ef022SMatthew Knepley       g0[c*Nc+d] = u_x[ c*Nc+d];
251649ef022SMatthew Knepley     }
252649ef022SMatthew Knepley   }
253649ef022SMatthew Knepley }
254649ef022SMatthew Knepley 
255649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
256649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
257649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
258649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
259649ef022SMatthew Knepley {
260649ef022SMatthew Knepley   PetscInt NcI = dim;
261649ef022SMatthew Knepley   PetscInt NcJ = dim;
262649ef022SMatthew Knepley   PetscInt c, d, e;
263649ef022SMatthew Knepley 
264649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
265649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
266649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
267649ef022SMatthew Knepley         if (c == d) {
268649ef022SMatthew Knepley           g1[(c*NcJ+d)*dim+e] = u[e];
269649ef022SMatthew Knepley         }
270649ef022SMatthew Knepley       }
271649ef022SMatthew Knepley     }
272649ef022SMatthew Knepley   }
273649ef022SMatthew Knepley }
274649ef022SMatthew Knepley 
275649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
276649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
277649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
278649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
279649ef022SMatthew Knepley {
280649ef022SMatthew Knepley   PetscInt d;
281649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
282649ef022SMatthew Knepley }
283649ef022SMatthew Knepley 
284649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
285649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
286649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
287649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
288649ef022SMatthew Knepley {
289649ef022SMatthew Knepley    const PetscReal      nu = PetscRealPart(constants[0]);
290649ef022SMatthew Knepley    const PetscInt         Nc = dim;
291649ef022SMatthew Knepley    PetscInt                     c, d;
292649ef022SMatthew Knepley 
293649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
294649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
295649ef022SMatthew Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU
296649ef022SMatthew Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose
297649ef022SMatthew Knepley     }
298649ef022SMatthew Knepley   }
299649ef022SMatthew Knepley }
300649ef022SMatthew Knepley 
301649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
302649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
303649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
304649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
305649ef022SMatthew Knepley {
306649ef022SMatthew Knepley   PetscInt d;
307649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
308649ef022SMatthew Knepley }
309649ef022SMatthew Knepley 
310649ef022SMatthew Knepley 
311649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
312649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
313649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
314649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
315649ef022SMatthew Knepley {
316649ef022SMatthew Knepley   PetscInt d;
317649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
318649ef022SMatthew Knepley }
319649ef022SMatthew Knepley 
320649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
321649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
322649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
323649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
324649ef022SMatthew Knepley {
325649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
326649ef022SMatthew Knepley   PetscInt        d;
327649ef022SMatthew Knepley 
328649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
329649ef022SMatthew Knepley }
330649ef022SMatthew Knepley 
331649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
332649ef022SMatthew Knepley {
333*30602db0SMatthew G. Knepley   PetscInt       sol;
334649ef022SMatthew Knepley   PetscErrorCode ierr;
335649ef022SMatthew Knepley 
336649ef022SMatthew Knepley 
337649ef022SMatthew Knepley   PetscFunctionBeginUser;
338649ef022SMatthew Knepley   options->solType   = SOL_QUADRATIC;
339649ef022SMatthew Knepley   options->showError = PETSC_FALSE;
340649ef022SMatthew Knepley 
341649ef022SMatthew Knepley   ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr);
342649ef022SMatthew Knepley   sol = options->solType;
343649ef022SMatthew Knepley   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
344649ef022SMatthew Knepley   options->solType = (SolType) sol;
345649ef022SMatthew Knepley   ierr = PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);CHKERRQ(ierr);
346649ef022SMatthew Knepley   ierr = PetscOptionsEnd();
347649ef022SMatthew Knepley   PetscFunctionReturn(0);
348649ef022SMatthew Knepley }
349649ef022SMatthew Knepley 
350649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user)
351649ef022SMatthew Knepley {
352649ef022SMatthew Knepley   PetscBag       bag;
353649ef022SMatthew Knepley   Parameter     *p;
354649ef022SMatthew Knepley   PetscErrorCode ierr;
355649ef022SMatthew Knepley 
356649ef022SMatthew Knepley   PetscFunctionBeginUser;
357649ef022SMatthew Knepley   /* setup PETSc parameter bag */
358649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
359649ef022SMatthew Knepley   ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr);
360649ef022SMatthew Knepley   bag  = user->bag;
361649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0,   "nu",      "Kinematic viscosity");CHKERRQ(ierr);
362649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0,   "alpha",   "Thermal diffusivity");CHKERRQ(ierr);
363649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->theta, 0.0,   "theta",   "Angle of pipe wall to x-axis");CHKERRQ(ierr);
364649ef022SMatthew Knepley   PetscFunctionReturn(0);
365649ef022SMatthew Knepley }
366649ef022SMatthew Knepley 
367649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
368649ef022SMatthew Knepley {
369649ef022SMatthew Knepley   PetscErrorCode ierr;
370649ef022SMatthew Knepley 
371649ef022SMatthew Knepley   PetscFunctionBeginUser;
372*30602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
373*30602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
374*30602db0SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
375649ef022SMatthew Knepley   {
376649ef022SMatthew Knepley     Parameter   *param;
377649ef022SMatthew Knepley     Vec          coordinates;
378649ef022SMatthew Knepley     PetscScalar *coords;
379649ef022SMatthew Knepley     PetscReal    theta;
380649ef022SMatthew Knepley     PetscInt     cdim, N, bs, i;
381649ef022SMatthew Knepley 
382649ef022SMatthew Knepley     ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr);
383649ef022SMatthew Knepley     ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr);
384649ef022SMatthew Knepley     ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
385649ef022SMatthew Knepley     ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr);
386649ef022SMatthew Knepley     if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim);
387649ef022SMatthew Knepley     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
388649ef022SMatthew Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
389649ef022SMatthew Knepley     theta = param->theta;
390649ef022SMatthew Knepley     for (i = 0; i < N; i += cdim) {
391649ef022SMatthew Knepley       PetscScalar x = coords[i+0];
392649ef022SMatthew Knepley       PetscScalar y = coords[i+1];
393649ef022SMatthew Knepley 
394649ef022SMatthew Knepley       coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y;
395649ef022SMatthew Knepley       coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y;
396649ef022SMatthew Knepley     }
397649ef022SMatthew Knepley     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
398649ef022SMatthew Knepley     ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr);
399649ef022SMatthew Knepley   }
400649ef022SMatthew Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
401649ef022SMatthew Knepley   PetscFunctionReturn(0);
402649ef022SMatthew Knepley }
403649ef022SMatthew Knepley 
404649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
405649ef022SMatthew Knepley {
406649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
407649ef022SMatthew Knepley   PetscDS          prob;
40845480ffeSMatthew G. Knepley   DMLabel          label;
409649ef022SMatthew Knepley   Parameter       *ctx;
410649ef022SMatthew Knepley   PetscInt         id;
411649ef022SMatthew Knepley   PetscErrorCode   ierr;
412649ef022SMatthew Knepley 
413649ef022SMatthew Knepley   PetscFunctionBeginUser;
41445480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
415649ef022SMatthew Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
416649ef022SMatthew Knepley   switch(user->solType){
417649ef022SMatthew Knepley   case SOL_QUADRATIC:
418649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr);
419649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr);
420649ef022SMatthew Knepley 
421649ef022SMatthew Knepley     exactFuncs[0] = quadratic_u;
422649ef022SMatthew Knepley     exactFuncs[1] = linear_p;
423649ef022SMatthew Knepley     exactFuncs[2] = linear_T;
424649ef022SMatthew Knepley     break;
425649ef022SMatthew Knepley   case SOL_CUBIC:
426649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr);
427649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr);
428649ef022SMatthew Knepley 
429649ef022SMatthew Knepley     exactFuncs[0] = cubic_u;
430649ef022SMatthew Knepley     exactFuncs[1] = quadratic_p;
431649ef022SMatthew Knepley     exactFuncs[2] = quadratic_T;
432649ef022SMatthew Knepley     break;
433649ef022SMatthew Knepley    default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
434649ef022SMatthew Knepley   }
435649ef022SMatthew Knepley 
436649ef022SMatthew Knepley   ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr);
437649ef022SMatthew Knepley 
438649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu);CHKERRQ(ierr);
439649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL);CHKERRQ(ierr);
440649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL);CHKERRQ(ierr);
441649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL);CHKERRQ(ierr);
442649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
443649ef022SMatthew Knepley   /* Setup constants */
444649ef022SMatthew Knepley   {
445649ef022SMatthew Knepley     Parameter  *param;
446649ef022SMatthew Knepley     PetscScalar constants[3];
447649ef022SMatthew Knepley 
448649ef022SMatthew Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
449649ef022SMatthew Knepley 
450649ef022SMatthew Knepley     constants[0] = param->nu;
451649ef022SMatthew Knepley     constants[1] = param->alpha;
452649ef022SMatthew Knepley     constants[2] = param->theta;
453649ef022SMatthew Knepley     ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr);
454649ef022SMatthew Knepley   }
455649ef022SMatthew Knepley   /* Setup Boundary Conditions */
456649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
457649ef022SMatthew Knepley   id   = 3;
45845480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity",    label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
459649ef022SMatthew Knepley   id   = 1;
46045480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
461649ef022SMatthew Knepley   id   = 2;
46245480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity",  label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
463649ef022SMatthew Knepley   id   = 4;
46445480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity",   label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
465649ef022SMatthew Knepley   id   = 3;
46645480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp",    label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
467649ef022SMatthew Knepley   id   = 1;
46845480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
469649ef022SMatthew Knepley   id   = 2;
47045480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp",  label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
471649ef022SMatthew Knepley   id   = 4;
47245480ffeSMatthew G. Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp",   label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
473649ef022SMatthew Knepley 
474649ef022SMatthew Knepley   /*setup exact solution.*/
475649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr);
476649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr);
477649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr);
478649ef022SMatthew Knepley   PetscFunctionReturn(0);
479649ef022SMatthew Knepley }
480649ef022SMatthew Knepley 
481649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
482649ef022SMatthew Knepley {
483649ef022SMatthew Knepley   DM              cdm   = dm;
484649ef022SMatthew Knepley   PetscFE         fe[3];
485649ef022SMatthew Knepley   Parameter      *param;
486649ef022SMatthew Knepley   MPI_Comm        comm;
487*30602db0SMatthew G. Knepley   PetscInt        dim;
488*30602db0SMatthew G. Knepley   PetscBool       simplex;
489649ef022SMatthew Knepley   PetscErrorCode  ierr;
490649ef022SMatthew Knepley 
491649ef022SMatthew Knepley   PetscFunctionBeginUser;
492*30602db0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
493*30602db0SMatthew G. Knepley   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
494649ef022SMatthew Knepley   /* Create finite element */
495649ef022SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
496*30602db0SMatthew G. Knepley   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
497649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
498649ef022SMatthew Knepley 
499*30602db0SMatthew G. Knepley   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
500649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
501649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
502649ef022SMatthew Knepley 
503*30602db0SMatthew G. Knepley   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
504649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
505649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
506649ef022SMatthew Knepley 
507649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
508649ef022SMatthew Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
509649ef022SMatthew Knepley   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
510649ef022SMatthew Knepley   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
511649ef022SMatthew Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
512649ef022SMatthew Knepley   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
513649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
514649ef022SMatthew Knepley   while (cdm) {
515649ef022SMatthew Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
516649ef022SMatthew Knepley     ierr = DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0);CHKERRQ(ierr);
517649ef022SMatthew Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
518649ef022SMatthew Knepley   }
519649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
520649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
521649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr);
522649ef022SMatthew Knepley   PetscFunctionReturn(0);
523649ef022SMatthew Knepley }
524649ef022SMatthew Knepley 
525649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
526649ef022SMatthew Knepley {
527649ef022SMatthew Knepley   Vec              vec;
528649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
529649ef022SMatthew Knepley   PetscErrorCode   ierr;
530649ef022SMatthew Knepley 
531649ef022SMatthew Knepley   PetscFunctionBeginUser;
532649ef022SMatthew Knepley   if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield);
533649ef022SMatthew Knepley   funcs[nfield] = constant;
534649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
535649ef022SMatthew Knepley   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
536649ef022SMatthew Knepley   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
537649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
538649ef022SMatthew Knepley   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
539649ef022SMatthew Knepley   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
540649ef022SMatthew Knepley   ierr = VecDestroy(&vec);CHKERRQ(ierr);
541649ef022SMatthew Knepley   PetscFunctionReturn(0);
542649ef022SMatthew Knepley }
543649ef022SMatthew Knepley 
544649ef022SMatthew Knepley int main(int argc, char **argv)
545649ef022SMatthew Knepley {
546649ef022SMatthew Knepley   SNES            snes;                 /* nonlinear solver */
547649ef022SMatthew Knepley   DM              dm;                   /* problem definition */
548649ef022SMatthew Knepley   Vec             u, r;                 /* solution, residual vectors */
549649ef022SMatthew Knepley   AppCtx          user;                 /* user-defined work context */
550649ef022SMatthew Knepley   PetscErrorCode  ierr;
551649ef022SMatthew Knepley 
552649ef022SMatthew Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
553649ef022SMatthew Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
554649ef022SMatthew Knepley   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
555649ef022SMatthew Knepley   ierr = SetupParameters(&user);CHKERRQ(ierr);
556649ef022SMatthew Knepley   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
557649ef022SMatthew Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
558649ef022SMatthew Knepley   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
559649ef022SMatthew Knepley   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
560649ef022SMatthew Knepley   /* Setup problem */
561649ef022SMatthew Knepley   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
562649ef022SMatthew Knepley   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
563649ef022SMatthew Knepley 
564649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
565649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
566649ef022SMatthew Knepley   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
567649ef022SMatthew Knepley 
568649ef022SMatthew Knepley   ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
569649ef022SMatthew Knepley   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
570649ef022SMatthew Knepley 
571649ef022SMatthew Knepley   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
572649ef022SMatthew Knepley   {
573649ef022SMatthew Knepley     PetscDS          ds;
574649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
575649ef022SMatthew Knepley     void            *ctxs[3];
576649ef022SMatthew Knepley 
577649ef022SMatthew Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
578649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr);
579649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr);
580649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr);
581649ef022SMatthew Knepley     ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
582649ef022SMatthew Knepley     ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr);
583649ef022SMatthew Knepley     ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr);
584649ef022SMatthew Knepley   }
585649ef022SMatthew Knepley   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
586649ef022SMatthew Knepley   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
587649ef022SMatthew Knepley   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
588649ef022SMatthew Knepley 
589649ef022SMatthew Knepley   if (user.showError) {
590649ef022SMatthew Knepley     PetscDS          ds;
591649ef022SMatthew Knepley     Vec              r;
592649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
593649ef022SMatthew Knepley     void            *ctxs[3];
594649ef022SMatthew Knepley 
595649ef022SMatthew Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
596649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr);
597649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr);
598649ef022SMatthew Knepley     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr);
599649ef022SMatthew Knepley     ierr = DMGetGlobalVector(dm, &r);CHKERRQ(ierr);
600649ef022SMatthew Knepley     ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r);CHKERRQ(ierr);
601649ef022SMatthew Knepley     ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr);
602649ef022SMatthew Knepley     ierr = PetscObjectSetName((PetscObject) r, "Solution Error");CHKERRQ(ierr);
603649ef022SMatthew Knepley     ierr = VecViewFromOptions(r, NULL, "-error_vec_view");CHKERRQ(ierr);
604649ef022SMatthew Knepley     ierr = DMRestoreGlobalVector(dm, &r);CHKERRQ(ierr);
605649ef022SMatthew Knepley   }
606649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
607649ef022SMatthew Knepley   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
608649ef022SMatthew Knepley 
609649ef022SMatthew Knepley   ierr = VecDestroy(&u);CHKERRQ(ierr);
610649ef022SMatthew Knepley   ierr = VecDestroy(&r);CHKERRQ(ierr);
611649ef022SMatthew Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
612649ef022SMatthew Knepley   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
613649ef022SMatthew Knepley   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
614649ef022SMatthew Knepley   ierr = PetscFinalize();
615649ef022SMatthew Knepley   return ierr;
616649ef022SMatthew Knepley }
617649ef022SMatthew Knepley 
618649ef022SMatthew Knepley /*TEST
619649ef022SMatthew Knepley 
620649ef022SMatthew Knepley   test:
621649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1
622649ef022SMatthew Knepley     requires: triangle !single
623649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
624649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
625649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
626649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
627649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
628649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
629649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
630649ef022SMatthew Knepley 
631649ef022SMatthew Knepley   test:
632649ef022SMatthew Knepley     # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9]
633649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1_conv
634649ef022SMatthew Knepley     requires: triangle !single
635649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
636649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
637649ef022SMatthew Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \
638649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
639649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
640649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
641649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
642649ef022SMatthew Knepley 
643649ef022SMatthew Knepley   test:
644649ef022SMatthew Knepley     suffix: 2d_tri_p3_p2_p2
645649ef022SMatthew Knepley     requires: triangle !single
646649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
647649ef022SMatthew Knepley       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
648649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
649649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
650649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
651649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
652649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
653649ef022SMatthew Knepley 
654649ef022SMatthew Knepley TEST*/
655