xref: /petsc/src/ts/tutorials/ex76.c (revision 606d57d4c425a15306d39756a53181975473b8fc)
1649ef022SMatthew Knepley static char help[] = "Time-dependent 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 time-dependent 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, du/dt> + <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 <petscts.h>
25649ef022SMatthew Knepley #include <petscds.h>
26649ef022SMatthew Knepley #include <petscbag.h>
27649ef022SMatthew Knepley 
28*606d57d4SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, SOL_TAYLOR_GREEN, NUM_SOL_TYPES} SolType;
29*606d57d4SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "unknown"};
30649ef022SMatthew Knepley 
31649ef022SMatthew Knepley typedef struct {
32649ef022SMatthew Knepley   PetscReal nu;    /* Kinematic viscosity */
33649ef022SMatthew Knepley   PetscReal alpha; /* Thermal diffusivity */
34649ef022SMatthew Knepley   PetscReal T_in;  /* Inlet temperature*/
35649ef022SMatthew Knepley } Parameter;
36649ef022SMatthew Knepley 
37649ef022SMatthew Knepley typedef struct {
38649ef022SMatthew Knepley   /* Problem definition */
39649ef022SMatthew Knepley   PetscBag bag;     /* Holds problem parameters */
40649ef022SMatthew Knepley   SolType  solType; /* MMS solution type */
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 = t + x^2 + y^2
62649ef022SMatthew Knepley     v = t + 2x^2 - 2xy
63649ef022SMatthew Knepley     p = x + y - 1
64649ef022SMatthew Knepley     T = t + x + y
65649ef022SMatthew Knepley     f = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 -4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 -4\nu + 2>
66649ef022SMatthew Knepley     Q = 1 + 2t + 3x^2 - 2xy + y^2
67649ef022SMatthew Knepley 
68649ef022SMatthew Knepley   so that
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley     \nabla \cdot u = 2x - 2x = 0
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
73649ef022SMatthew Knepley     = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1>
74649ef022SMatthew Knepley     = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> + <-4 \nu + 2, -4\nu + 2>
75649ef022SMatthew Knepley     = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 - 4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 - 4\nu + 2>
76649ef022SMatthew Knepley 
77649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
78649ef022SMatthew Knepley     = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0
79649ef022SMatthew Knepley     = 1 + 2t + 3x^2 - 2xy + y^2
80649ef022SMatthew Knepley */
81649ef022SMatthew Knepley 
82649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
83649ef022SMatthew Knepley {
84649ef022SMatthew Knepley   u[0] = time + X[0]*X[0] + X[1]*X[1];
85649ef022SMatthew Knepley   u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
86649ef022SMatthew Knepley   return 0;
87649ef022SMatthew Knepley }
88649ef022SMatthew Knepley static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
89649ef022SMatthew Knepley {
90649ef022SMatthew Knepley   u[0] = 1.0;
91649ef022SMatthew Knepley   u[1] = 1.0;
92649ef022SMatthew Knepley   return 0;
93649ef022SMatthew Knepley }
94649ef022SMatthew Knepley 
95649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
96649ef022SMatthew Knepley {
97649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
98649ef022SMatthew Knepley   return 0;
99649ef022SMatthew Knepley }
100649ef022SMatthew Knepley 
101649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
102649ef022SMatthew Knepley {
103649ef022SMatthew Knepley   T[0] = time + X[0] + X[1];
104649ef022SMatthew Knepley   return 0;
105649ef022SMatthew Knepley }
106649ef022SMatthew Knepley static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
107649ef022SMatthew Knepley {
108649ef022SMatthew Knepley   T[0] = 1.0;
109649ef022SMatthew Knepley   return 0;
110649ef022SMatthew Knepley }
111649ef022SMatthew Knepley 
112649ef022SMatthew Knepley /* f0_v = du/dt - f */
113649ef022SMatthew Knepley static void f0_quadratic_v(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   const PetscReal nu = PetscRealPart(constants[0]);
119649ef022SMatthew Knepley   PetscInt        Nc = dim;
120649ef022SMatthew Knepley   PetscInt        c, d;
121649ef022SMatthew Knepley 
122649ef022SMatthew Knepley   for (d = 0; d<dim; ++d) f0[d] = u_t[uOff[0]+d];
123649ef022SMatthew Knepley 
124649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
125649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
126649ef022SMatthew Knepley   }
127649ef022SMatthew Knepley   f0[0] -= (t*(2*X[0] + 2*X[1]) + 2*X[0]*X[0]*X[0] + 4*X[0]*X[0]*X[1] - 2*X[0]*X[1]*X[1] - 4.0*nu + 2);
128649ef022SMatthew Knepley   f0[1] -= (t*(2*X[0] - 2*X[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 + 2);
129649ef022SMatthew Knepley }
130649ef022SMatthew Knepley 
131649ef022SMatthew Knepley /* f0_w = dT/dt + u.grad(T) - Q */
132649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
136649ef022SMatthew Knepley {
137649ef022SMatthew Knepley   PetscInt d;
138649ef022SMatthew Knepley   f0[0] = 0;
139649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
140649ef022SMatthew Knepley   f0[0] += u_t[uOff[2]] - (2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]);
141649ef022SMatthew Knepley }
142649ef022SMatthew Knepley 
143649ef022SMatthew Knepley /*
144649ef022SMatthew Knepley   CASE: cubic
145649ef022SMatthew Knepley   In 2D we use exact solution:
146649ef022SMatthew Knepley 
147649ef022SMatthew Knepley     u = t + x^3 + y^3
148649ef022SMatthew Knepley     v = t + 2x^3 - 3x^2y
149649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
150649ef022SMatthew Knepley     T = t + 1/2 x^2 + 1/2 y^2
151649ef022SMatthew Knepley     f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1,
152649ef022SMatthew Knepley           t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1>
153649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1
154649ef022SMatthew Knepley 
155649ef022SMatthew Knepley   so that
156649ef022SMatthew Knepley 
157649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
158649ef022SMatthew Knepley 
159649ef022SMatthew Knepley   du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f
160649ef022SMatthew Knepley   = <1,1> + <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1, t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1>  = 0
161649ef022SMatthew Knepley 
162649ef022SMatthew Knepley   dT/dt + u \cdot \nabla T - \alpha \Delta T - Q = 1 + (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2*\alpha +1)   = 0
163649ef022SMatthew Knepley */
164649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
165649ef022SMatthew Knepley {
166649ef022SMatthew Knepley   u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
167649ef022SMatthew Knepley   u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
168649ef022SMatthew Knepley   return 0;
169649ef022SMatthew Knepley }
170649ef022SMatthew Knepley static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
171649ef022SMatthew Knepley {
172649ef022SMatthew Knepley   u[0] = 1.0;
173649ef022SMatthew Knepley   u[1] = 1.0;
174649ef022SMatthew Knepley   return 0;
175649ef022SMatthew Knepley }
176649ef022SMatthew Knepley 
177649ef022SMatthew Knepley static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
178649ef022SMatthew Knepley {
179649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
180649ef022SMatthew Knepley   return 0;
181649ef022SMatthew Knepley }
182649ef022SMatthew Knepley 
183649ef022SMatthew Knepley static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
184649ef022SMatthew Knepley {
185649ef022SMatthew Knepley   T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
186649ef022SMatthew Knepley   return 0;
187649ef022SMatthew Knepley }
188649ef022SMatthew Knepley static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
189649ef022SMatthew Knepley {
190649ef022SMatthew Knepley   T[0] = 1.0;
191649ef022SMatthew Knepley   return 0;
192649ef022SMatthew Knepley }
193649ef022SMatthew Knepley 
194649ef022SMatthew Knepley 
195649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
196649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
197649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
198649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
199649ef022SMatthew Knepley {
200649ef022SMatthew Knepley   PetscInt                   c, d;
201649ef022SMatthew Knepley   PetscInt                   Nc = dim;
202649ef022SMatthew Knepley   const PetscReal            nu = PetscRealPart(constants[0]);
203649ef022SMatthew Knepley 
204649ef022SMatthew Knepley   for (d=0; d<dim; ++d) f0[d] = u_t[uOff[0]+d];
205649ef022SMatthew Knepley 
206649ef022SMatthew Knepley   for (c=0; c<Nc; ++c) {
207649ef022SMatthew Knepley     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
208649ef022SMatthew Knepley   }
209649ef022SMatthew Knepley   f0[0] -= (t*(3*X[0]*X[0] + 3*X[1]*X[1]) + 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] + 1);
210649ef022SMatthew Knepley   f0[1] -= (t*(3*X[0]*X[0] - 6*X[0]*X[1]) + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1] + 1);
211649ef022SMatthew Knepley }
212649ef022SMatthew Knepley 
213649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
214649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
215649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
216649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
217649ef022SMatthew Knepley {
218649ef022SMatthew Knepley   PetscInt              d;
219649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
220649ef022SMatthew Knepley 
221649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
222649ef022SMatthew Knepley   f0[0] += u_t[uOff[2]] - (X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] + X[0]*t + X[1]*t - 2.0*alpha + 1);
223649ef022SMatthew Knepley }
224649ef022SMatthew Knepley 
225649ef022SMatthew Knepley /*
226649ef022SMatthew Knepley   CASE: cubic-trigonometric
227649ef022SMatthew Knepley   In 2D we use exact solution:
228649ef022SMatthew Knepley 
229649ef022SMatthew Knepley     u = beta cos t + x^3 + y^3
230649ef022SMatthew Knepley     v = beta sin t + 2x^3 - 3x^2y
231649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
232649ef022SMatthew Knepley     T = 20 cos t + 1/2 x^2 + 1/2 y^2
233649ef022SMatthew Knepley     f = < beta cos t 3x^2         + beta sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y)  + 3x,
234649ef022SMatthew Knepley           beta cos t (6x^2 - 6xy) - beta sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu(12x - 6y) + 3y>
235649ef022SMatthew Knepley     Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha
236649ef022SMatthew Knepley 
237649ef022SMatthew Knepley   so that
238649ef022SMatthew Knepley 
239649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
240649ef022SMatthew Knepley 
241649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
242649ef022SMatthew Knepley     = <-sin t, cos t> + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> <<3x^2, 6x^2 - 6xy>, <3y^2, -3x^2>> - \nu <6x + 6y, 12x - 6y> + <3x, 3y>
243649ef022SMatthew Knepley     = <-sin t, cos t> + <cos t 3x^2 + 3x^5 + 3x^2y^3 + sin t 3y^2 + 6x^3y^2 - 9x^2y^3, cos t (6x^2 - 6xy) + 6x^5 - 6x^4y + 6x^2y^3 - 6xy^4 + sin t (-3x^2) - 6x^5 + 9x^4y> - \nu <6x + 6y, 12x - 6y> + <3x, 3y>
244649ef022SMatthew Knepley     = <cos t (3x^2)       + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y)  + 3x,
245649ef022SMatthew Knepley        cos t (6x^2 - 6xy) - sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu (12x - 6y) + 3y>
246649ef022SMatthew Knepley 
247649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
248649ef022SMatthew Knepley     = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha
249649ef022SMatthew Knepley     = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha
250649ef022SMatthew Knepley     = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha)
251649ef022SMatthew Knepley */
252649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
253649ef022SMatthew Knepley {
254649ef022SMatthew Knepley   u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
255649ef022SMatthew Knepley   u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
256649ef022SMatthew Knepley   return 0;
257649ef022SMatthew Knepley }
258649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
259649ef022SMatthew Knepley {
260649ef022SMatthew Knepley   u[0] = -100.*PetscSinReal(time);
261649ef022SMatthew Knepley   u[1] =  100.*PetscCosReal(time);
262649ef022SMatthew Knepley   return 0;
263649ef022SMatthew Knepley }
264649ef022SMatthew Knepley 
265649ef022SMatthew Knepley static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
266649ef022SMatthew Knepley {
267649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
268649ef022SMatthew Knepley   return 0;
269649ef022SMatthew Knepley }
270649ef022SMatthew Knepley 
271649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
272649ef022SMatthew Knepley {
273649ef022SMatthew Knepley   T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
274649ef022SMatthew Knepley   return 0;
275649ef022SMatthew Knepley }
276649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
277649ef022SMatthew Knepley {
278649ef022SMatthew Knepley   T[0] = -100.*PetscSinReal(time);
279649ef022SMatthew Knepley   return 0;
280649ef022SMatthew Knepley }
281649ef022SMatthew Knepley 
282649ef022SMatthew Knepley static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
283649ef022SMatthew Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
284649ef022SMatthew Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
285649ef022SMatthew Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
286649ef022SMatthew Knepley {
287649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
288649ef022SMatthew Knepley   PetscInt        Nc = dim;
289649ef022SMatthew Knepley   PetscInt        c, d;
290649ef022SMatthew Knepley 
291649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d];
292649ef022SMatthew Knepley 
293649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
294649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
295649ef022SMatthew Knepley   }
296649ef022SMatthew Knepley   f0[0] -= 100.*PetscCosReal(t)*(3*X[0]*X[0])               + 100.*PetscSinReal(t)*(3*X[1]*X[1] - 1.) + 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];
297649ef022SMatthew Knepley   f0[1] -= 100.*PetscCosReal(t)*(6*X[0]*X[0] - 6*X[0]*X[1]) - 100.*PetscSinReal(t)*(3*X[0]*X[0])      + 3*X[0]*X[0]*X[0]*X[0]*X[1] + 6*X[0]*X[0]*X[1]*X[1]*X[1] - 6*X[0]*X[1]*X[1]*X[1]*X[1] - (12*X[0] - 6*X[1])*nu + 3*X[1];
298649ef022SMatthew Knepley }
299649ef022SMatthew Knepley 
300649ef022SMatthew Knepley static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
301649ef022SMatthew Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
302649ef022SMatthew Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
303649ef022SMatthew Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
304649ef022SMatthew Knepley {
305649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
306649ef022SMatthew Knepley   PetscInt        d;
307649ef022SMatthew Knepley 
308649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
309649ef022SMatthew Knepley   f0[0] += u_t[uOff[2]] - (100.*PetscCosReal(t)*X[0] + 100.*PetscSinReal(t)*(X[1] - 1.) + X[0]*X[0]*X[0]*X[0] + 2.0*X[0]*X[0]*X[0]*X[1] - 3.0*X[0]*X[0]*X[1]*X[1] + X[0]*X[1]*X[1]*X[1] - 2.0*alpha);
310649ef022SMatthew Knepley }
311649ef022SMatthew Knepley 
312*606d57d4SMatthew G. Knepley /*
313*606d57d4SMatthew G. Knepley   CASE: taylor-green vortex
314*606d57d4SMatthew G. Knepley   In 2D we use exact solution:
315*606d57d4SMatthew G. Knepley 
316*606d57d4SMatthew G. Knepley     u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)
317*606d57d4SMatthew G. Knepley     v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)
318*606d57d4SMatthew G. Knepley     p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t)
319*606d57d4SMatthew G. Knepley     T = t + x + y
320*606d57d4SMatthew G. Knepley     f = <\nu \pi^2 exp(-2\nu \pi^2 t) cos(\pi(x-t)) sin(\pi(y-t)), -\nu \pi^2 exp(-2\nu \pi^2 t) sin(\pi(x-t)) cos(\pi(y-t))  >
321*606d57d4SMatthew G. Knepley     Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t)
322*606d57d4SMatthew G. Knepley 
323*606d57d4SMatthew G. Knepley   so that
324*606d57d4SMatthew G. Knepley 
325*606d57d4SMatthew G. Knepley   \nabla \cdot u = \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) - \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t) = 0
326*606d57d4SMatthew G. Knepley 
327*606d57d4SMatthew G. Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
328*606d57d4SMatthew G. Knepley     = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t),
329*606d57d4SMatthew G. Knepley         \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)>
330*606d57d4SMatthew G. Knepley     + < \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
331*606d57d4SMatthew G. Knepley         \pi (1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
332*606d57d4SMatthew G. Knepley     + <-\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
333*606d57d4SMatthew G. Knepley        -\pi (1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)) sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
334*606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
335*606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
336*606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
337*606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
338*606d57d4SMatthew G. Knepley     = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi cos(\pi(x - t)) sin(\pi(y - t))) exp(-2 \pi^2 \nu t),
339*606d57d4SMatthew G. Knepley         \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t)) - 2\pi sin(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)>
340*606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
341*606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
342*606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
343*606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
344*606d57d4SMatthew G. Knepley     + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
345*606d57d4SMatthew G. Knepley        -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
346*606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
347*606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
348*606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
349*606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
350*606d57d4SMatthew G. Knepley     = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t),
351*606d57d4SMatthew G. Knepley         \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)>
352*606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
353*606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
354*606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
355*606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
356*606d57d4SMatthew G. Knepley     = < \pi cos(\pi(x - t)) cos(\pi(y - t)),
357*606d57d4SMatthew G. Knepley         \pi sin(\pi(x - t)) sin(\pi(y - t))>
358*606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)),
359*606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0
360*606d57d4SMatthew G. Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
361*606d57d4SMatthew G. Knepley     = 1 + u \cdot <1, 1> - 0
362*606d57d4SMatthew G. Knepley     = 1 + u + v
363*606d57d4SMatthew G. Knepley */
364*606d57d4SMatthew G. Knepley 
365*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
366*606d57d4SMatthew G. Knepley {
367*606d57d4SMatthew G. Knepley   u[0] = 1 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
368*606d57d4SMatthew G. Knepley   u[1] = 1 + PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
369*606d57d4SMatthew G. Knepley   return 0;
370*606d57d4SMatthew G. Knepley }
371*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
372*606d57d4SMatthew G. Knepley {
373*606d57d4SMatthew G. Knepley   u[0] = -PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
374*606d57d4SMatthew G. Knepley                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
375*606d57d4SMatthew G. Knepley                   - 2*PETSC_PI*PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
376*606d57d4SMatthew G. Knepley   u[1] =  PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
377*606d57d4SMatthew G. Knepley                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
378*606d57d4SMatthew G. Knepley                   - 2*PETSC_PI*PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
379*606d57d4SMatthew G. Knepley   return 0;
380*606d57d4SMatthew G. Knepley }
381*606d57d4SMatthew G. Knepley 
382*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
383*606d57d4SMatthew G. Knepley {
384*606d57d4SMatthew G. Knepley   p[0] = -0.25*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time)))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time);
385*606d57d4SMatthew G. Knepley   return 0;
386*606d57d4SMatthew G. Knepley }
387*606d57d4SMatthew G. Knepley 
388*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
389*606d57d4SMatthew G. Knepley {
390*606d57d4SMatthew G. Knepley   p[0] = PETSC_PI*(0.5*(PetscSinReal(2*PETSC_PI*(X[0]-time)) + PetscSinReal(2*PETSC_PI*(X[1]-time)))
391*606d57d4SMatthew G. Knepley                  + PETSC_PI*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time))))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time);
392*606d57d4SMatthew G. Knepley   return 0;
393*606d57d4SMatthew G. Knepley }
394*606d57d4SMatthew G. Knepley 
395*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
396*606d57d4SMatthew G. Knepley {
397*606d57d4SMatthew G. Knepley   T[0] = time + X[0] + X[1];
398*606d57d4SMatthew G. Knepley   return 0;
399*606d57d4SMatthew G. Knepley }
400*606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
401*606d57d4SMatthew G. Knepley {
402*606d57d4SMatthew G. Knepley   T[0] = 1.0;
403*606d57d4SMatthew G. Knepley   return 0;
404*606d57d4SMatthew G. Knepley }
405*606d57d4SMatthew G. Knepley 
406*606d57d4SMatthew G. Knepley static void f0_taylor_green_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
407*606d57d4SMatthew G. Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
408*606d57d4SMatthew G. Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
409*606d57d4SMatthew G. Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
410*606d57d4SMatthew G. Knepley {
411*606d57d4SMatthew G. Knepley   PetscInt        Nc = dim;
412*606d57d4SMatthew G. Knepley   PetscInt        c, d;
413*606d57d4SMatthew G. Knepley 
414*606d57d4SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d];
415*606d57d4SMatthew G. Knepley 
416*606d57d4SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
417*606d57d4SMatthew G. Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
418*606d57d4SMatthew G. Knepley   }
419*606d57d4SMatthew G. Knepley }
420*606d57d4SMatthew G. Knepley 
421*606d57d4SMatthew G. Knepley static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
422*606d57d4SMatthew G. Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
423*606d57d4SMatthew G. Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
424*606d57d4SMatthew G. Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
425*606d57d4SMatthew G. Knepley {
426*606d57d4SMatthew G. Knepley   PetscScalar vel[2];
427*606d57d4SMatthew G. Knepley   PetscInt    d;
428*606d57d4SMatthew G. Knepley 
429*606d57d4SMatthew G. Knepley   taylor_green_u(dim, t, X, Nf, vel, NULL);
430*606d57d4SMatthew G. Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
431*606d57d4SMatthew G. Knepley   f0[0] += u_t[uOff[2]] - (1.0 + vel[0] + vel[1]);
432*606d57d4SMatthew G. Knepley }
433*606d57d4SMatthew G. Knepley 
434649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
435649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
436649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
437649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
438649ef022SMatthew Knepley {
439649ef022SMatthew Knepley   PetscInt d;
440649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
441649ef022SMatthew Knepley }
442649ef022SMatthew Knepley 
443649ef022SMatthew Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
444649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
445649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
446649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
447649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
448649ef022SMatthew Knepley {
449649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
450649ef022SMatthew Knepley   const PetscInt    Nc = dim;
451649ef022SMatthew Knepley   PetscInt        c, d;
452649ef022SMatthew Knepley 
453649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
454649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
455649ef022SMatthew Knepley       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
456649ef022SMatthew Knepley     }
457649ef022SMatthew Knepley     f1[c*dim+c] -= u[uOff[1]];
458649ef022SMatthew Knepley   }
459649ef022SMatthew Knepley }
460649ef022SMatthew Knepley 
461649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
462649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
463649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
464649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
465649ef022SMatthew Knepley {
466649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
467649ef022SMatthew Knepley   PetscInt d;
468649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
469649ef022SMatthew Knepley }
470649ef022SMatthew Knepley 
471649ef022SMatthew Knepley /*Jacobians*/
472649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
473649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
474649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
475649ef022SMatthew Knepley                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
476649ef022SMatthew Knepley {
477649ef022SMatthew Knepley   PetscInt d;
478649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
479649ef022SMatthew Knepley }
480649ef022SMatthew Knepley 
481649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
482649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
483649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
484649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
485649ef022SMatthew Knepley {
486649ef022SMatthew Knepley   PetscInt c, d;
487649ef022SMatthew Knepley   const PetscInt  Nc = dim;
488649ef022SMatthew Knepley 
489649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift;
490649ef022SMatthew Knepley 
491649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
492649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
493649ef022SMatthew Knepley       g0[c*Nc+d] += u_x[ c*Nc+d];
494649ef022SMatthew Knepley     }
495649ef022SMatthew Knepley   }
496649ef022SMatthew Knepley }
497649ef022SMatthew Knepley 
498649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
499649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
500649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
501649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
502649ef022SMatthew Knepley {
503649ef022SMatthew Knepley   PetscInt NcI = dim;
504649ef022SMatthew Knepley   PetscInt NcJ = dim;
505649ef022SMatthew Knepley   PetscInt c, d, e;
506649ef022SMatthew Knepley 
507649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
508649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
509649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
510649ef022SMatthew Knepley         if (c == d) {
511649ef022SMatthew Knepley           g1[(c*NcJ+d)*dim+e] += u[e];
512649ef022SMatthew Knepley         }
513649ef022SMatthew Knepley       }
514649ef022SMatthew Knepley     }
515649ef022SMatthew Knepley   }
516649ef022SMatthew Knepley }
517649ef022SMatthew Knepley 
518649ef022SMatthew Knepley 
519649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
520649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
521649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
522649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
523649ef022SMatthew Knepley {
524649ef022SMatthew Knepley   PetscInt d;
525649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
526649ef022SMatthew Knepley }
527649ef022SMatthew Knepley 
528649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
529649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
530649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
531649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
532649ef022SMatthew Knepley {
533649ef022SMatthew Knepley    const PetscReal nu = PetscRealPart(constants[0]);
534649ef022SMatthew Knepley    const PetscInt  Nc = dim;
535649ef022SMatthew Knepley    PetscInt        c, d;
536649ef022SMatthew Knepley 
537649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
538649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
539*606d57d4SMatthew G. Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += nu;
540*606d57d4SMatthew G. Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += nu;
541649ef022SMatthew Knepley     }
542649ef022SMatthew Knepley   }
543649ef022SMatthew Knepley }
544649ef022SMatthew Knepley 
545649ef022SMatthew Knepley static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
546649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
547649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
548649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
549649ef022SMatthew Knepley {
550649ef022SMatthew Knepley   PetscInt d;
551649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_tShift;
552649ef022SMatthew Knepley }
553649ef022SMatthew Knepley 
554649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
555649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
556649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
557649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
558649ef022SMatthew Knepley {
559649ef022SMatthew Knepley   PetscInt d;
560649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
561649ef022SMatthew Knepley }
562649ef022SMatthew Knepley 
563649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
564649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
565649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
566649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
567649ef022SMatthew Knepley {
568649ef022SMatthew Knepley   PetscInt d;
569649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
570649ef022SMatthew Knepley }
571649ef022SMatthew Knepley 
572649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
573649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
574649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
575649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
576649ef022SMatthew Knepley {
577649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
578649ef022SMatthew Knepley   PetscInt               d;
579649ef022SMatthew Knepley 
580649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
581649ef022SMatthew Knepley }
582649ef022SMatthew Knepley 
583649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
584649ef022SMatthew Knepley {
585649ef022SMatthew Knepley   PetscInt       sol;
586649ef022SMatthew Knepley   PetscErrorCode ierr;
587649ef022SMatthew Knepley 
588649ef022SMatthew Knepley 
589649ef022SMatthew Knepley   PetscFunctionBeginUser;
590649ef022SMatthew Knepley   options->solType = SOL_QUADRATIC;
591649ef022SMatthew Knepley 
592649ef022SMatthew Knepley   ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr);
593649ef022SMatthew Knepley   sol = options->solType;
594649ef022SMatthew Knepley   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
595649ef022SMatthew Knepley   options->solType = (SolType) sol;
596649ef022SMatthew Knepley   ierr = PetscOptionsEnd();
597649ef022SMatthew Knepley   PetscFunctionReturn(0);
598649ef022SMatthew Knepley }
599649ef022SMatthew Knepley 
600649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user)
601649ef022SMatthew Knepley {
602649ef022SMatthew Knepley   PetscBag       bag;
603649ef022SMatthew Knepley   Parameter     *p;
604649ef022SMatthew Knepley   PetscErrorCode ierr;
605649ef022SMatthew Knepley 
606649ef022SMatthew Knepley   PetscFunctionBeginUser;
607649ef022SMatthew Knepley   /* setup PETSc parameter bag */
608649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
609649ef022SMatthew Knepley   ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr);
610649ef022SMatthew Knepley   bag  = user->bag;
611649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
612649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr);
613649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->T_in,  1.0, "T_in",  "Inlet temperature");CHKERRQ(ierr);
614649ef022SMatthew Knepley   PetscFunctionReturn(0);
615649ef022SMatthew Knepley }
616649ef022SMatthew Knepley 
617649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
618649ef022SMatthew Knepley {
619649ef022SMatthew Knepley   PetscErrorCode ierr;
620649ef022SMatthew Knepley 
621649ef022SMatthew Knepley   PetscFunctionBeginUser;
622649ef022SMatthew Knepley   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
623649ef022SMatthew Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
624649ef022SMatthew Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
625649ef022SMatthew Knepley   PetscFunctionReturn(0);
626649ef022SMatthew Knepley }
627649ef022SMatthew Knepley 
628649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
629649ef022SMatthew Knepley {
630649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
631649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
632649ef022SMatthew Knepley   PetscDS          prob;
633649ef022SMatthew Knepley   Parameter       *ctx;
634649ef022SMatthew Knepley   PetscInt         id;
635649ef022SMatthew Knepley   PetscErrorCode   ierr;
636649ef022SMatthew Knepley 
637649ef022SMatthew Knepley   PetscFunctionBeginUser;
638649ef022SMatthew Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
639649ef022SMatthew Knepley   switch(user->solType){
640649ef022SMatthew Knepley   case SOL_QUADRATIC:
641649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr);
642649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr);
643649ef022SMatthew Knepley 
644649ef022SMatthew Knepley     exactFuncs[0]   = quadratic_u;
645649ef022SMatthew Knepley     exactFuncs[1]   = quadratic_p;
646649ef022SMatthew Knepley     exactFuncs[2]   = quadratic_T;
647649ef022SMatthew Knepley     exactFuncs_t[0] = quadratic_u_t;
648649ef022SMatthew Knepley     exactFuncs_t[1] = NULL;
649649ef022SMatthew Knepley     exactFuncs_t[2] = quadratic_T_t;
650649ef022SMatthew Knepley     break;
651649ef022SMatthew Knepley   case SOL_CUBIC:
652649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr);
653649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr);
654649ef022SMatthew Knepley 
655649ef022SMatthew Knepley     exactFuncs[0]   = cubic_u;
656649ef022SMatthew Knepley     exactFuncs[1]   = cubic_p;
657649ef022SMatthew Knepley     exactFuncs[2]   = cubic_T;
658649ef022SMatthew Knepley     exactFuncs_t[0] = cubic_u_t;
659649ef022SMatthew Knepley     exactFuncs_t[1] = NULL;
660649ef022SMatthew Knepley     exactFuncs_t[2] = cubic_T_t;
661649ef022SMatthew Knepley     break;
662649ef022SMatthew Knepley   case SOL_CUBIC_TRIG:
663649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 0, f0_cubic_trig_v, f1_v);CHKERRQ(ierr);
664649ef022SMatthew Knepley     ierr = PetscDSSetResidual(prob, 2, f0_cubic_trig_w, f1_w);CHKERRQ(ierr);
665649ef022SMatthew Knepley 
666649ef022SMatthew Knepley     exactFuncs[0]   = cubic_trig_u;
667649ef022SMatthew Knepley     exactFuncs[1]   = cubic_trig_p;
668649ef022SMatthew Knepley     exactFuncs[2]   = cubic_trig_T;
669649ef022SMatthew Knepley     exactFuncs_t[0] = cubic_trig_u_t;
670649ef022SMatthew Knepley     exactFuncs_t[1] = NULL;
671649ef022SMatthew Knepley     exactFuncs_t[2] = cubic_trig_T_t;
672649ef022SMatthew Knepley     break;
673*606d57d4SMatthew G. Knepley   case SOL_TAYLOR_GREEN:
674*606d57d4SMatthew G. Knepley     ierr = PetscDSSetResidual(prob, 0, f0_taylor_green_v, f1_v);CHKERRQ(ierr);
675*606d57d4SMatthew G. Knepley     ierr = PetscDSSetResidual(prob, 2, f0_taylor_green_w, f1_w);CHKERRQ(ierr);
676*606d57d4SMatthew G. Knepley 
677*606d57d4SMatthew G. Knepley     exactFuncs[0]   = taylor_green_u;
678*606d57d4SMatthew G. Knepley     exactFuncs[1]   = taylor_green_p;
679*606d57d4SMatthew G. Knepley     exactFuncs[2]   = taylor_green_T;
680*606d57d4SMatthew G. Knepley     exactFuncs_t[0] = taylor_green_u_t;
681*606d57d4SMatthew G. Knepley     exactFuncs_t[1] = taylor_green_p_t;
682*606d57d4SMatthew G. Knepley     exactFuncs_t[2] = taylor_green_T_t;
683*606d57d4SMatthew G. Knepley     break;
684649ef022SMatthew Knepley    default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
685649ef022SMatthew Knepley   }
686649ef022SMatthew Knepley 
687649ef022SMatthew Knepley   ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr);
688649ef022SMatthew Knepley 
689649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu);CHKERRQ(ierr);
690649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL);CHKERRQ(ierr);
691649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL);CHKERRQ(ierr);
692649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL);CHKERRQ(ierr);
693649ef022SMatthew Knepley   ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
694649ef022SMatthew Knepley   /* Setup constants */
695649ef022SMatthew Knepley   {
696649ef022SMatthew Knepley     Parameter  *param;
697649ef022SMatthew Knepley     PetscScalar constants[3];
698649ef022SMatthew Knepley 
699649ef022SMatthew Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
700649ef022SMatthew Knepley 
701649ef022SMatthew Knepley     constants[0] = param->nu;
702649ef022SMatthew Knepley     constants[1] = param->alpha;
703649ef022SMatthew Knepley     constants[2] = param->T_in;
704649ef022SMatthew Knepley     ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr);
705649ef022SMatthew Knepley   }
706649ef022SMatthew Knepley   /* Setup Boundary Conditions */
707649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
708649ef022SMatthew Knepley   id   = 3;
709649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity",    "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
710649ef022SMatthew Knepley   id   = 1;
711649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
712649ef022SMatthew Knepley   id   = 2;
713649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity",  "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
714649ef022SMatthew Knepley   id   = 4;
715649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity",   "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], (void (*)(void)) exactFuncs_t[0], 1, &id, ctx);CHKERRQ(ierr);
716649ef022SMatthew Knepley   id   = 3;
717649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp",    "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
718649ef022SMatthew Knepley   id   = 1;
719649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
720649ef022SMatthew Knepley   id   = 2;
721649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp",  "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
722649ef022SMatthew Knepley   id   = 4;
723649ef022SMatthew Knepley   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp",   "marker", 2, 0, NULL, (void (*)(void)) exactFuncs[2], (void (*)(void)) exactFuncs_t[2], 1, &id, ctx);CHKERRQ(ierr);
724649ef022SMatthew Knepley 
725649ef022SMatthew Knepley   /*setup exact solution.*/
726649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr);
727649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr);
728649ef022SMatthew Knepley   ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr);
729649ef022SMatthew Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr);
730649ef022SMatthew Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr);
731649ef022SMatthew Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr);
732649ef022SMatthew Knepley   PetscFunctionReturn(0);
733649ef022SMatthew Knepley }
734649ef022SMatthew Knepley 
735649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
736649ef022SMatthew Knepley {
737649ef022SMatthew Knepley   DM              cdm   = dm;
738649ef022SMatthew Knepley   PetscFE         fe[3];
739649ef022SMatthew Knepley   Parameter      *param;
740649ef022SMatthew Knepley   MPI_Comm        comm;
741649ef022SMatthew Knepley   DMPolytopeType  ct;
742649ef022SMatthew Knepley   PetscInt        dim, cStart;
743649ef022SMatthew Knepley   PetscBool       simplex;
744649ef022SMatthew Knepley   PetscErrorCode  ierr;
745649ef022SMatthew Knepley 
746649ef022SMatthew Knepley   PetscFunctionBeginUser;
747649ef022SMatthew Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
748649ef022SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
749649ef022SMatthew Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
750649ef022SMatthew Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
751649ef022SMatthew Knepley   /* Create finite element */
752649ef022SMatthew Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
753649ef022SMatthew Knepley   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
754649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
755649ef022SMatthew Knepley 
756649ef022SMatthew Knepley   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
757649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
758649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
759649ef022SMatthew Knepley 
760649ef022SMatthew Knepley   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
761649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
762649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
763649ef022SMatthew Knepley 
764649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
765649ef022SMatthew Knepley   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
766649ef022SMatthew Knepley   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
767649ef022SMatthew Knepley   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
768649ef022SMatthew Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
769649ef022SMatthew Knepley   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
770649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
771649ef022SMatthew Knepley   while (cdm) {
772649ef022SMatthew Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
773649ef022SMatthew Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
774649ef022SMatthew Knepley   }
775649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
776649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
777649ef022SMatthew Knepley   ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr);
778649ef022SMatthew Knepley 
779649ef022SMatthew Knepley   {
780649ef022SMatthew Knepley     PetscObject  pressure;
781649ef022SMatthew Knepley     MatNullSpace nullspacePres;
782649ef022SMatthew Knepley 
783649ef022SMatthew Knepley     ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr);
784649ef022SMatthew Knepley     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
785649ef022SMatthew Knepley     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
786649ef022SMatthew Knepley     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
787649ef022SMatthew Knepley   }
788649ef022SMatthew Knepley 
789649ef022SMatthew Knepley   PetscFunctionReturn(0);
790649ef022SMatthew Knepley }
791649ef022SMatthew Knepley 
792649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
793649ef022SMatthew Knepley {
794649ef022SMatthew Knepley   Vec              vec;
795649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
796649ef022SMatthew Knepley   PetscErrorCode   ierr;
797649ef022SMatthew Knepley 
798649ef022SMatthew Knepley   PetscFunctionBeginUser;
799649ef022SMatthew Knepley   if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield);
800649ef022SMatthew Knepley   funcs[nfield] = constant;
801649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
802649ef022SMatthew Knepley   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
803649ef022SMatthew Knepley   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
804649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
805649ef022SMatthew Knepley   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
806649ef022SMatthew Knepley   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
807649ef022SMatthew Knepley   ierr = VecDestroy(&vec);CHKERRQ(ierr);
808649ef022SMatthew Knepley   PetscFunctionReturn(0);
809649ef022SMatthew Knepley }
810649ef022SMatthew Knepley 
811649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u)
812649ef022SMatthew Knepley {
813649ef022SMatthew Knepley   DM             dm;
814649ef022SMatthew Knepley   MatNullSpace   nullsp;
815649ef022SMatthew Knepley   PetscErrorCode ierr;
816649ef022SMatthew Knepley 
817649ef022SMatthew Knepley   PetscFunctionBegin;
818649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
819649ef022SMatthew Knepley   ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr);
820649ef022SMatthew Knepley   ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr);
821649ef022SMatthew Knepley   ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
822649ef022SMatthew Knepley   PetscFunctionReturn(0);
823649ef022SMatthew Knepley }
824649ef022SMatthew Knepley 
825649ef022SMatthew Knepley /* Make the discrete pressure discretely divergence free */
826649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts)
827649ef022SMatthew Knepley {
828649ef022SMatthew Knepley   Vec            u;
829649ef022SMatthew Knepley   PetscErrorCode ierr;
830649ef022SMatthew Knepley 
831649ef022SMatthew Knepley   PetscFunctionBegin;
832649ef022SMatthew Knepley   ierr = TSGetSolution(ts, &u);CHKERRQ(ierr);
833649ef022SMatthew Knepley   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
834649ef022SMatthew Knepley   PetscFunctionReturn(0);
835649ef022SMatthew Knepley }
836649ef022SMatthew Knepley 
837649ef022SMatthew Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u)
838649ef022SMatthew Knepley {
839649ef022SMatthew Knepley   DM             dm;
840649ef022SMatthew Knepley   PetscReal      t;
841649ef022SMatthew Knepley   PetscErrorCode ierr;
842649ef022SMatthew Knepley 
843649ef022SMatthew Knepley   PetscFunctionBegin;
844649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
845649ef022SMatthew Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
846649ef022SMatthew Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
847649ef022SMatthew Knepley   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
848649ef022SMatthew Knepley   PetscFunctionReturn(0);
849649ef022SMatthew Knepley }
850649ef022SMatthew Knepley 
851649ef022SMatthew Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
852649ef022SMatthew Knepley {
853649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
854649ef022SMatthew Knepley   void            *ctxs[3];
855649ef022SMatthew Knepley   DM               dm;
856649ef022SMatthew Knepley   PetscDS          ds;
857649ef022SMatthew Knepley   Vec              v;
858649ef022SMatthew Knepley   PetscReal        ferrors[3];
859649ef022SMatthew Knepley   PetscInt         f;
860649ef022SMatthew Knepley   PetscErrorCode   ierr;
861649ef022SMatthew Knepley 
862649ef022SMatthew Knepley   PetscFunctionBeginUser;
863649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
864649ef022SMatthew Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
865649ef022SMatthew Knepley 
866649ef022SMatthew Knepley   for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);}
867649ef022SMatthew Knepley   ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr);
868649ef022SMatthew Knepley   ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2]);CHKERRQ(ierr);
869649ef022SMatthew Knepley 
870649ef022SMatthew Knepley   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
871649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
872649ef022SMatthew Knepley   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
873649ef022SMatthew Knepley   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
874649ef022SMatthew Knepley 
875649ef022SMatthew Knepley   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
876649ef022SMatthew Knepley   ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr);
877649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr);
878649ef022SMatthew Knepley   ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr);
879649ef022SMatthew Knepley   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
880649ef022SMatthew Knepley 
881649ef022SMatthew Knepley   PetscFunctionReturn(0);
882649ef022SMatthew Knepley }
883649ef022SMatthew Knepley 
884649ef022SMatthew Knepley int main(int argc, char **argv)
885649ef022SMatthew Knepley {
886649ef022SMatthew Knepley   DM              dm;   /* problem definition */
887649ef022SMatthew Knepley   TS              ts;   /* timestepper */
888649ef022SMatthew Knepley   Vec             u;    /* solution */
889649ef022SMatthew Knepley   AppCtx          user; /* user-defined work context */
890649ef022SMatthew Knepley   PetscReal       t;
891649ef022SMatthew Knepley   PetscErrorCode  ierr;
892649ef022SMatthew Knepley 
893649ef022SMatthew Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
894649ef022SMatthew Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
895649ef022SMatthew Knepley   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
896649ef022SMatthew Knepley   ierr = SetupParameters(&user);CHKERRQ(ierr);
897649ef022SMatthew Knepley   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
898649ef022SMatthew Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
899649ef022SMatthew Knepley   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
900649ef022SMatthew Knepley   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
901649ef022SMatthew Knepley   /* Setup problem */
902649ef022SMatthew Knepley   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
903649ef022SMatthew Knepley   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
904649ef022SMatthew Knepley 
905649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
906649ef022SMatthew Knepley   ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
907649ef022SMatthew Knepley 
908649ef022SMatthew Knepley   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr);
909649ef022SMatthew Knepley   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr);
910649ef022SMatthew Knepley   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr);
911649ef022SMatthew Knepley   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
912649ef022SMatthew Knepley   ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr);
913649ef022SMatthew Knepley   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
914649ef022SMatthew Knepley 
915649ef022SMatthew Knepley   ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */
916649ef022SMatthew Knepley   ierr = SetInitialConditions(ts, u);CHKERRQ(ierr);
917649ef022SMatthew Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
918649ef022SMatthew Knepley   ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr);
919649ef022SMatthew Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
920649ef022SMatthew Knepley   ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
921649ef022SMatthew Knepley 
922649ef022SMatthew Knepley   ierr = TSSolve(ts, u);CHKERRQ(ierr);
923649ef022SMatthew Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
924649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
925649ef022SMatthew Knepley 
926649ef022SMatthew Knepley   ierr = VecDestroy(&u);CHKERRQ(ierr);
927649ef022SMatthew Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
928649ef022SMatthew Knepley   ierr = TSDestroy(&ts);CHKERRQ(ierr);
929649ef022SMatthew Knepley   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
930649ef022SMatthew Knepley   ierr = PetscFinalize();
931649ef022SMatthew Knepley   return ierr;
932649ef022SMatthew Knepley }
933649ef022SMatthew Knepley 
934649ef022SMatthew Knepley /*TEST
935649ef022SMatthew Knepley 
936649ef022SMatthew Knepley   test:
937649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1
938649ef022SMatthew Knepley     requires: triangle !single
939649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
940649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
941649ef022SMatthew Knepley       -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
942649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
943649ef022SMatthew 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 \
944649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
945649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
946649ef022SMatthew Knepley 
947649ef022SMatthew Knepley   # TODO Need trig t for convergence in time, also need to refine in space
948649ef022SMatthew Knepley   test:
949649ef022SMatthew Knepley     # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0]
950649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1_tconv
951649ef022SMatthew Knepley     requires: triangle !single
952649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic_trig -dm_refine 0 \
953649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
954649ef022SMatthew Knepley       -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \
955649ef022SMatthew Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
956649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
957649ef022SMatthew 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 \
958649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
959649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
960649ef022SMatthew Knepley 
961649ef022SMatthew Knepley   test:
962649ef022SMatthew Knepley     # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9]
963649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1_sconv
964649ef022SMatthew Knepley     requires: triangle !single
965649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
966649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
967649ef022SMatthew Knepley       -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
968649ef022SMatthew Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
969649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
970649ef022SMatthew 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 \
971649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
972649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
973649ef022SMatthew Knepley 
974649ef022SMatthew Knepley   test:
975649ef022SMatthew Knepley     suffix: 2d_tri_p3_p2_p2
976649ef022SMatthew Knepley     requires: triangle !single
977649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
978649ef022SMatthew Knepley       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
979649ef022SMatthew Knepley       -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
980649ef022SMatthew Knepley       -snes_convergence_test correct_pressure \
981649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
982649ef022SMatthew 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 \
983649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
984649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
985649ef022SMatthew Knepley 
986*606d57d4SMatthew G. Knepley   test:
987*606d57d4SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1]
988*606d57d4SMatthew G. Knepley     suffix: 2d_tri_p2_p1_p1_tg_sconv
989*606d57d4SMatthew G. Knepley     requires: triangle !single
990*606d57d4SMatthew G. Knepley     args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \
991*606d57d4SMatthew G. Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
992*606d57d4SMatthew G. Knepley       -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
993*606d57d4SMatthew G. Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
994*606d57d4SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
995*606d57d4SMatthew G. 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 \
996*606d57d4SMatthew G. Knepley         -fieldsplit_0_pc_type lu \
997*606d57d4SMatthew G. Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
998*606d57d4SMatthew G. Knepley 
999*606d57d4SMatthew G. Knepley   test:
1000*606d57d4SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2]
1001*606d57d4SMatthew G. Knepley     suffix: 2d_tri_p2_p1_p1_tg_tconv
1002*606d57d4SMatthew G. Knepley     requires: triangle !single
1003*606d57d4SMatthew G. Knepley     args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \
1004*606d57d4SMatthew G. Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1005*606d57d4SMatthew G. Knepley       -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \
1006*606d57d4SMatthew G. Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
1007*606d57d4SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
1008*606d57d4SMatthew G. 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 \
1009*606d57d4SMatthew G. Knepley         -fieldsplit_0_pc_type lu \
1010*606d57d4SMatthew G. Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1011*606d57d4SMatthew G. Knepley 
1012649ef022SMatthew Knepley TEST*/
1013