xref: /petsc/src/ts/tutorials/ex76.c (revision efa12513287cff49a2b9648ae83199dcbfaad71a)
1 static char help[] = "Time-dependent Low Mach Flow in 2d and 3d channels with finite elements.\n\
2 We solve the Low Mach flow problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4 
5 /*F
6 This Low Mach flow is time-dependent isoviscous Navier-Stokes flow. We discretize using the
7 finite element method on an unstructured mesh. The weak form equations are
8 
9 \begin{align*}
10     < q, \nabla\cdot u > = 0
11     <v, du/dt> + <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
12     < w, u \cdot \nabla T > + < \nabla w, \alpha \nabla T > - < w, Q > = 0
13 \end{align*}
14 
15 where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
16 
17 For visualization, use
18 
19   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
20 F*/
21 
22 #include <petscdmplex.h>
23 #include <petscsnes.h>
24 #include <petscts.h>
25 #include <petscds.h>
26 #include <petscbag.h>
27 
28 typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, SOL_TAYLOR_GREEN, NUM_SOL_TYPES} SolType;
29 const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "unknown"};
30 
31 typedef struct {
32   PetscReal nu;    /* Kinematic viscosity */
33   PetscReal alpha; /* Thermal diffusivity */
34   PetscReal T_in;  /* Inlet temperature*/
35 } Parameter;
36 
37 typedef struct {
38   /* Problem definition */
39   PetscBag bag;     /* Holds problem parameters */
40   SolType  solType; /* MMS solution type */
41 } AppCtx;
42 
43 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
44 {
45   PetscInt d;
46   for (d = 0; d < Nc; ++d) u[d] = 0.0;
47   return 0;
48 }
49 
50 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51 {
52   PetscInt d;
53   for (d = 0; d < Nc; ++d) u[d] = 1.0;
54   return 0;
55 }
56 
57 /*
58   CASE: quadratic
59   In 2D we use exact solution:
60 
61     u = t + x^2 + y^2
62     v = t + 2x^2 - 2xy
63     p = x + y - 1
64     T = t + x + y
65     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>
66     Q = 1 + 2t + 3x^2 - 2xy + y^2
67 
68   so that
69 
70     \nabla \cdot u = 2x - 2x = 0
71 
72   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
73     = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1>
74     = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> + <-4 \nu + 2, -4\nu + 2>
75     = <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2 - 4\nu + 2, t (2x - 2y) + 4xy^2 + 2x^2y - 2y^3 - 4\nu + 2>
76 
77   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
78     = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0
79     = 1 + 2t + 3x^2 - 2xy + y^2
80 */
81 
82 static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
83 {
84   u[0] = time + X[0]*X[0] + X[1]*X[1];
85   u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
86   return 0;
87 }
88 static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
89 {
90   u[0] = 1.0;
91   u[1] = 1.0;
92   return 0;
93 }
94 
95 static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
96 {
97   p[0] = X[0] + X[1] - 1.0;
98   return 0;
99 }
100 
101 static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
102 {
103   T[0] = time + X[0] + X[1];
104   return 0;
105 }
106 static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
107 {
108   T[0] = 1.0;
109   return 0;
110 }
111 
112 /* f0_v = du/dt - f */
113 static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
117 {
118   const PetscReal nu = PetscRealPart(constants[0]);
119   PetscInt        Nc = dim;
120   PetscInt        c, d;
121 
122   for (d = 0; d<dim; ++d) f0[d] = u_t[uOff[0]+d];
123 
124   for (c = 0; c < Nc; ++c) {
125     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
126   }
127   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);
128   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);
129 }
130 
131 /* f0_w = dT/dt + u.grad(T) - Q */
132 static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
136 {
137   PetscInt d;
138   f0[0] = 0;
139   for (d = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
140   f0[0] += u_t[uOff[2]] - (2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1]);
141 }
142 
143 /*
144   CASE: cubic
145   In 2D we use exact solution:
146 
147     u = t + x^3 + y^3
148     v = t + 2x^3 - 3x^2y
149     p = 3/2 x^2 + 3/2 y^2 - 1
150     T = t + 1/2 x^2 + 1/2 y^2
151     f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1,
152           t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1>
153     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1
154 
155   so that
156 
157     \nabla \cdot u = 3x^2 - 3x^2 = 0
158 
159   du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f
160   = <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
161 
162   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
163 */
164 static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
165 {
166   u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
167   u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
168   return 0;
169 }
170 static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
171 {
172   u[0] = 1.0;
173   u[1] = 1.0;
174   return 0;
175 }
176 
177 static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
178 {
179   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
180   return 0;
181 }
182 
183 static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
184 {
185   T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
186   return 0;
187 }
188 static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
189 {
190   T[0] = 1.0;
191   return 0;
192 }
193 
194 
195 static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
196                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
197                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
198                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
199 {
200   PetscInt                   c, d;
201   PetscInt                   Nc = dim;
202   const PetscReal            nu = PetscRealPart(constants[0]);
203 
204   for (d=0; d<dim; ++d) f0[d] = u_t[uOff[0]+d];
205 
206   for (c=0; c<Nc; ++c) {
207     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
208   }
209   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);
210   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);
211 }
212 
213 static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
214                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
215                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
216                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
217 {
218   PetscInt              d;
219   const PetscReal alpha = PetscRealPart(constants[1]);
220 
221   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
222   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);
223 }
224 
225 /*
226   CASE: cubic-trigonometric
227   In 2D we use exact solution:
228 
229     u = beta cos t + x^3 + y^3
230     v = beta sin t + 2x^3 - 3x^2y
231     p = 3/2 x^2 + 3/2 y^2 - 1
232     T = 20 cos t + 1/2 x^2 + 1/2 y^2
233     f = < beta cos t 3x^2         + beta sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y)  + 3x,
234           beta cos t (6x^2 - 6xy) - beta sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu(12x - 6y) + 3y>
235     Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha
236 
237   so that
238 
239     \nabla \cdot u = 3x^2 - 3x^2 = 0
240 
241   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
242     = <-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>
243     = <-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>
244     = <cos t (3x^2)       + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y)  + 3x,
245        cos t (6x^2 - 6xy) - sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu (12x - 6y) + 3y>
246 
247   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
248     = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha
249     = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha
250     = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha)
251 */
252 static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
253 {
254   u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
255   u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
256   return 0;
257 }
258 static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
259 {
260   u[0] = -100.*PetscSinReal(time);
261   u[1] =  100.*PetscCosReal(time);
262   return 0;
263 }
264 
265 static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
266 {
267   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
268   return 0;
269 }
270 
271 static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
272 {
273   T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
274   return 0;
275 }
276 static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
277 {
278   T[0] = -100.*PetscSinReal(time);
279   return 0;
280 }
281 
282 static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
283                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
284                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
285                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
286 {
287   const PetscReal nu = PetscRealPart(constants[0]);
288   PetscInt        Nc = dim;
289   PetscInt        c, d;
290 
291   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d];
292 
293   for (c = 0; c < Nc; ++c) {
294     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
295   }
296   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];
297   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];
298 }
299 
300 static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
301                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
302                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
303                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
304 {
305   const PetscReal alpha = PetscRealPart(constants[1]);
306   PetscInt        d;
307 
308   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
309   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);
310 }
311 
312 /*
313   CASE: taylor-green vortex
314   In 2D we use exact solution:
315 
316     u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)
317     v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)
318     p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t)
319     T = t + x + y
320     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     Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t)
322 
323   so that
324 
325   \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 
327   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
328     = <-\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         \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     + < \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         \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     + <-\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        -\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     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
335         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
336     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
337         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
338     = <-\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         \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     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
341         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
342     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
343        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
344     + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
345        -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
346     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
347         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
348     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
349         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
350     = <-\pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t),
351         \pi (sin(\pi(x - t)) sin(\pi(y - t)) - cos(\pi(x - t)) cos(\pi(y - t))) exp(-2 \pi^2 \nu t)>
352     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
353         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
354     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
355        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
356     = < \pi cos(\pi(x - t)) cos(\pi(y - t)),
357         \pi sin(\pi(x - t)) sin(\pi(y - t))>
358     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)),
359        -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0
360   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
361     = 1 + u \cdot <1, 1> - 0
362     = 1 + u + v
363 */
364 
365 static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
366 {
367   u[0] = 1 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
368   u[1] = 1 + PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
369   return 0;
370 }
371 static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
372 {
373   u[0] = -PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
374                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
375                   - 2*PETSC_PI*PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
376   u[1] =  PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
377                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
378                   - 2*PETSC_PI*PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
379   return 0;
380 }
381 
382 static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
383 {
384   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   return 0;
386 }
387 
388 static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
389 {
390   p[0] = PETSC_PI*(0.5*(PetscSinReal(2*PETSC_PI*(X[0]-time)) + PetscSinReal(2*PETSC_PI*(X[1]-time)))
391                  + PETSC_PI*(PetscCosReal(2*PETSC_PI*(X[0]-time)) + PetscCosReal(2*PETSC_PI*(X[1]-time))))*PetscExpReal(-4*PETSC_PI*PETSC_PI*time);
392   return 0;
393 }
394 
395 static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
396 {
397   T[0] = time + X[0] + X[1];
398   return 0;
399 }
400 static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
401 {
402   T[0] = 1.0;
403   return 0;
404 }
405 
406 static void f0_taylor_green_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
407                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
408                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
409                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
410 {
411   PetscInt        Nc = dim;
412   PetscInt        c, d;
413 
414   for (d = 0; d < dim; ++d) f0[d] = u_t[uOff[0]+d];
415 
416   for (c = 0; c < Nc; ++c) {
417     for (d = 0; d < dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
418   }
419 }
420 
421 static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
422                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
423                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
424                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
425 {
426   PetscScalar vel[2];
427   PetscInt    d;
428 
429   taylor_green_u(dim, t, X, Nf, vel, NULL);
430   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
431   f0[0] += u_t[uOff[2]] - (1.0 + vel[0] + vel[1]);
432 }
433 
434 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
435                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
436                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
437                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
438 {
439   PetscInt d;
440   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
441 }
442 
443 /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
444 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
445                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
446                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
447                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
448 {
449   const PetscReal nu = PetscRealPart(constants[0]);
450   const PetscInt    Nc = dim;
451   PetscInt        c, d;
452 
453   for (c = 0; c < Nc; ++c) {
454     for (d = 0; d < dim; ++d) {
455       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
456     }
457     f1[c*dim+c] -= u[uOff[1]];
458   }
459 }
460 
461 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
462                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
463                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
464                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
465 {
466   const PetscReal alpha = PetscRealPart(constants[1]);
467   PetscInt d;
468   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
469 }
470 
471 /*Jacobians*/
472 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
473                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
474                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
475                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
476 {
477   PetscInt d;
478   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
479 }
480 
481 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
482                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
483                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
484                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
485 {
486   PetscInt c, d;
487   const PetscInt  Nc = dim;
488 
489   for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift;
490 
491   for (c = 0; c < Nc; ++c) {
492     for (d = 0; d < dim; ++d) {
493       g0[c*Nc+d] += u_x[ c*Nc+d];
494     }
495   }
496 }
497 
498 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
499                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
500                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
501                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
502 {
503   PetscInt NcI = dim;
504   PetscInt NcJ = dim;
505   PetscInt c, d, e;
506 
507   for (c = 0; c < NcI; ++c) {
508     for (d = 0; d < NcJ; ++d) {
509       for (e = 0; e < dim; ++e) {
510         if (c == d) {
511           g1[(c*NcJ+d)*dim+e] += u[e];
512         }
513       }
514     }
515   }
516 }
517 
518 
519 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
520                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
521                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
522                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
523 {
524   PetscInt d;
525   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
526 }
527 
528 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
529                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
530                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
531                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
532 {
533    const PetscReal nu = PetscRealPart(constants[0]);
534    const PetscInt  Nc = dim;
535    PetscInt        c, d;
536 
537   for (c = 0; c < Nc; ++c) {
538     for (d = 0; d < dim; ++d) {
539       g3[((c*Nc+c)*dim+d)*dim+d] += nu;
540       g3[((c*Nc+d)*dim+d)*dim+c] += nu;
541     }
542   }
543 }
544 
545 static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
546                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
547                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
548                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
549 {
550   PetscInt d;
551   for (d = 0; d < dim; ++d) g0[d] = u_tShift;
552 }
553 
554 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
555                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
556                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
557                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
558 {
559   PetscInt d;
560   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
561 }
562 
563 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
564                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
565                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
566                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
567 {
568   PetscInt d;
569   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
570 }
571 
572 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
573                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
574                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
575                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
576 {
577   const PetscReal alpha = PetscRealPart(constants[1]);
578   PetscInt               d;
579 
580   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
581 }
582 
583 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
584 {
585   PetscInt       sol;
586   PetscErrorCode ierr;
587 
588 
589   PetscFunctionBeginUser;
590   options->solType = SOL_QUADRATIC;
591 
592   ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr);
593   sol = options->solType;
594   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
595   options->solType = (SolType) sol;
596   ierr = PetscOptionsEnd();
597   PetscFunctionReturn(0);
598 }
599 
600 static PetscErrorCode SetupParameters(AppCtx *user)
601 {
602   PetscBag       bag;
603   Parameter     *p;
604   PetscErrorCode ierr;
605 
606   PetscFunctionBeginUser;
607   /* setup PETSc parameter bag */
608   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
609   ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr);
610   bag  = user->bag;
611   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
612   ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr);
613   ierr = PetscBagRegisterReal(bag, &p->T_in,  1.0, "T_in",  "Inlet temperature");CHKERRQ(ierr);
614   PetscFunctionReturn(0);
615 }
616 
617 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
618 {
619   PetscErrorCode ierr;
620 
621   PetscFunctionBeginUser;
622   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
623   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
624   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
625   PetscFunctionReturn(0);
626 }
627 
628 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
629 {
630   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
631   PetscErrorCode (*exactFuncs_t[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
632   PetscDS          prob;
633   Parameter       *ctx;
634   PetscInt         id;
635   PetscErrorCode   ierr;
636 
637   PetscFunctionBeginUser;
638   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
639   switch(user->solType){
640   case SOL_QUADRATIC:
641     ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr);
642     ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr);
643 
644     exactFuncs[0]   = quadratic_u;
645     exactFuncs[1]   = quadratic_p;
646     exactFuncs[2]   = quadratic_T;
647     exactFuncs_t[0] = quadratic_u_t;
648     exactFuncs_t[1] = NULL;
649     exactFuncs_t[2] = quadratic_T_t;
650     break;
651   case SOL_CUBIC:
652     ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr);
653     ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr);
654 
655     exactFuncs[0]   = cubic_u;
656     exactFuncs[1]   = cubic_p;
657     exactFuncs[2]   = cubic_T;
658     exactFuncs_t[0] = cubic_u_t;
659     exactFuncs_t[1] = NULL;
660     exactFuncs_t[2] = cubic_T_t;
661     break;
662   case SOL_CUBIC_TRIG:
663     ierr = PetscDSSetResidual(prob, 0, f0_cubic_trig_v, f1_v);CHKERRQ(ierr);
664     ierr = PetscDSSetResidual(prob, 2, f0_cubic_trig_w, f1_w);CHKERRQ(ierr);
665 
666     exactFuncs[0]   = cubic_trig_u;
667     exactFuncs[1]   = cubic_trig_p;
668     exactFuncs[2]   = cubic_trig_T;
669     exactFuncs_t[0] = cubic_trig_u_t;
670     exactFuncs_t[1] = NULL;
671     exactFuncs_t[2] = cubic_trig_T_t;
672     break;
673   case SOL_TAYLOR_GREEN:
674     ierr = PetscDSSetResidual(prob, 0, f0_taylor_green_v, f1_v);CHKERRQ(ierr);
675     ierr = PetscDSSetResidual(prob, 2, f0_taylor_green_w, f1_w);CHKERRQ(ierr);
676 
677     exactFuncs[0]   = taylor_green_u;
678     exactFuncs[1]   = taylor_green_p;
679     exactFuncs[2]   = taylor_green_T;
680     exactFuncs_t[0] = taylor_green_u_t;
681     exactFuncs_t[1] = taylor_green_p_t;
682     exactFuncs_t[2] = taylor_green_T_t;
683     break;
684    default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
685   }
686 
687   ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr);
688 
689   ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu);CHKERRQ(ierr);
690   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL);CHKERRQ(ierr);
691   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL);CHKERRQ(ierr);
692   ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL);CHKERRQ(ierr);
693   ierr = PetscDSSetJacobian(prob, 2, 2, g0_wT, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
694   /* Setup constants */
695   {
696     Parameter  *param;
697     PetscScalar constants[3];
698 
699     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
700 
701     constants[0] = param->nu;
702     constants[1] = param->alpha;
703     constants[2] = param->T_in;
704     ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr);
705   }
706   /* Setup Boundary Conditions */
707   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
708   id   = 3;
709   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);
710   id   = 1;
711   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);
712   id   = 2;
713   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);
714   id   = 4;
715   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);
716   id   = 3;
717   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);
718   id   = 1;
719   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);
720   id   = 2;
721   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);
722   id   = 4;
723   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);
724 
725   /*setup exact solution.*/
726   ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr);
727   ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr);
728   ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr);
729   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 0, exactFuncs_t[0], ctx);CHKERRQ(ierr);
730   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 1, exactFuncs_t[1], ctx);CHKERRQ(ierr);
731   ierr = PetscDSSetExactSolutionTimeDerivative(prob, 2, exactFuncs_t[2], ctx);CHKERRQ(ierr);
732   PetscFunctionReturn(0);
733 }
734 
735 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
736 {
737   DM              cdm   = dm;
738   PetscFE         fe[3];
739   Parameter      *param;
740   MPI_Comm        comm;
741   DMPolytopeType  ct;
742   PetscInt        dim, cStart;
743   PetscBool       simplex;
744   PetscErrorCode  ierr;
745 
746   PetscFunctionBeginUser;
747   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
748   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
749   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
750   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
751   /* Create finite element */
752   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
753   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
754   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
755 
756   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
757   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
758   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
759 
760   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
761   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
762   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
763 
764   /* Set discretization and boundary conditions for each mesh */
765   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
766   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
767   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
768   ierr = DMCreateDS(dm);CHKERRQ(ierr);
769   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
770   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
771   while (cdm) {
772     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
773     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
774   }
775   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
776   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
777   ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr);
778 
779   {
780     PetscObject  pressure;
781     MatNullSpace nullspacePres;
782 
783     ierr = DMGetField(dm, 1, NULL, &pressure);CHKERRQ(ierr);
784     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
785     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
786     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
787   }
788 
789   PetscFunctionReturn(0);
790 }
791 
792 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
793 {
794   Vec              vec;
795   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
796   PetscErrorCode   ierr;
797 
798   PetscFunctionBeginUser;
799   if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield);
800   funcs[nfield] = constant;
801   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
802   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
803   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
804   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
805   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
806   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
807   ierr = VecDestroy(&vec);CHKERRQ(ierr);
808   PetscFunctionReturn(0);
809 }
810 
811 static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u)
812 {
813   DM             dm;
814   MatNullSpace   nullsp;
815   PetscErrorCode ierr;
816 
817   PetscFunctionBegin;
818   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
819   ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr);
820   ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr);
821   ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
822   PetscFunctionReturn(0);
823 }
824 
825 /* Make the discrete pressure discretely divergence free */
826 static PetscErrorCode RemoveDiscretePressureNullspace(TS ts)
827 {
828   Vec            u;
829   PetscErrorCode ierr;
830 
831   PetscFunctionBegin;
832   ierr = TSGetSolution(ts, &u);CHKERRQ(ierr);
833   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
834   PetscFunctionReturn(0);
835 }
836 
837 static PetscErrorCode SetInitialConditions(TS ts, Vec u)
838 {
839   DM             dm;
840   PetscReal      t;
841   PetscErrorCode ierr;
842 
843   PetscFunctionBegin;
844   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
845   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
846   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
847   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
848   PetscFunctionReturn(0);
849 }
850 
851 static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
852 {
853   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
854   void            *ctxs[3];
855   DM               dm;
856   PetscDS          ds;
857   Vec              v;
858   PetscReal        ferrors[3];
859   PetscInt         f;
860   PetscErrorCode   ierr;
861 
862   PetscFunctionBeginUser;
863   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
864   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
865 
866   for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);}
867   ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr);
868   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);
869 
870   ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr);
871   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
872   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
873   ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr);
874 
875   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
876   ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr);
877   ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr);
878   ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr);
879   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
880 
881   PetscFunctionReturn(0);
882 }
883 
884 int main(int argc, char **argv)
885 {
886   DM              dm;   /* problem definition */
887   TS              ts;   /* timestepper */
888   Vec             u;    /* solution */
889   AppCtx          user; /* user-defined work context */
890   PetscReal       t;
891   PetscErrorCode  ierr;
892 
893   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
894   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
895   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
896   ierr = SetupParameters(&user);CHKERRQ(ierr);
897   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
898   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
899   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
900   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
901   /* Setup problem */
902   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
903   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
904 
905   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
906   ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
907 
908   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr);
909   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr);
910   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr);
911   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
912   ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr);
913   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
914 
915   ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */
916   ierr = SetInitialConditions(ts, u);CHKERRQ(ierr);
917   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
918   ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr);
919   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
920   ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
921 
922   ierr = TSSolve(ts, u);CHKERRQ(ierr);
923   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
924   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
925 
926   ierr = VecDestroy(&u);CHKERRQ(ierr);
927   ierr = DMDestroy(&dm);CHKERRQ(ierr);
928   ierr = TSDestroy(&ts);CHKERRQ(ierr);
929   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
930   ierr = PetscFinalize();
931   return ierr;
932 }
933 
934 /*TEST
935 
936   test:
937     suffix: 2d_tri_p2_p1_p1
938     requires: triangle !single
939     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
940       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
941       -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
942       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
943       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
944         -fieldsplit_0_pc_type lu \
945         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
946 
947   # TODO Need trig t for convergence in time, also need to refine in space
948   test:
949     # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0]
950     suffix: 2d_tri_p2_p1_p1_tconv
951     requires: triangle !single
952     args: -dm_plex_separate_marker -sol_type cubic_trig -dm_refine 0 \
953       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
954       -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \
955       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
956       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
957       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
958         -fieldsplit_0_pc_type lu \
959         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
960 
961   test:
962     # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9]
963     suffix: 2d_tri_p2_p1_p1_sconv
964     requires: triangle !single
965     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
966       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
967       -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
968       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
969       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
970       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
971         -fieldsplit_0_pc_type lu \
972         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
973 
974   test:
975     suffix: 2d_tri_p3_p2_p2
976     requires: triangle !single
977     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
978       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
979       -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
980       -snes_convergence_test correct_pressure \
981       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
982       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
983         -fieldsplit_0_pc_type lu \
984         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
985 
986   test:
987     # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1]
988     suffix: 2d_tri_p2_p1_p1_tg_sconv
989     requires: triangle !single
990     args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \
991       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
992       -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
993       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
994       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
995       -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         -fieldsplit_0_pc_type lu \
997         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
998 
999   test:
1000     # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2]
1001     suffix: 2d_tri_p2_p1_p1_tg_tconv
1002     requires: triangle !single
1003     args: -dm_plex_separate_marker -sol_type taylor_green -dm_refine 0 \
1004       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1005       -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1 \
1006       -snes_error_if_not_converged -snes_convergence_test correct_pressure \
1007       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_atol 1e-16 -ksp_error_if_not_converged \
1008       -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         -fieldsplit_0_pc_type lu \
1010         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1011 
1012 TEST*/
1013