xref: /petsc/src/ts/tutorials/ex76.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1649ef022SMatthew Knepley static char help[] = "Time-dependent Low Mach Flow in 2d and 3d channels with finite elements.\n\
2444129b9SMatthew G. Knepley We solve the Low Mach flow problem for both conducting and non-conducting fluids,\n\
3444129b9SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4649ef022SMatthew Knepley 
5649ef022SMatthew Knepley /*F
6444129b9SMatthew G. Knepley The non-conducting 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 
17444129b9SMatthew G. Knepley The conducting form is given in the ABLATE documentation [1,2] and derived in Principe and Codina [2].
18444129b9SMatthew G. Knepley 
19649ef022SMatthew Knepley For visualization, use
20649ef022SMatthew Knepley 
21649ef022SMatthew Knepley   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
22444129b9SMatthew G. Knepley 
23444129b9SMatthew G. Knepley [1] https://ubchrest.github.io/ablate/content/formulations/lowMachFlow/
24444129b9SMatthew G. Knepley [2] https://github.com/UBCHREST/ablate/blob/main/ablateCore/flow/lowMachFlow.c
25444129b9SMatthew G. Knepley [3] J. Principe and R. Codina, "Mathematical models for thermally coupled low speed flows", Adv. in Theo. and App. Mech., 2(1), pp.93--112, 2009.
26649ef022SMatthew Knepley F*/
27649ef022SMatthew Knepley 
28649ef022SMatthew Knepley #include <petscdmplex.h>
29649ef022SMatthew Knepley #include <petscsnes.h>
30649ef022SMatthew Knepley #include <petscts.h>
31649ef022SMatthew Knepley #include <petscds.h>
32649ef022SMatthew Knepley #include <petscbag.h>
33649ef022SMatthew Knepley 
34444129b9SMatthew G. Knepley typedef enum {MOD_INCOMPRESSIBLE, MOD_CONDUCTING, NUM_MOD_TYPES} ModType;
35444129b9SMatthew G. Knepley const char *modTypes[NUM_MOD_TYPES+1] = {"incompressible", "conducting", "unknown"};
36444129b9SMatthew G. Knepley 
37444129b9SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, SOL_CUBIC_TRIG, SOL_TAYLOR_GREEN, SOL_PIPE, NUM_SOL_TYPES} SolType;
38444129b9SMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "pipe", "unknown"};
39444129b9SMatthew G. Knepley 
40444129b9SMatthew G. Knepley /* Fields */
41444129b9SMatthew G. Knepley const PetscInt VEL      = 0;
42444129b9SMatthew G. Knepley const PetscInt PRES     = 1;
43444129b9SMatthew G. Knepley const PetscInt TEMP     = 2;
44444129b9SMatthew G. Knepley /* Sources */
45444129b9SMatthew G. Knepley const PetscInt MOMENTUM = 0;
46444129b9SMatthew G. Knepley const PetscInt MASS     = 1;
47444129b9SMatthew G. Knepley const PetscInt ENERGY   = 2;
48444129b9SMatthew G. Knepley /* Constants */
49444129b9SMatthew G. Knepley const PetscInt STROUHAL = 0;
50444129b9SMatthew G. Knepley const PetscInt FROUDE   = 1;
51444129b9SMatthew G. Knepley const PetscInt REYNOLDS = 2;
52444129b9SMatthew G. Knepley const PetscInt PECLET   = 3;
53444129b9SMatthew G. Knepley const PetscInt P_TH     = 4;
54444129b9SMatthew G. Knepley const PetscInt MU       = 5;
55444129b9SMatthew G. Knepley const PetscInt NU       = 6;
56444129b9SMatthew G. Knepley const PetscInt C_P      = 7;
57444129b9SMatthew G. Knepley const PetscInt K        = 8;
58444129b9SMatthew G. Knepley const PetscInt ALPHA    = 9;
59444129b9SMatthew G. Knepley const PetscInt T_IN     = 10;
60444129b9SMatthew G. Knepley const PetscInt G_DIR    = 11;
61649ef022SMatthew Knepley 
62649ef022SMatthew Knepley typedef struct {
63444129b9SMatthew G. Knepley   PetscReal Strouhal; /* Strouhal number */
64444129b9SMatthew G. Knepley   PetscReal Froude;   /* Froude number */
65444129b9SMatthew G. Knepley   PetscReal Reynolds; /* Reynolds number */
66444129b9SMatthew G. Knepley   PetscReal Peclet;   /* Peclet number */
67444129b9SMatthew G. Knepley   PetscReal p_th;     /* Thermodynamic pressure */
68444129b9SMatthew G. Knepley   PetscReal mu;       /* Dynamic viscosity */
69649ef022SMatthew Knepley   PetscReal nu;       /* Kinematic viscosity */
70444129b9SMatthew G. Knepley   PetscReal c_p;      /* Specific heat at constant pressure */
71444129b9SMatthew G. Knepley   PetscReal k;        /* Thermal conductivity */
72649ef022SMatthew Knepley   PetscReal alpha;    /* Thermal diffusivity */
73649ef022SMatthew Knepley   PetscReal T_in;     /* Inlet temperature */
74444129b9SMatthew G. Knepley   PetscReal g_dir;    /* Gravity direction */
75649ef022SMatthew Knepley } Parameter;
76649ef022SMatthew Knepley 
77649ef022SMatthew Knepley typedef struct {
78649ef022SMatthew Knepley   /* Problem definition */
79649ef022SMatthew Knepley   PetscBag  bag;          /* Holds problem parameters */
80444129b9SMatthew G. Knepley   ModType   modType;      /* Model type */
81649ef022SMatthew Knepley   SolType   solType;      /* MMS solution type */
82444129b9SMatthew G. Knepley   PetscBool hasNullSpace; /* Problem has the constant null space for pressure */
83a712f3bbSMatthew G. Knepley   /* Flow diagnostics */
84a712f3bbSMatthew G. Knepley   DM        dmCell;       /* A DM with piecewise constant discretization */
85649ef022SMatthew Knepley } AppCtx;
86649ef022SMatthew Knepley 
87649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
88649ef022SMatthew Knepley {
89649ef022SMatthew Knepley   PetscInt d;
90649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
91649ef022SMatthew Knepley   return 0;
92649ef022SMatthew Knepley }
93649ef022SMatthew Knepley 
94649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
95649ef022SMatthew Knepley {
96649ef022SMatthew Knepley   PetscInt d;
97649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
98649ef022SMatthew Knepley   return 0;
99649ef022SMatthew Knepley }
100649ef022SMatthew Knepley 
101649ef022SMatthew Knepley /*
102649ef022SMatthew Knepley   CASE: quadratic
103649ef022SMatthew Knepley   In 2D we use exact solution:
104649ef022SMatthew Knepley 
105649ef022SMatthew Knepley     u = t + x^2 + y^2
106649ef022SMatthew Knepley     v = t + 2x^2 - 2xy
107649ef022SMatthew Knepley     p = x + y - 1
108444129b9SMatthew G. Knepley     T = t + x + y + 1
109649ef022SMatthew 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>
110649ef022SMatthew Knepley     Q = 1 + 2t + 3x^2 - 2xy + y^2
111649ef022SMatthew Knepley 
112649ef022SMatthew Knepley   so that
113649ef022SMatthew Knepley 
114649ef022SMatthew Knepley     \nabla \cdot u = 2x - 2x = 0
115649ef022SMatthew Knepley 
116649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
117649ef022SMatthew Knepley     = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1>
118649ef022SMatthew 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>
119649ef022SMatthew 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>
120649ef022SMatthew Knepley 
121649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
122649ef022SMatthew Knepley     = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0
123649ef022SMatthew Knepley     = 1 + 2t + 3x^2 - 2xy + y^2
124649ef022SMatthew Knepley */
125649ef022SMatthew Knepley 
126649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
127649ef022SMatthew Knepley {
128649ef022SMatthew Knepley   u[0] = time + X[0]*X[0] + X[1]*X[1];
129649ef022SMatthew Knepley   u[1] = time + 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
130649ef022SMatthew Knepley   return 0;
131649ef022SMatthew Knepley }
132649ef022SMatthew Knepley static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
133649ef022SMatthew Knepley {
134649ef022SMatthew Knepley   u[0] = 1.0;
135649ef022SMatthew Knepley   u[1] = 1.0;
136649ef022SMatthew Knepley   return 0;
137649ef022SMatthew Knepley }
138649ef022SMatthew Knepley 
139649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
140649ef022SMatthew Knepley {
141649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
142649ef022SMatthew Knepley   return 0;
143649ef022SMatthew Knepley }
144649ef022SMatthew Knepley 
145649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
146649ef022SMatthew Knepley {
147444129b9SMatthew G. Knepley   T[0] = time + X[0] + X[1] + 1.0;
148649ef022SMatthew Knepley   return 0;
149649ef022SMatthew Knepley }
150649ef022SMatthew Knepley static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
151649ef022SMatthew Knepley {
152649ef022SMatthew Knepley   T[0] = 1.0;
153649ef022SMatthew Knepley   return 0;
154649ef022SMatthew Knepley }
155649ef022SMatthew Knepley 
156649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
160649ef022SMatthew Knepley {
161444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
162649ef022SMatthew Knepley 
163444129b9SMatthew G. 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;
164444129b9SMatthew G. 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;
165649ef022SMatthew Knepley }
166649ef022SMatthew Knepley 
167649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
168649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
169649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
170649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
171649ef022SMatthew Knepley {
172444129b9SMatthew G. Knepley   f0[0] -= 2*t + 1 + 3*X[0]*X[0] - 2*X[0]*X[1] + X[1]*X[1];
173444129b9SMatthew G. Knepley }
174444129b9SMatthew G. Knepley 
175444129b9SMatthew G. Knepley /*
176444129b9SMatthew G. Knepley   CASE: quadratic
177444129b9SMatthew G. Knepley   In 2D we use exact solution:
178444129b9SMatthew G. Knepley 
179444129b9SMatthew G. Knepley     u = t + x^2 + y^2
180444129b9SMatthew G. Knepley     v = t + 2x^2 - 2xy
181444129b9SMatthew G. Knepley     p = x + y - 1
182444129b9SMatthew G. Knepley     T = t + x + y + 1
183444129b9SMatthew G. Knepley   rho = p^{th} / T
184444129b9SMatthew G. Knepley 
185444129b9SMatthew G. Knepley   so that
186444129b9SMatthew G. Knepley 
187444129b9SMatthew G. Knepley     \nabla \cdot u = 2x - 2x = 0
188444129b9SMatthew G. Knepley     grad u = <<2 x, 4x - 2y>, <2 y, -2x>>
189444129b9SMatthew G. Knepley     epsilon(u) = 1/2 (grad u + grad u^T) = <<2x, 2x>, <2x, -2x>>
190444129b9SMatthew G. Knepley     epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u)
191444129b9SMatthew G. Knepley     div epsilon'(u) = <2, 2>
192444129b9SMatthew G. Knepley 
193444129b9SMatthew G. Knepley   f = rho S du/dt + rho u \cdot \nabla u - 2\mu/Re div \epsilon'(u) + \nabla p + rho / F^2 \hat y
194444129b9SMatthew G. Knepley     = rho S <1, 1> + rho <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - 2\mu/Re <2, 2> + <1, 1> + rho/F^2 <0, 1>
195444129b9SMatthew G. Knepley     = rho S <1, 1> + rho <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t (2x - 2y) + 2x^2y + 4xy^2 - 2y^3> - mu/Re <4, 4> + <1, 1> + rho/F^2 <0, 1>
196444129b9SMatthew G. Knepley 
197444129b9SMatthew G. Knepley   g = S rho_t + div (rho u)
198444129b9SMatthew G. Knepley     = -S pth T_t/T^2 + rho div (u) + u . grad rho
199444129b9SMatthew G. Knepley     = -S pth 1/T^2 - pth u . grad T / T^2
200444129b9SMatthew G. Knepley     = -pth / T^2 (S + 2t + 3 x^2 - 2xy + y^2)
201444129b9SMatthew G. Knepley 
202444129b9SMatthew G. Knepley   Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T
203444129b9SMatthew G. Knepley     = c_p S pth / T + c_p pth (2t + 3 x^2 - 2xy + y^2) / T - k/Pe 0
204444129b9SMatthew G. Knepley     = c_p pth / T (S + 2t + 3 x^2 - 2xy + y^2)
205444129b9SMatthew G. Knepley */
206444129b9SMatthew G. Knepley static void f0_conduct_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
207444129b9SMatthew G. Knepley                                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
208444129b9SMatthew G. Knepley                                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
209444129b9SMatthew G. Knepley                                    PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
210444129b9SMatthew G. Knepley {
211444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
212444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
213444129b9SMatthew G. Knepley   const PetscReal Re   = PetscRealPart(constants[REYNOLDS]);
214444129b9SMatthew G. Knepley   const PetscReal mu   = PetscRealPart(constants[MU]);
215444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
216444129b9SMatthew G. Knepley   const PetscReal rho  = p_th / (t + X[0] + X[1] + 1.);
217444129b9SMatthew G. Knepley   const PetscInt  gd   = (PetscInt) PetscRealPart(constants[G_DIR]);
218444129b9SMatthew G. Knepley 
219444129b9SMatthew G. Knepley   f0[0]  -= rho * S + rho * (2.*t*(X[0] + X[1]) + 2.*X[0]*X[0]*X[0] + 4.*X[0]*X[0]*X[1] - 2.*X[0]*X[1]*X[1]) - 4.*mu/Re + 1.;
220444129b9SMatthew G. Knepley   f0[1]  -= rho * S + rho * (2.*t*(X[0] - X[1]) + 2.*X[0]*X[0]*X[1] + 4.*X[0]*X[1]*X[1] - 2.*X[1]*X[1]*X[1]) - 4.*mu/Re + 1.;
221444129b9SMatthew G. Knepley   f0[gd] -= rho/PetscSqr(F);
222444129b9SMatthew G. Knepley }
223444129b9SMatthew G. Knepley 
224444129b9SMatthew G. Knepley static void f0_conduct_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
225444129b9SMatthew G. Knepley                                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
226444129b9SMatthew G. Knepley                                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
227444129b9SMatthew G. Knepley                                    PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
228444129b9SMatthew G. Knepley {
229444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
230444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
231444129b9SMatthew G. Knepley 
232444129b9SMatthew G. Knepley   f0[0] += p_th * (S + 2.*t + 3.*X[0]*X[0] - 2.*X[0]*X[1] + X[1]*X[1]) / PetscSqr(t + X[0] + X[1] + 1.);
233444129b9SMatthew G. Knepley }
234444129b9SMatthew G. Knepley 
235444129b9SMatthew G. Knepley static void f0_conduct_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
236444129b9SMatthew G. Knepley                                    const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
237444129b9SMatthew G. Knepley                                    const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
238444129b9SMatthew G. Knepley                                    PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
239444129b9SMatthew G. Knepley {
240444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
241444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
242444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
243444129b9SMatthew G. Knepley 
244444129b9SMatthew G. Knepley   f0[0] -= c_p*p_th * (S + 2.*t + 3.*X[0]*X[0] - 2.*X[0]*X[1] + X[1]*X[1]) / (t + X[0] + X[1] + 1.);
245649ef022SMatthew Knepley }
246649ef022SMatthew Knepley 
247649ef022SMatthew Knepley /*
248649ef022SMatthew Knepley   CASE: cubic
249649ef022SMatthew Knepley   In 2D we use exact solution:
250649ef022SMatthew Knepley 
251649ef022SMatthew Knepley     u = t + x^3 + y^3
252649ef022SMatthew Knepley     v = t + 2x^3 - 3x^2y
253649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
254649ef022SMatthew Knepley     T = t + 1/2 x^2 + 1/2 y^2
255649ef022SMatthew Knepley     f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1,
256649ef022SMatthew Knepley           t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1>
257649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1
258649ef022SMatthew Knepley 
259649ef022SMatthew Knepley   so that
260649ef022SMatthew Knepley 
261649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
262649ef022SMatthew Knepley 
263649ef022SMatthew Knepley   du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f
264649ef022SMatthew 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
265649ef022SMatthew Knepley 
266649ef022SMatthew 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
267649ef022SMatthew Knepley */
268649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
269649ef022SMatthew Knepley {
270649ef022SMatthew Knepley   u[0] = time + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
271649ef022SMatthew Knepley   u[1] = time + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
272649ef022SMatthew Knepley   return 0;
273649ef022SMatthew Knepley }
274649ef022SMatthew Knepley static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
275649ef022SMatthew Knepley {
276649ef022SMatthew Knepley   u[0] = 1.0;
277649ef022SMatthew Knepley   u[1] = 1.0;
278649ef022SMatthew Knepley   return 0;
279649ef022SMatthew Knepley }
280649ef022SMatthew Knepley 
281649ef022SMatthew Knepley static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
282649ef022SMatthew Knepley {
283649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
284649ef022SMatthew Knepley   return 0;
285649ef022SMatthew Knepley }
286649ef022SMatthew Knepley 
287649ef022SMatthew Knepley static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
288649ef022SMatthew Knepley {
289649ef022SMatthew Knepley   T[0] = time + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
290649ef022SMatthew Knepley   return 0;
291649ef022SMatthew Knepley }
292649ef022SMatthew Knepley static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
293649ef022SMatthew Knepley {
294649ef022SMatthew Knepley   T[0] = 1.0;
295649ef022SMatthew Knepley   return 0;
296649ef022SMatthew Knepley }
297649ef022SMatthew Knepley 
298649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
299649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
300649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
301649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
302649ef022SMatthew Knepley {
303444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
304649ef022SMatthew Knepley 
305649ef022SMatthew 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);
306649ef022SMatthew 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);
307649ef022SMatthew Knepley }
308649ef022SMatthew Knepley 
309649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
310649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
311649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
312649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
313649ef022SMatthew Knepley {
314444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
315649ef022SMatthew Knepley 
316444129b9SMatthew G. Knepley   f0[0] -= 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;
317649ef022SMatthew Knepley }
318649ef022SMatthew Knepley 
319649ef022SMatthew Knepley /*
320649ef022SMatthew Knepley   CASE: cubic-trigonometric
321649ef022SMatthew Knepley   In 2D we use exact solution:
322649ef022SMatthew Knepley 
323649ef022SMatthew Knepley     u = beta cos t + x^3 + y^3
324649ef022SMatthew Knepley     v = beta sin t + 2x^3 - 3x^2y
325649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
326649ef022SMatthew Knepley     T = 20 cos t + 1/2 x^2 + 1/2 y^2
327649ef022SMatthew 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,
328649ef022SMatthew Knepley           beta cos t (6x^2 - 6xy) - beta sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu(12x - 6y) + 3y>
329649ef022SMatthew Knepley     Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha
330649ef022SMatthew Knepley 
331649ef022SMatthew Knepley   so that
332649ef022SMatthew Knepley 
333649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
334649ef022SMatthew Knepley 
335649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
336649ef022SMatthew 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>
337649ef022SMatthew 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>
338649ef022SMatthew Knepley     = <cos t (3x^2)       + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y)  + 3x,
339649ef022SMatthew Knepley        cos t (6x^2 - 6xy) - sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu (12x - 6y) + 3y>
340649ef022SMatthew Knepley 
341649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
342649ef022SMatthew Knepley     = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha
343649ef022SMatthew Knepley     = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha
344649ef022SMatthew Knepley     = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha)
345649ef022SMatthew Knepley */
346649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
347649ef022SMatthew Knepley {
348649ef022SMatthew Knepley   u[0] = 100.*PetscCosReal(time) + X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
349649ef022SMatthew Knepley   u[1] = 100.*PetscSinReal(time) + 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
350649ef022SMatthew Knepley   return 0;
351649ef022SMatthew Knepley }
352649ef022SMatthew Knepley static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
353649ef022SMatthew Knepley {
354649ef022SMatthew Knepley   u[0] = -100.*PetscSinReal(time);
355649ef022SMatthew Knepley   u[1] =  100.*PetscCosReal(time);
356649ef022SMatthew Knepley   return 0;
357649ef022SMatthew Knepley }
358649ef022SMatthew Knepley 
359649ef022SMatthew Knepley static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
360649ef022SMatthew Knepley {
361649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
362649ef022SMatthew Knepley   return 0;
363649ef022SMatthew Knepley }
364649ef022SMatthew Knepley 
365649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
366649ef022SMatthew Knepley {
367649ef022SMatthew Knepley   T[0] = 100.*PetscCosReal(time) + X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
368649ef022SMatthew Knepley   return 0;
369649ef022SMatthew Knepley }
370649ef022SMatthew Knepley static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
371649ef022SMatthew Knepley {
372649ef022SMatthew Knepley   T[0] = -100.*PetscSinReal(time);
373649ef022SMatthew Knepley   return 0;
374649ef022SMatthew Knepley }
375649ef022SMatthew Knepley 
376649ef022SMatthew Knepley static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
377649ef022SMatthew Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
378649ef022SMatthew Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
379649ef022SMatthew Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
380649ef022SMatthew Knepley {
381444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
382649ef022SMatthew Knepley 
383649ef022SMatthew 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];
384649ef022SMatthew 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];
385649ef022SMatthew Knepley }
386649ef022SMatthew Knepley 
387649ef022SMatthew Knepley static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
388649ef022SMatthew Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
389649ef022SMatthew Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
390649ef022SMatthew Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
391649ef022SMatthew Knepley {
392444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
393649ef022SMatthew Knepley 
394444129b9SMatthew G. Knepley   f0[0] -= 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;
395649ef022SMatthew Knepley }
396649ef022SMatthew Knepley 
397606d57d4SMatthew G. Knepley /*
398444129b9SMatthew G. Knepley   CASE: Taylor-Green vortex
399606d57d4SMatthew G. Knepley   In 2D we use exact solution:
400606d57d4SMatthew G. Knepley 
401606d57d4SMatthew G. Knepley     u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)
402606d57d4SMatthew G. Knepley     v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)
403606d57d4SMatthew G. Knepley     p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t)
404606d57d4SMatthew G. Knepley     T = t + x + y
405606d57d4SMatthew 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))  >
406606d57d4SMatthew G. Knepley     Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t)
407606d57d4SMatthew G. Knepley 
408606d57d4SMatthew G. Knepley   so that
409606d57d4SMatthew G. Knepley 
410606d57d4SMatthew 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
411606d57d4SMatthew G. Knepley 
412606d57d4SMatthew G. Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
413606d57d4SMatthew 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),
414606d57d4SMatthew 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)>
415606d57d4SMatthew 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),
416606d57d4SMatthew 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)>
417606d57d4SMatthew 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),
418606d57d4SMatthew 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)>
419606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
420606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
421606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
422606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
423606d57d4SMatthew 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),
424606d57d4SMatthew 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)>
425606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
426606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
427606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
428606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
429606d57d4SMatthew G. Knepley     + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
430606d57d4SMatthew G. Knepley        -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
431606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
432606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
433606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
434606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
435606d57d4SMatthew 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),
436606d57d4SMatthew 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)>
437606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
438606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
439606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
440606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
441606d57d4SMatthew G. Knepley     = < \pi cos(\pi(x - t)) cos(\pi(y - t)),
442606d57d4SMatthew G. Knepley         \pi sin(\pi(x - t)) sin(\pi(y - t))>
443606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)),
444606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0
445606d57d4SMatthew G. Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
446606d57d4SMatthew G. Knepley     = 1 + u \cdot <1, 1> - 0
447606d57d4SMatthew G. Knepley     = 1 + u + v
448606d57d4SMatthew G. Knepley */
449606d57d4SMatthew G. Knepley 
450606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
451606d57d4SMatthew G. Knepley {
452606d57d4SMatthew G. Knepley   u[0] = 1 - PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
453606d57d4SMatthew G. Knepley   u[1] = 1 + PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
454606d57d4SMatthew G. Knepley   return 0;
455606d57d4SMatthew G. Knepley }
456606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
457606d57d4SMatthew G. Knepley {
458606d57d4SMatthew G. Knepley   u[0] = -PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
459606d57d4SMatthew G. Knepley                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
460606d57d4SMatthew G. Knepley                   - 2*PETSC_PI*PetscCosReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
461606d57d4SMatthew G. Knepley   u[1] =  PETSC_PI*(PetscSinReal(PETSC_PI*(X[0]-time))*PetscSinReal(PETSC_PI*(X[1]-time))
462606d57d4SMatthew G. Knepley                   - PetscCosReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time))
463606d57d4SMatthew G. Knepley                   - 2*PETSC_PI*PetscSinReal(PETSC_PI*(X[0]-time))*PetscCosReal(PETSC_PI*(X[1]-time)))*PetscExpReal(-2*PETSC_PI*PETSC_PI*time);
464606d57d4SMatthew G. Knepley   return 0;
465606d57d4SMatthew G. Knepley }
466606d57d4SMatthew G. Knepley 
467606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
468606d57d4SMatthew G. Knepley {
469606d57d4SMatthew 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);
470606d57d4SMatthew G. Knepley   return 0;
471606d57d4SMatthew G. Knepley }
472606d57d4SMatthew G. Knepley 
473606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
474606d57d4SMatthew G. Knepley {
475606d57d4SMatthew G. Knepley   p[0] = PETSC_PI*(0.5*(PetscSinReal(2*PETSC_PI*(X[0]-time)) + PetscSinReal(2*PETSC_PI*(X[1]-time)))
476606d57d4SMatthew 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);
477606d57d4SMatthew G. Knepley   return 0;
478606d57d4SMatthew G. Knepley }
479606d57d4SMatthew G. Knepley 
480606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
481606d57d4SMatthew G. Knepley {
482606d57d4SMatthew G. Knepley   T[0] = time + X[0] + X[1];
483606d57d4SMatthew G. Knepley   return 0;
484606d57d4SMatthew G. Knepley }
485606d57d4SMatthew G. Knepley static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
486606d57d4SMatthew G. Knepley {
487606d57d4SMatthew G. Knepley   T[0] = 1.0;
488606d57d4SMatthew G. Knepley   return 0;
489606d57d4SMatthew G. Knepley }
490606d57d4SMatthew G. Knepley 
491606d57d4SMatthew G. Knepley static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
492606d57d4SMatthew G. Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
493606d57d4SMatthew G. Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
494606d57d4SMatthew G. Knepley                             PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
495606d57d4SMatthew G. Knepley {
496606d57d4SMatthew G. Knepley   PetscScalar vel[2];
497606d57d4SMatthew G. Knepley 
498606d57d4SMatthew G. Knepley   taylor_green_u(dim, t, X, Nf, vel, NULL);
499444129b9SMatthew G. Knepley   f0[0] -= 1.0 + vel[0] + vel[1];
500606d57d4SMatthew G. Knepley }
501606d57d4SMatthew G. Knepley 
502444129b9SMatthew G. Knepley /*
503444129b9SMatthew G. Knepley   CASE: Pipe flow
504444129b9SMatthew G. Knepley   Poiseuille flow, with the incoming fluid having a parabolic temperature profile and the side walls being held at T_in
505444129b9SMatthew G. Knepley 
506444129b9SMatthew G. Knepley     u = \Delta Re/(2 mu) y (1 - y)
507444129b9SMatthew G. Knepley     v = 0
508444129b9SMatthew G. Knepley     p = -\Delta x
509444129b9SMatthew G. Knepley     T = y (1 - y) + T_in
510444129b9SMatthew G. Knepley   rho = p^{th} / T
511444129b9SMatthew G. Knepley 
512444129b9SMatthew G. Knepley   so that
513444129b9SMatthew G. Knepley 
514444129b9SMatthew G. Knepley     \nabla \cdot u = 0 - 0 = 0
515444129b9SMatthew G. Knepley     grad u = \Delta Re/(2 mu) <<0, 0>, <1 - 2y, 0>>
516444129b9SMatthew G. Knepley     epsilon(u) = 1/2 (grad u + grad u^T) = \Delta Re/(4 mu) <<0, 1 - 2y>, <<1 - 2y, 0>>
517444129b9SMatthew G. Knepley     epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u)
518444129b9SMatthew G. Knepley     div epsilon'(u) = -\Delta Re/(2 mu) <1, 0>
519444129b9SMatthew G. Knepley 
520444129b9SMatthew G. Knepley   f = rho S du/dt + rho u \cdot \nabla u - 2\mu/Re div \epsilon'(u) + \nabla p + rho / F^2 \hat y
521444129b9SMatthew G. Knepley     = 0 + 0 - div (2\mu/Re \epsilon'(u) - pI) + rho / F^2 \hat y
522444129b9SMatthew G. Knepley     = -\Delta div <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> + rho / F^2 \hat y
523444129b9SMatthew G. Knepley     = \Delta <1, 0> - \Delta <1, 0> + rho/F^2 <0, 1>
524444129b9SMatthew G. Knepley     = rho/F^2 <0, 1>
525444129b9SMatthew G. Knepley 
526444129b9SMatthew G. Knepley   g = S rho_t + div (rho u)
527444129b9SMatthew G. Knepley     = 0 + rho div (u) + u . grad rho
528444129b9SMatthew G. Knepley     = 0 + 0 - pth u . grad T / T^2
529444129b9SMatthew G. Knepley     = 0
530444129b9SMatthew G. Knepley 
531444129b9SMatthew G. Knepley   Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T
532444129b9SMatthew G. Knepley     = 0 + c_p pth / T 0 + 2 k/Pe
533444129b9SMatthew G. Knepley     = 2 k/Pe
534444129b9SMatthew G. Knepley 
535444129b9SMatthew G. Knepley   The boundary conditions on the top and bottom are zero velocity and T_in temperature. The boundary term is
536444129b9SMatthew G. Knepley 
537444129b9SMatthew G. Knepley     (2\mu/Re \epsilon'(u) - p I) . n = \Delta <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> . n
538444129b9SMatthew G. Knepley 
539444129b9SMatthew G. Knepley   so that
540444129b9SMatthew G. Knepley 
541444129b9SMatthew G. Knepley     x = 0: \Delta <<0, (1 - 2y)/2>, <<(1 - 2y)/2, 0>> . <-1, 0> = <0, (2y - 1)/2>
542444129b9SMatthew G. Knepley     x = 1: \Delta <<1, (1 - 2y)/2>, <<(1 - 2y)/2, 1>> . <1, 0> = <1, (1 - 2y)/2>
543444129b9SMatthew G. Knepley */
544444129b9SMatthew G. Knepley static PetscErrorCode pipe_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
545444129b9SMatthew G. Knepley {
546444129b9SMatthew G. Knepley   Parameter *param = (Parameter *) ctx;
547444129b9SMatthew G. Knepley 
548444129b9SMatthew G. Knepley   u[0] = (0.5*param->Reynolds / param->mu) * X[1]*(1.0 - X[1]);
549444129b9SMatthew G. Knepley   u[1] = 0.0;
550444129b9SMatthew G. Knepley   return 0;
551444129b9SMatthew G. Knepley }
552444129b9SMatthew G. Knepley static PetscErrorCode pipe_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
553444129b9SMatthew G. Knepley {
554444129b9SMatthew G. Knepley   u[0] = 0.0;
555444129b9SMatthew G. Knepley   u[1] = 0.0;
556444129b9SMatthew G. Knepley   return 0;
557444129b9SMatthew G. Knepley }
558444129b9SMatthew G. Knepley 
559444129b9SMatthew G. Knepley static PetscErrorCode pipe_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
560444129b9SMatthew G. Knepley {
561444129b9SMatthew G. Knepley   p[0] = -X[0];
562444129b9SMatthew G. Knepley   return 0;
563444129b9SMatthew G. Knepley }
564444129b9SMatthew G. Knepley static PetscErrorCode pipe_p_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
565444129b9SMatthew G. Knepley {
566444129b9SMatthew G. Knepley   p[0] = 0.0;
567444129b9SMatthew G. Knepley   return 0;
568444129b9SMatthew G. Knepley }
569444129b9SMatthew G. Knepley 
570444129b9SMatthew G. Knepley static PetscErrorCode pipe_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
571444129b9SMatthew G. Knepley {
572444129b9SMatthew G. Knepley   Parameter *param = (Parameter *) ctx;
573444129b9SMatthew G. Knepley 
574444129b9SMatthew G. Knepley   T[0] = X[1]*(1.0 - X[1]) + param->T_in;
575444129b9SMatthew G. Knepley   return 0;
576444129b9SMatthew G. Knepley }
577444129b9SMatthew G. Knepley static PetscErrorCode pipe_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
578444129b9SMatthew G. Knepley {
579444129b9SMatthew G. Knepley   T[0] = 0.0;
580444129b9SMatthew G. Knepley   return 0;
581444129b9SMatthew G. Knepley }
582444129b9SMatthew G. Knepley 
583444129b9SMatthew G. Knepley static void f0_conduct_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
584444129b9SMatthew G. Knepley                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
585444129b9SMatthew G. Knepley                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
586444129b9SMatthew G. Knepley                               PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
587444129b9SMatthew G. Knepley {
588444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
589444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
590444129b9SMatthew G. Knepley   const PetscReal T_in = PetscRealPart(constants[T_IN]);
591444129b9SMatthew G. Knepley   const PetscReal rho  = p_th / (X[1]*(1. - X[1]) + T_in);
592444129b9SMatthew G. Knepley   const PetscInt  gd   = (PetscInt) PetscRealPart(constants[G_DIR]);
593444129b9SMatthew G. Knepley 
594444129b9SMatthew G. Knepley   f0[gd] -= rho/PetscSqr(F);
595444129b9SMatthew G. Knepley }
596444129b9SMatthew G. Knepley 
597444129b9SMatthew G. Knepley static void f0_conduct_bd_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
598444129b9SMatthew G. Knepley                                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
599444129b9SMatthew G. Knepley                                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
600444129b9SMatthew G. Knepley                                  PetscReal t, const PetscReal X[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
601444129b9SMatthew G. Knepley {
602444129b9SMatthew G. Knepley   PetscReal sigma[4] = {X[0], 0.5*(1. - 2.*X[1]), 0.5*(1. - 2.*X[1]), X[0]};
603444129b9SMatthew G. Knepley   PetscInt  d, e;
604444129b9SMatthew G. Knepley 
605444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
606444129b9SMatthew G. Knepley     for (e = 0; e < dim; ++e) {
607444129b9SMatthew G. Knepley       f0[d] -= sigma[d*dim + e] * n[e];
608444129b9SMatthew G. Knepley     }
609444129b9SMatthew G. Knepley   }
610444129b9SMatthew G. Knepley }
611444129b9SMatthew G. Knepley 
612444129b9SMatthew G. Knepley static void f0_conduct_pipe_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
613444129b9SMatthew G. Knepley                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
614444129b9SMatthew G. Knepley                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
615444129b9SMatthew G. Knepley                               PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
616444129b9SMatthew G. Knepley {
617444129b9SMatthew G. Knepley   f0[0] += 0.0;
618444129b9SMatthew G. Knepley }
619444129b9SMatthew G. Knepley 
620444129b9SMatthew G. Knepley static void f0_conduct_pipe_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
621444129b9SMatthew G. Knepley                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
622444129b9SMatthew G. Knepley                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
623444129b9SMatthew G. Knepley                               PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
624444129b9SMatthew G. Knepley {
625444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
626444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
627444129b9SMatthew G. Knepley 
628444129b9SMatthew G. Knepley   f0[0] -= 2*k/Pe;
629444129b9SMatthew G. Knepley }
630444129b9SMatthew G. Knepley 
631444129b9SMatthew G. Knepley /*      Physics Kernels      */
632444129b9SMatthew G. Knepley 
633649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
634649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
635649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
636649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
637649ef022SMatthew Knepley {
638649ef022SMatthew Knepley   PetscInt d;
639649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
640649ef022SMatthew Knepley }
641649ef022SMatthew Knepley 
642444129b9SMatthew G. Knepley /* -\frac{Sp^{th}}{T^2} \frac{\partial T}{\partial t} + \frac{p^{th}}{T} \nabla \cdot \vb{u} - \frac{p^{th}}{T^2} \vb{u} \cdot \nabla T */
643444129b9SMatthew G. Knepley static void f0_conduct_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
644444129b9SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
645444129b9SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
646444129b9SMatthew G. Knepley                          PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
647444129b9SMatthew G. Knepley {
648444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
649444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
650444129b9SMatthew G. Knepley   PetscInt        d;
651444129b9SMatthew G. Knepley 
652444129b9SMatthew G. Knepley   // -\frac{S p^{th}}{T^2} \frac{\partial T}{\partial t}
653444129b9SMatthew G. Knepley   f0[0] += -u_t[uOff[TEMP]] * S * p_th / PetscSqr(u[uOff[TEMP]]);
654444129b9SMatthew G. Knepley 
655444129b9SMatthew G. Knepley   // \frac{p^{th}}{T} \nabla \cdot \vb{u}
656444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
657444129b9SMatthew G. Knepley     f0[0] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + d*dim + d];
658444129b9SMatthew G. Knepley   }
659444129b9SMatthew G. Knepley 
660444129b9SMatthew G. Knepley   // - \frac{p^{th}}{T^2} \vb{u} \cdot \nabla T
661444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
662444129b9SMatthew G. Knepley     f0[0] -= p_th / (u[uOff[TEMP]] * u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d];
663444129b9SMatthew G. Knepley   }
664444129b9SMatthew G. Knepley 
665444129b9SMatthew G. Knepley   // Add in any fixed source term
666444129b9SMatthew G. Knepley   if (NfAux > 0) {
667444129b9SMatthew G. Knepley     f0[0] += a[aOff[MASS]];
668444129b9SMatthew G. Knepley   }
669444129b9SMatthew G. Knepley }
670444129b9SMatthew G. Knepley 
671444129b9SMatthew G. Knepley /* \vb{u}_t + \vb{u} \cdot \nabla\vb{u} */
672444129b9SMatthew G. Knepley static void f0_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
673444129b9SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
674444129b9SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
675444129b9SMatthew G. Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
676444129b9SMatthew G. Knepley {
677444129b9SMatthew G. Knepley   const PetscInt Nc = dim;
678444129b9SMatthew G. Knepley   PetscInt       c, d;
679444129b9SMatthew G. Knepley 
680444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
681444129b9SMatthew G. Knepley     /* \vb{u}_t */
682444129b9SMatthew G. Knepley     f0[c] += u_t[uOff[VEL] + c];
683444129b9SMatthew G. Knepley     /* \vb{u} \cdot \nabla\vb{u} */
684444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) f0[c] += u[uOff[VEL] + d]*u_x[uOff_x[VEL] + c*dim + d];
685444129b9SMatthew G. Knepley   }
686444129b9SMatthew G. Knepley }
687444129b9SMatthew G. Knepley 
688444129b9SMatthew G. Knepley /* \rho S \frac{\partial \vb{u}}{\partial t} + \rho \vb{u} \cdot \nabla \vb{u} + \rho \frac{\hat{\vb{z}}}{F^2} */
689444129b9SMatthew G. Knepley static void f0_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
690444129b9SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
691444129b9SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
692444129b9SMatthew G. Knepley                          PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
693444129b9SMatthew G. Knepley {
694444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
695444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
696444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
697444129b9SMatthew G. Knepley   const PetscReal rho  = p_th / PetscRealPart(u[uOff[TEMP]]);
698444129b9SMatthew G. Knepley   const PetscInt  gdir = (PetscInt) PetscRealPart(constants[G_DIR]);
699444129b9SMatthew G. Knepley   PetscInt        Nc   = dim;
700444129b9SMatthew G. Knepley   PetscInt        c, d;
701444129b9SMatthew G. Knepley 
702444129b9SMatthew G. Knepley   // \rho S \frac{\partial \vb{u}}{\partial t}
703444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
704444129b9SMatthew G. Knepley     f0[d] = rho * S * u_t[uOff[VEL] + d];
705444129b9SMatthew G. Knepley   }
706444129b9SMatthew G. Knepley 
707444129b9SMatthew G. Knepley   // \rho \vb{u} \cdot \nabla \vb{u}
708444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
709444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
710444129b9SMatthew G. Knepley       f0[c] += rho * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c*dim + d];
711444129b9SMatthew G. Knepley     }
712444129b9SMatthew G. Knepley   }
713444129b9SMatthew G. Knepley 
714444129b9SMatthew G. Knepley   // rho \hat{z}/F^2
715444129b9SMatthew G. Knepley   f0[gdir] += rho / (F*F);
716444129b9SMatthew G. Knepley 
717444129b9SMatthew G. Knepley   // Add in any fixed source term
718444129b9SMatthew G. Knepley   if (NfAux > 0) {
719444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
720444129b9SMatthew G. Knepley       f0[d] += a[aOff[MOMENTUM] + d];
721444129b9SMatthew G. Knepley     }
722444129b9SMatthew G. Knepley   }
723444129b9SMatthew G. Knepley }
724444129b9SMatthew G. Knepley 
725649ef022SMatthew Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
726649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
727649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
728649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
729649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
730649ef022SMatthew Knepley {
731444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
732649ef022SMatthew Knepley   const PetscInt  Nc = dim;
733649ef022SMatthew Knepley   PetscInt        c, d;
734649ef022SMatthew Knepley 
735649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
736649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
737649ef022SMatthew Knepley       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
738649ef022SMatthew Knepley     }
739649ef022SMatthew Knepley     f1[c*dim+c] -= u[uOff[1]];
740649ef022SMatthew Knepley   }
741649ef022SMatthew Knepley }
742649ef022SMatthew Knepley 
743444129b9SMatthew G. Knepley /* 2 \mu/Re (1/2 (\nabla \vb{u} + \nabla \vb{u}^T) - 1/3 (\nabla \cdot \vb{u}) I) - p I */
744444129b9SMatthew G. Knepley static void f1_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
745444129b9SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
746444129b9SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
747444129b9SMatthew G. Knepley                          PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
748444129b9SMatthew G. Knepley {
749444129b9SMatthew G. Knepley   const PetscReal Re    = PetscRealPart(constants[REYNOLDS]);
750444129b9SMatthew G. Knepley   const PetscReal mu    = PetscRealPart(constants[MU]);
751444129b9SMatthew G. Knepley   const PetscReal coef  = mu / Re;
752444129b9SMatthew G. Knepley   PetscReal       u_div = 0.0;
753444129b9SMatthew G. Knepley   const PetscInt  Nc    = dim;
754444129b9SMatthew G. Knepley   PetscInt        c, d;
755444129b9SMatthew G. Knepley 
756444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
757444129b9SMatthew G. Knepley     u_div += PetscRealPart(u_x[uOff_x[VEL] + c*dim + c]);
758444129b9SMatthew G. Knepley   }
759444129b9SMatthew G. Knepley 
760444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
761444129b9SMatthew G. Knepley     // 2 \mu/Re 1/2 (\nabla \vb{u} + \nabla \vb{u}^T
762444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
763444129b9SMatthew G. Knepley       f1[c*dim + d] += coef * (u_x[uOff_x[VEL] + c*dim + d] + u_x[uOff_x[VEL] + d*dim + c]);
764444129b9SMatthew G. Knepley     }
765444129b9SMatthew G. Knepley     // -2/3 \mu/Re (\nabla \cdot \vb{u}) I
766444129b9SMatthew G. Knepley     f1[c * dim + c] -= 2.0 * coef / 3.0 * u_div;
767444129b9SMatthew G. Knepley   }
768444129b9SMatthew G. Knepley 
769444129b9SMatthew G. Knepley   // -p I
770444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
771444129b9SMatthew G. Knepley     f1[c*dim + c] -= u[uOff[PRES]];
772444129b9SMatthew G. Knepley   }
773444129b9SMatthew G. Knepley }
774444129b9SMatthew G. Knepley 
775444129b9SMatthew G. Knepley /* T_t + \vb{u} \cdot \nabla T */
776444129b9SMatthew G. Knepley static void f0_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
777444129b9SMatthew G. Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
778444129b9SMatthew G. Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
779444129b9SMatthew G. Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
780444129b9SMatthew G. Knepley {
781444129b9SMatthew G. Knepley   PetscInt d;
782444129b9SMatthew G. Knepley 
783444129b9SMatthew G. Knepley   /* T_t */
784444129b9SMatthew G. Knepley   f0[0] += u_t[uOff[TEMP]];
785444129b9SMatthew G. Knepley   /* \vb{u} \cdot \nabla T */
786444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d];
787444129b9SMatthew G. Knepley }
788444129b9SMatthew G. Knepley 
789444129b9SMatthew G. Knepley /* \frac{C_p S p^{th}}{T} \frac{\partial T}{\partial t} + \frac{C_p p^{th}}{T} \vb{u} \cdot \nabla T */
790444129b9SMatthew G. Knepley static void f0_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
791444129b9SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
792444129b9SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
793444129b9SMatthew G. Knepley                          PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
794444129b9SMatthew G. Knepley {
795444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
796444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
797444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
798444129b9SMatthew G. Knepley   PetscInt        d;
799444129b9SMatthew G. Knepley 
800444129b9SMatthew G. Knepley   // \frac{C_p S p^{th}}{T} \frac{\partial T}{\partial t}
801444129b9SMatthew G. Knepley   f0[0] = c_p * S * p_th / u[uOff[TEMP]] * u_t[uOff[TEMP]];
802444129b9SMatthew G. Knepley 
803444129b9SMatthew G. Knepley   // \frac{C_p p^{th}}{T} \vb{u} \cdot \nabla T
804444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
805444129b9SMatthew G. Knepley     f0[0] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d];
806444129b9SMatthew G. Knepley   }
807444129b9SMatthew G. Knepley 
808444129b9SMatthew G. Knepley   // Add in any fixed source term
809444129b9SMatthew G. Knepley   if (NfAux > 0) {
810444129b9SMatthew G. Knepley     f0[0] += a[aOff[ENERGY]];
811444129b9SMatthew G. Knepley   }
812444129b9SMatthew G. Knepley }
813444129b9SMatthew G. Knepley 
814649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
815649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
816649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
817649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
818649ef022SMatthew Knepley {
819444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
820649ef022SMatthew Knepley   PetscInt        d;
821444129b9SMatthew G. Knepley 
822649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
823649ef022SMatthew Knepley }
824649ef022SMatthew Knepley 
825444129b9SMatthew G. Knepley /* \frac{k}{Pe} \nabla T */
826444129b9SMatthew G. Knepley static void f1_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
827444129b9SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
828444129b9SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
829444129b9SMatthew G. Knepley                          PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
830444129b9SMatthew G. Knepley {
831444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
832444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
833444129b9SMatthew G. Knepley   PetscInt        d;
834444129b9SMatthew G. Knepley 
835444129b9SMatthew G. Knepley   // \frac{k}{Pe} \nabla T
836444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
837444129b9SMatthew G. Knepley     f1[d] = k / Pe * u_x[uOff_x[TEMP] + d];
838444129b9SMatthew G. Knepley   }
839444129b9SMatthew G. Knepley }
840444129b9SMatthew G. Knepley 
841649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
842649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
843649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
844649ef022SMatthew Knepley                  PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
845649ef022SMatthew Knepley {
846649ef022SMatthew Knepley   PetscInt d;
847649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
848649ef022SMatthew Knepley }
849649ef022SMatthew Knepley 
850649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
851649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
852649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
853649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
854649ef022SMatthew Knepley {
855649ef022SMatthew Knepley   PetscInt c, d;
856649ef022SMatthew Knepley   const PetscInt  Nc = dim;
857649ef022SMatthew Knepley 
858649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d*dim+d] = u_tShift;
859649ef022SMatthew Knepley 
860649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
861649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
862649ef022SMatthew Knepley       g0[c*Nc+d] += u_x[ c*Nc+d];
863649ef022SMatthew Knepley     }
864649ef022SMatthew Knepley   }
865649ef022SMatthew Knepley }
866649ef022SMatthew Knepley 
867649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
868649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
869649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
870649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
871649ef022SMatthew Knepley {
872649ef022SMatthew Knepley   PetscInt NcI = dim;
873649ef022SMatthew Knepley   PetscInt NcJ = dim;
874649ef022SMatthew Knepley   PetscInt c, d, e;
875649ef022SMatthew Knepley 
876649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
877649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
878649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
879649ef022SMatthew Knepley         if (c == d) {
880649ef022SMatthew Knepley           g1[(c*NcJ+d)*dim+e] += u[e];
881649ef022SMatthew Knepley         }
882649ef022SMatthew Knepley       }
883649ef022SMatthew Knepley     }
884649ef022SMatthew Knepley   }
885649ef022SMatthew Knepley }
886649ef022SMatthew Knepley 
887444129b9SMatthew G. Knepley static void g0_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
888444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
889444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
890444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
891444129b9SMatthew G. Knepley {
892444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
893444129b9SMatthew G. Knepley   PetscInt        d;
894444129b9SMatthew G. Knepley 
895444129b9SMatthew G. Knepley   // - \phi_i \frac{p^{th}}{T^2} \frac{\partial T}{\partial x_c} \psi_{j, u_c}
896444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
897444129b9SMatthew G. Knepley     g0[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u_x[uOff_x[TEMP] + d];
898444129b9SMatthew G. Knepley   }
899444129b9SMatthew G. Knepley }
900444129b9SMatthew G. Knepley 
901444129b9SMatthew G. Knepley static void g1_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
902444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
903444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
904444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
905444129b9SMatthew G. Knepley {
906444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
907444129b9SMatthew G. Knepley   PetscInt        d;
908444129b9SMatthew G. Knepley 
909444129b9SMatthew G. Knepley   // \phi_i \frac{p^{th}}{T} \frac{\partial \psi_{u_c,j}}{\partial x_c}
910444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
911444129b9SMatthew G. Knepley     g1[d * dim + d] = p_th / u[uOff[TEMP]];
912444129b9SMatthew G. Knepley   }
913444129b9SMatthew G. Knepley }
914444129b9SMatthew G. Knepley 
915444129b9SMatthew G. Knepley static void g0_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
916444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
917444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
918444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
919444129b9SMatthew G. Knepley {
920444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
921444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
922444129b9SMatthew G. Knepley   PetscInt        d;
923444129b9SMatthew G. Knepley 
924444129b9SMatthew G. Knepley   // - \phi_i \frac{S p^{th}}{T^2} \psi_j
925444129b9SMatthew G. Knepley   g0[0] -= S * p_th / PetscSqr(u[uOff[TEMP]]) * u_tShift;
926444129b9SMatthew G. Knepley   // \phi_i 2 \frac{S p^{th}}{T^3} T_t \psi_j
927444129b9SMatthew G. Knepley   g0[0] += 2.0 * S * p_th / PetscPowScalarInt(u[uOff[TEMP]], 3) * u_t[uOff[TEMP]];
928444129b9SMatthew G. Knepley   // \phi_i \frac{p^{th}}{T^2} \left( - \nabla \cdot \vb{u} \psi_j + \frac{2}{T} \vb{u} \cdot \nabla T \psi_j \right)
929444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
930444129b9SMatthew G. Knepley     g0[0] += p_th / PetscSqr(u[uOff[TEMP]]) * (-u_x[uOff_x[VEL] + d * dim + d] + 2.0 / u[uOff[TEMP]] * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]);
931444129b9SMatthew G. Knepley   }
932444129b9SMatthew G. Knepley }
933444129b9SMatthew G. Knepley 
934444129b9SMatthew G. Knepley static void g1_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
935444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
936444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
937444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
938444129b9SMatthew G. Knepley {
939444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
940444129b9SMatthew G. Knepley   PetscInt        d;
941444129b9SMatthew G. Knepley 
942444129b9SMatthew G. Knepley   // - \phi_i \frac{p^{th}}{T^2} \vb{u} \cdot \nabla \psi_j
943444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
944444129b9SMatthew G. Knepley     g1[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d];
945444129b9SMatthew G. Knepley   }
946444129b9SMatthew G. Knepley }
947444129b9SMatthew G. Knepley 
948649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
949649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
950649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
951649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
952649ef022SMatthew Knepley {
953649ef022SMatthew Knepley   PetscInt d;
954649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
955649ef022SMatthew Knepley }
956649ef022SMatthew Knepley 
957649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
958649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
959649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
960649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
961649ef022SMatthew Knepley {
962444129b9SMatthew G. Knepley    const PetscReal nu = PetscRealPart(constants[NU]);
963649ef022SMatthew Knepley    const PetscInt  Nc = dim;
964649ef022SMatthew Knepley    PetscInt        c, d;
965649ef022SMatthew Knepley 
966649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
967649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
968606d57d4SMatthew G. Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += nu;
969606d57d4SMatthew G. Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += nu;
970649ef022SMatthew Knepley     }
971649ef022SMatthew Knepley   }
972649ef022SMatthew Knepley }
973649ef022SMatthew Knepley 
974444129b9SMatthew G. Knepley static void g0_conduct_vT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
975444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
976444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
977444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
978444129b9SMatthew G. Knepley {
979444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
980444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
981444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
982444129b9SMatthew G. Knepley   const PetscInt  gdir = (PetscInt) PetscRealPart(constants[G_DIR]);
983444129b9SMatthew G. Knepley   const PetscInt  Nc = dim;
984444129b9SMatthew G. Knepley   PetscInt        c, d;
985444129b9SMatthew G. Knepley 
986444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vb{u}_t \frac{p^{th} S}{T^2} \psi_j
987444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
988444129b9SMatthew G. Knepley     g0[d] -= p_th * S / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[VEL] + d];
989444129b9SMatthew G. Knepley   }
990444129b9SMatthew G. Knepley 
991444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vb{u} \cdot \nabla \vb{u} \frac{p^{th}}{T^2} \psi_j
992444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
993444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
994444129b9SMatthew G. Knepley       g0[c] -= p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c * dim + d];
995444129b9SMatthew G. Knepley     }
996444129b9SMatthew G. Knepley   }
997444129b9SMatthew G. Knepley 
998444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vu{z} \frac{p^{th}}{T^2 F^2} \psi_j
999444129b9SMatthew G. Knepley   g0[gdir] -= p_th / PetscSqr(u[uOff[TEMP]] * F);
1000444129b9SMatthew G. Knepley }
1001444129b9SMatthew G. Knepley 
1002444129b9SMatthew G. Knepley static void g0_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1003444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1004444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1005444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1006444129b9SMatthew G. Knepley {
1007444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
1008444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1009444129b9SMatthew G. Knepley   const PetscInt  Nc   = dim;
1010444129b9SMatthew G. Knepley   PetscInt        c, d;
1011444129b9SMatthew G. Knepley 
1012444129b9SMatthew G. Knepley   // \vb{\phi}_i \cdot S \rho \psi_j
1013444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1014444129b9SMatthew G. Knepley     g0[d * dim + d] = S * p_th / u[uOff[TEMP]] * u_tShift;
1015444129b9SMatthew G. Knepley   }
1016444129b9SMatthew G. Knepley 
1017444129b9SMatthew G. Knepley   // \phi^c_i \cdot \rho \frac{\partial u^c}{\partial x^d} \psi^d_j
1018444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1019444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1020444129b9SMatthew G. Knepley       g0[c * Nc + d] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + c * Nc + d];
1021444129b9SMatthew G. Knepley     }
1022444129b9SMatthew G. Knepley   }
1023444129b9SMatthew G. Knepley }
1024444129b9SMatthew G. Knepley 
1025444129b9SMatthew G. Knepley static void g1_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1026444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1027444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1028444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
1029444129b9SMatthew G. Knepley {
1030444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1031444129b9SMatthew G. Knepley   const PetscInt  NcI = dim;
1032444129b9SMatthew G. Knepley   const PetscInt  NcJ = dim;
1033444129b9SMatthew G. Knepley   PetscInt        c, d, e;
1034444129b9SMatthew G. Knepley 
1035444129b9SMatthew G. Knepley   // \phi^c_i \rho u^e \frac{\partial \psi^d_j}{\partial x^e}
1036444129b9SMatthew G. Knepley   for (c = 0; c < NcI; ++c) {
1037444129b9SMatthew G. Knepley     for (d = 0; d < NcJ; ++d) {
1038444129b9SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
1039444129b9SMatthew G. Knepley         if (c == d) {
1040444129b9SMatthew G. Knepley           g1[(c * NcJ + d) * dim + e] += p_th / u[uOff[TEMP]] * u[uOff[VEL] + e];
1041444129b9SMatthew G. Knepley         }
1042444129b9SMatthew G. Knepley       }
1043444129b9SMatthew G. Knepley     }
1044444129b9SMatthew G. Knepley   }
1045444129b9SMatthew G. Knepley }
1046444129b9SMatthew G. Knepley 
1047444129b9SMatthew G. Knepley static void g3_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1048444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1049444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1050444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
1051444129b9SMatthew G. Knepley {
1052444129b9SMatthew G. Knepley   const PetscReal Re = PetscRealPart(constants[REYNOLDS]);
1053444129b9SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[MU]);
1054444129b9SMatthew G. Knepley   const PetscInt  Nc = dim;
1055444129b9SMatthew G. Knepley   PetscInt        c, d;
1056444129b9SMatthew G. Knepley 
1057444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
1058444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
1059444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^c_i}{\partial x^d}
1060444129b9SMatthew G. Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += mu / Re;  // gradU
1061444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^d_i}{\partial x^c}
1062444129b9SMatthew G. Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += mu / Re;  // gradU transpose
1063444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} -2/3 \mu/Re \frac{\partial \psi^d_i}{\partial x^c}
1064444129b9SMatthew G. Knepley       g3[((c * Nc + d) * dim + c) * dim + d] -= 2.0 / 3.0 * mu / Re;
1065444129b9SMatthew G. Knepley     }
1066444129b9SMatthew G. Knepley   }
1067444129b9SMatthew G. Knepley }
1068444129b9SMatthew G. Knepley 
1069444129b9SMatthew G. Knepley static void g2_conduct_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1070444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1071444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1072444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
1073444129b9SMatthew G. Knepley {
1074444129b9SMatthew G. Knepley   PetscInt d;
1075444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1076444129b9SMatthew G. Knepley     g2[d * dim + d] = -1.0;
1077444129b9SMatthew G. Knepley   }
1078444129b9SMatthew G. Knepley }
1079444129b9SMatthew G. Knepley 
1080649ef022SMatthew Knepley static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1081649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1082649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1083649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1084649ef022SMatthew Knepley {
1085a712f3bbSMatthew G. Knepley   g0[0] = u_tShift;
1086649ef022SMatthew Knepley }
1087649ef022SMatthew Knepley 
1088649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1089649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1090649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1091649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1092649ef022SMatthew Knepley {
1093649ef022SMatthew Knepley   PetscInt d;
1094649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
1095649ef022SMatthew Knepley }
1096649ef022SMatthew Knepley 
1097649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1098649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1099649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1100649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
1101649ef022SMatthew Knepley {
1102649ef022SMatthew Knepley   PetscInt d;
1103649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
1104649ef022SMatthew Knepley }
1105649ef022SMatthew Knepley 
1106649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1107649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1108649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1109649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
1110649ef022SMatthew Knepley {
1111444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
1112649ef022SMatthew Knepley   PetscInt        d;
1113649ef022SMatthew Knepley 
1114649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
1115649ef022SMatthew Knepley }
1116649ef022SMatthew Knepley 
1117444129b9SMatthew G. Knepley static void g0_conduct_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1118444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1119444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1120444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1121444129b9SMatthew G. Knepley {
1122444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1123444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
1124444129b9SMatthew G. Knepley   PetscInt        d;
1125444129b9SMatthew G. Knepley 
1126444129b9SMatthew G. Knepley   // \phi_i \frac{C_p p^{th}}{T} \nabla T \cdot \psi_j
1127444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1128444129b9SMatthew G. Knepley     g0[d] = c_p * p_th / u[uOff[TEMP]] * u_x[uOff_x[TEMP] + d];
1129444129b9SMatthew G. Knepley   }
1130444129b9SMatthew G. Knepley }
1131444129b9SMatthew G. Knepley 
1132444129b9SMatthew G. Knepley static void g0_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1133444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1134444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1135444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
1136444129b9SMatthew G. Knepley {
1137444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
1138444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1139444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
1140444129b9SMatthew G. Knepley   PetscInt        d;
1141444129b9SMatthew G. Knepley 
1142444129b9SMatthew G. Knepley   // \psi_i C_p S p^{th}\T \psi_{j}
1143444129b9SMatthew G. Knepley   g0[0] += c_p * S * p_th / u[uOff[TEMP]] * u_tShift;
1144444129b9SMatthew G. Knepley   // - \phi_i C_p S p^{th}/T^2 T_t \psi_j
1145444129b9SMatthew G. Knepley   g0[0] -= c_p * S * p_th / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[TEMP]];
1146444129b9SMatthew G. Knepley   // - \phi_i C_p p^{th}/T^2 \vb{u} \cdot \nabla T \psi_j
1147444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1148444129b9SMatthew G. Knepley     g0[0] -= c_p * p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d];
1149444129b9SMatthew G. Knepley   }
1150444129b9SMatthew G. Knepley }
1151444129b9SMatthew G. Knepley 
1152444129b9SMatthew G. Knepley static void g1_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1153444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1154444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1155444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
1156444129b9SMatthew G. Knepley {
1157444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1158444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
1159444129b9SMatthew G. Knepley   PetscInt        d;
1160444129b9SMatthew G. Knepley 
1161444129b9SMatthew G. Knepley   // \phi_i C_p p^{th}/T \vb{u} \cdot \nabla \psi_j
1162444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1163444129b9SMatthew G. Knepley     g1[d] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d];
1164444129b9SMatthew G. Knepley   }
1165444129b9SMatthew G. Knepley }
1166444129b9SMatthew G. Knepley 
1167444129b9SMatthew G. Knepley static void g3_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1168444129b9SMatthew G. Knepley                           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1169444129b9SMatthew G. Knepley                           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1170444129b9SMatthew G. Knepley                           PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
1171444129b9SMatthew G. Knepley {
1172444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
1173444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
1174444129b9SMatthew G. Knepley   PetscInt        d;
1175444129b9SMatthew G. Knepley 
1176444129b9SMatthew G. Knepley   // \nabla \phi_i \frac{k}{Pe} \nabla \phi_j
1177444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
1178444129b9SMatthew G. Knepley     g3[d * dim + d] = k / Pe;
1179444129b9SMatthew G. Knepley   }
1180444129b9SMatthew G. Knepley }
1181444129b9SMatthew G. Knepley 
1182649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
1183649ef022SMatthew Knepley {
1184444129b9SMatthew G. Knepley   PetscInt       mod, sol;
1185649ef022SMatthew Knepley   PetscErrorCode ierr;
1186649ef022SMatthew Knepley 
1187649ef022SMatthew Knepley   PetscFunctionBeginUser;
1188444129b9SMatthew G. Knepley   options->modType      = MOD_INCOMPRESSIBLE;
1189649ef022SMatthew Knepley   options->solType      = SOL_QUADRATIC;
1190444129b9SMatthew G. Knepley   options->hasNullSpace = PETSC_TRUE;
1191649ef022SMatthew Knepley 
1192649ef022SMatthew Knepley   ierr = PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");CHKERRQ(ierr);
1193444129b9SMatthew G. Knepley   mod = options->modType;
1194444129b9SMatthew G. Knepley   ierr = PetscOptionsEList("-mod_type", "The model type", "ex76.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL);CHKERRQ(ierr);
1195444129b9SMatthew G. Knepley   options->modType = (ModType) mod;
1196649ef022SMatthew Knepley   sol = options->solType;
1197444129b9SMatthew G. Knepley   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex76.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
1198649ef022SMatthew Knepley   options->solType = (SolType) sol;
1199*1e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1200649ef022SMatthew Knepley   PetscFunctionReturn(0);
1201649ef022SMatthew Knepley }
1202649ef022SMatthew Knepley 
1203444129b9SMatthew G. Knepley static PetscErrorCode SetupParameters(DM dm, AppCtx *user)
1204649ef022SMatthew Knepley {
1205649ef022SMatthew Knepley   PetscBag       bag;
1206649ef022SMatthew Knepley   Parameter     *p;
1207444129b9SMatthew G. Knepley   PetscReal      dir;
1208444129b9SMatthew G. Knepley   PetscInt       dim;
1209649ef022SMatthew Knepley   PetscErrorCode ierr;
1210649ef022SMatthew Knepley 
1211649ef022SMatthew Knepley   PetscFunctionBeginUser;
1212444129b9SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1213444129b9SMatthew G. Knepley   dir  = (PetscReal) (dim-1);
1214649ef022SMatthew Knepley   /* setup PETSc parameter bag */
1215649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
1216649ef022SMatthew Knepley   ierr = PetscBagSetName(user->bag, "par", "Low Mach flow parameters");CHKERRQ(ierr);
1217649ef022SMatthew Knepley   bag  = user->bag;
1218444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->Strouhal, 1.0, "S",     "Strouhal number");CHKERRQ(ierr);
1219444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->Froude,   1.0, "Fr",    "Froude number");CHKERRQ(ierr);
1220444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->Reynolds, 1.0, "Re",    "Reynolds number");CHKERRQ(ierr);
1221444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->Peclet,   1.0, "Pe",    "Peclet number");CHKERRQ(ierr);
1222444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->p_th,     1.0, "p_th",  "Thermodynamic pressure");CHKERRQ(ierr);
1223444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->mu,       1.0, "mu",    "Dynamic viscosity");CHKERRQ(ierr);
1224649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->nu,       1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
1225444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->c_p,      1.0, "c_p",   "Specific heat at constant pressure");CHKERRQ(ierr);
1226444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->k,        1.0, "k",     "Thermal conductivity");CHKERRQ(ierr);
1227649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->alpha,    1.0, "alpha", "Thermal diffusivity");CHKERRQ(ierr);
1228649ef022SMatthew Knepley   ierr = PetscBagRegisterReal(bag, &p->T_in,     1.0, "T_in",  "Inlet temperature");CHKERRQ(ierr);
1229444129b9SMatthew G. Knepley   ierr = PetscBagRegisterReal(bag, &p->g_dir,    dir, "g_dir", "Gravity direction");CHKERRQ(ierr);
1230649ef022SMatthew Knepley   PetscFunctionReturn(0);
1231649ef022SMatthew Knepley }
1232649ef022SMatthew Knepley 
1233649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
1234649ef022SMatthew Knepley {
1235649ef022SMatthew Knepley   PetscErrorCode ierr;
1236649ef022SMatthew Knepley 
1237649ef022SMatthew Knepley   PetscFunctionBeginUser;
123830602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
123930602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1240649ef022SMatthew Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
1241649ef022SMatthew Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
1242649ef022SMatthew Knepley   PetscFunctionReturn(0);
1243649ef022SMatthew Knepley }
1244649ef022SMatthew Knepley 
1245444129b9SMatthew G. Knepley static PetscErrorCode UniformBoundaryConditions(DM dm, DMLabel label, PetscSimplePointFunc exactFuncs[], PetscSimplePointFunc exactFuncs_t[], AppCtx *user)
1246444129b9SMatthew G. Knepley {
1247444129b9SMatthew G. Knepley   PetscDS        ds;
1248444129b9SMatthew G. Knepley   PetscInt       id;
1249444129b9SMatthew G. Knepley   void          *ctx;
1250444129b9SMatthew G. Knepley   PetscErrorCode ierr;
1251444129b9SMatthew G. Knepley 
1252444129b9SMatthew G. Knepley   PetscFunctionBeginUser;
1253444129b9SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1254444129b9SMatthew G. Knepley   ierr = PetscBagGetData(user->bag, &ctx);CHKERRQ(ierr);
1255444129b9SMatthew G. Knepley   id   = 3;
1256444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall velocity",    label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1257444129b9SMatthew G. Knepley   id   = 1;
1258444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1259444129b9SMatthew G. Knepley   id   = 2;
1260444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "right wall velocity",  label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1261444129b9SMatthew G. Knepley   id   = 4;
1262444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall velocity",   label, 1, &id, VEL, 0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1263444129b9SMatthew G. Knepley   id   = 3;
1264444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall temp",    label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1265444129b9SMatthew G. Knepley   id   = 1;
1266444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1267444129b9SMatthew G. Knepley   id   = 2;
1268444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "right wall temp",  label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1269444129b9SMatthew G. Knepley   id   = 4;
1270444129b9SMatthew G. Knepley   ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall temp",   label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1271444129b9SMatthew G. Knepley   PetscFunctionReturn(0);
1272444129b9SMatthew G. Knepley }
1273444129b9SMatthew G. Knepley 
1274649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
1275649ef022SMatthew Knepley {
127645480ffeSMatthew G. Knepley   PetscSimplePointFunc exactFuncs[3];
127745480ffeSMatthew G. Knepley   PetscSimplePointFunc exactFuncs_t[3];
1278444129b9SMatthew G. Knepley   PetscDS              ds;
1279444129b9SMatthew G. Knepley   PetscWeakForm        wf;
128045480ffeSMatthew G. Knepley   DMLabel              label;
1281649ef022SMatthew Knepley   Parameter           *ctx;
1282444129b9SMatthew G. Knepley   PetscInt             id, bd;
1283649ef022SMatthew Knepley   PetscErrorCode       ierr;
1284649ef022SMatthew Knepley 
1285649ef022SMatthew Knepley   PetscFunctionBeginUser;
128645480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
1287444129b9SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1288444129b9SMatthew G. Knepley   ierr = PetscDSGetWeakForm(ds, &wf);CHKERRQ(ierr);
1289444129b9SMatthew G. Knepley 
1290444129b9SMatthew G. Knepley   switch(user->modType) {
1291444129b9SMatthew G. Knepley     case MOD_INCOMPRESSIBLE:
1292444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, VEL,  f0_v, f1_v);CHKERRQ(ierr);
1293444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, PRES, f0_q, NULL);CHKERRQ(ierr);
1294444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, TEMP, f0_w, f1_w);CHKERRQ(ierr);
1295444129b9SMatthew G. Knepley 
1296444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, VEL,  VEL,  g0_vu, g1_vu, NULL,  g3_vu);CHKERRQ(ierr);
1297444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, VEL,  PRES, NULL,  NULL,  g2_vp, NULL);CHKERRQ(ierr);
1298444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, PRES, VEL,  NULL,  g1_qu, NULL,  NULL);CHKERRQ(ierr);
1299444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, TEMP, VEL,  g0_wu, NULL,  NULL,  NULL);CHKERRQ(ierr);
1300444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, TEMP, TEMP, g0_wT, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
1301444129b9SMatthew G. Knepley 
1302649ef022SMatthew Knepley       switch(user->solType) {
1303649ef022SMatthew Knepley       case SOL_QUADRATIC:
1304444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL,  0, 1, f0_quadratic_v, 0, NULL);CHKERRQ(ierr);
1305444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_quadratic_w, 0, NULL);CHKERRQ(ierr);
1306649ef022SMatthew Knepley 
1307444129b9SMatthew G. Knepley         exactFuncs[VEL]    = quadratic_u;
1308444129b9SMatthew G. Knepley         exactFuncs[PRES]   = quadratic_p;
1309444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = quadratic_T;
1310444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = quadratic_u_t;
1311444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = NULL;
1312444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = quadratic_T_t;
1313444129b9SMatthew G. Knepley 
1314444129b9SMatthew G. Knepley         ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr);
1315649ef022SMatthew Knepley         break;
1316649ef022SMatthew Knepley       case SOL_CUBIC:
1317444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL,  0, 1, f0_cubic_v, 0, NULL);CHKERRQ(ierr);
1318444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_w, 0, NULL);CHKERRQ(ierr);
1319649ef022SMatthew Knepley 
1320444129b9SMatthew G. Knepley         exactFuncs[VEL]    = cubic_u;
1321444129b9SMatthew G. Knepley         exactFuncs[PRES]   = cubic_p;
1322444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = cubic_T;
1323444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = cubic_u_t;
1324444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = NULL;
1325444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = cubic_T_t;
1326444129b9SMatthew G. Knepley 
1327444129b9SMatthew G. Knepley         ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr);
1328649ef022SMatthew Knepley         break;
1329649ef022SMatthew Knepley       case SOL_CUBIC_TRIG:
1330444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL,  0, 1, f0_cubic_trig_v, 0, NULL);CHKERRQ(ierr);
1331444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_trig_w, 0, NULL);CHKERRQ(ierr);
1332649ef022SMatthew Knepley 
1333444129b9SMatthew G. Knepley         exactFuncs[VEL]    = cubic_trig_u;
1334444129b9SMatthew G. Knepley         exactFuncs[PRES]   = cubic_trig_p;
1335444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = cubic_trig_T;
1336444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = cubic_trig_u_t;
1337444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = NULL;
1338444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = cubic_trig_T_t;
1339444129b9SMatthew G. Knepley 
1340444129b9SMatthew G. Knepley         ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr);
1341649ef022SMatthew Knepley         break;
1342606d57d4SMatthew G. Knepley       case SOL_TAYLOR_GREEN:
1343444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_taylor_green_w, 0, NULL);CHKERRQ(ierr);
1344606d57d4SMatthew G. Knepley 
1345444129b9SMatthew G. Knepley         exactFuncs[VEL]    = taylor_green_u;
1346444129b9SMatthew G. Knepley         exactFuncs[PRES]   = taylor_green_p;
1347444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = taylor_green_T;
1348444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = taylor_green_u_t;
1349444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = taylor_green_p_t;
1350444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = taylor_green_T_t;
1351444129b9SMatthew G. Knepley 
1352444129b9SMatthew G. Knepley         ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr);
1353606d57d4SMatthew G. Knepley         break;
1354444129b9SMatthew G. Knepley        default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
1355649ef022SMatthew Knepley       }
1356444129b9SMatthew G. Knepley       break;
1357444129b9SMatthew G. Knepley     case MOD_CONDUCTING:
1358444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, VEL,  f0_conduct_v, f1_conduct_v);CHKERRQ(ierr);
1359444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, PRES, f0_conduct_q, NULL);CHKERRQ(ierr);
1360444129b9SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, TEMP, f0_conduct_w, f1_conduct_w);CHKERRQ(ierr);
1361649ef022SMatthew Knepley 
1362444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, VEL,  VEL,  g0_conduct_vu, g1_conduct_vu, NULL,          g3_conduct_vu);CHKERRQ(ierr);
1363444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, VEL,  PRES, NULL,          NULL,          g2_conduct_vp, NULL);CHKERRQ(ierr);
1364444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, VEL,  TEMP, g0_conduct_vT, NULL,          NULL,          NULL);CHKERRQ(ierr);
1365444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, PRES, VEL,  g0_conduct_qu, g1_conduct_qu, NULL,          NULL);CHKERRQ(ierr);
1366444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, PRES, TEMP, g0_conduct_qT, g1_conduct_qT, NULL,          NULL);CHKERRQ(ierr);
1367444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, TEMP, VEL,  g0_conduct_wu, NULL,          NULL,          NULL);CHKERRQ(ierr);
1368444129b9SMatthew G. Knepley       ierr = PetscDSSetJacobian(ds, TEMP, TEMP, g0_conduct_wT, g1_conduct_wT, NULL,          g3_conduct_wT);CHKERRQ(ierr);
1369649ef022SMatthew Knepley 
1370444129b9SMatthew G. Knepley       switch(user->solType) {
1371444129b9SMatthew G. Knepley       case SOL_QUADRATIC:
1372444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL,  0, 1, f0_conduct_quadratic_v, 0, NULL);CHKERRQ(ierr);
1373444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_quadratic_q, 0, NULL);CHKERRQ(ierr);
1374444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_quadratic_w, 0, NULL);CHKERRQ(ierr);
1375444129b9SMatthew G. Knepley 
1376444129b9SMatthew G. Knepley         exactFuncs[VEL]    = quadratic_u;
1377444129b9SMatthew G. Knepley         exactFuncs[PRES]   = quadratic_p;
1378444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = quadratic_T;
1379444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = quadratic_u_t;
1380444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = NULL;
1381444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = quadratic_T_t;
1382444129b9SMatthew G. Knepley 
1383444129b9SMatthew G. Knepley         ierr = UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user);CHKERRQ(ierr);
1384444129b9SMatthew G. Knepley         break;
1385444129b9SMatthew G. Knepley       case SOL_PIPE:
1386444129b9SMatthew G. Knepley         user->hasNullSpace = PETSC_FALSE;
1387444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL,  0, 1, f0_conduct_pipe_v, 0, NULL);CHKERRQ(ierr);
1388444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_pipe_q, 0, NULL);CHKERRQ(ierr);
1389444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_pipe_w, 0, NULL);CHKERRQ(ierr);
1390444129b9SMatthew G. Knepley 
1391444129b9SMatthew G. Knepley         exactFuncs[VEL]    = pipe_u;
1392444129b9SMatthew G. Knepley         exactFuncs[PRES]   = pipe_p;
1393444129b9SMatthew G. Knepley         exactFuncs[TEMP]   = pipe_T;
1394444129b9SMatthew G. Knepley         exactFuncs_t[VEL]  = pipe_u_t;
1395444129b9SMatthew G. Knepley         exactFuncs_t[PRES] = pipe_p_t;
1396444129b9SMatthew G. Knepley         exactFuncs_t[TEMP] = pipe_T_t;
1397444129b9SMatthew G. Knepley 
1398444129b9SMatthew G. Knepley         ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
1399444129b9SMatthew G. Knepley         id   = 2;
1400444129b9SMatthew G. Knepley         ierr = DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr);
1401444129b9SMatthew G. Knepley         ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
1402444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL);CHKERRQ(ierr);
1403444129b9SMatthew G. Knepley         id   = 4;
1404444129b9SMatthew G. Knepley         ierr = DMAddBoundary(dm, DM_BC_NATURAL, "left wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd);CHKERRQ(ierr);
1405444129b9SMatthew G. Knepley         ierr = PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
1406444129b9SMatthew G. Knepley         ierr = PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL);CHKERRQ(ierr);
1407444129b9SMatthew G. Knepley         id   = 4;
1408444129b9SMatthew G. Knepley         ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "left wall temperature",   label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1409444129b9SMatthew G. Knepley         id   = 3;
1410444129b9SMatthew G. Knepley         ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall velocity",       label, 1, &id, VEL,  0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1411444129b9SMatthew G. Knepley         ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "top wall temperature",    label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1412444129b9SMatthew G. Knepley         id   = 1;
1413444129b9SMatthew G. Knepley         ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall velocity",    label, 1, &id, VEL,  0, NULL, (void (*)(void)) exactFuncs[VEL], (void (*)(void)) exactFuncs_t[VEL], ctx, NULL);CHKERRQ(ierr);
1414444129b9SMatthew G. Knepley         ierr = PetscDSAddBoundary(ds, DM_BC_ESSENTIAL, "bottom wall temperature", label, 1, &id, TEMP, 0, NULL, (void (*)(void)) exactFuncs[TEMP], (void (*)(void)) exactFuncs_t[TEMP], ctx, NULL);CHKERRQ(ierr);
1415444129b9SMatthew G. Knepley         break;
1416444129b9SMatthew G. Knepley        default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
1417444129b9SMatthew G. Knepley       }
1418444129b9SMatthew G. Knepley       break;
1419444129b9SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%D)", solTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType);
1420444129b9SMatthew G. Knepley   }
1421649ef022SMatthew Knepley   /* Setup constants */
1422649ef022SMatthew Knepley   {
1423649ef022SMatthew Knepley     Parameter  *param;
1424444129b9SMatthew G. Knepley     PetscScalar constants[12];
1425649ef022SMatthew Knepley 
1426649ef022SMatthew Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
1427649ef022SMatthew Knepley 
1428444129b9SMatthew G. Knepley     constants[STROUHAL] = param->Strouhal;
1429444129b9SMatthew G. Knepley     constants[FROUDE]   = param->Froude;
1430444129b9SMatthew G. Knepley     constants[REYNOLDS] = param->Reynolds;
1431444129b9SMatthew G. Knepley     constants[PECLET]   = param->Peclet;
1432444129b9SMatthew G. Knepley     constants[P_TH]     = param->p_th;
1433444129b9SMatthew G. Knepley     constants[MU]       = param->mu;
1434444129b9SMatthew G. Knepley     constants[NU]       = param->nu;
1435444129b9SMatthew G. Knepley     constants[C_P]      = param->c_p;
1436444129b9SMatthew G. Knepley     constants[K]        = param->k;
1437444129b9SMatthew G. Knepley     constants[ALPHA]    = param->alpha;
1438444129b9SMatthew G. Knepley     constants[T_IN]     = param->T_in;
1439444129b9SMatthew G. Knepley     constants[G_DIR]    = param->g_dir;
1440444129b9SMatthew G. Knepley     ierr = PetscDSSetConstants(ds, 12, constants);CHKERRQ(ierr);
1441649ef022SMatthew Knepley   }
1442649ef022SMatthew Knepley 
1443444129b9SMatthew G. Knepley   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
1444444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, VEL,  exactFuncs[VEL],  ctx);CHKERRQ(ierr);
1445444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, PRES, exactFuncs[PRES], ctx);CHKERRQ(ierr);
1446444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, TEMP, exactFuncs[TEMP], ctx);CHKERRQ(ierr);
1447444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(ds, VEL,  exactFuncs_t[VEL],  ctx);CHKERRQ(ierr);
1448444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(ds, PRES, exactFuncs_t[PRES], ctx);CHKERRQ(ierr);
1449444129b9SMatthew G. Knepley   ierr = PetscDSSetExactSolutionTimeDerivative(ds, TEMP, exactFuncs_t[TEMP], ctx);CHKERRQ(ierr);
1450649ef022SMatthew Knepley   PetscFunctionReturn(0);
1451649ef022SMatthew Knepley }
1452649ef022SMatthew Knepley 
1453649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
1454649ef022SMatthew Knepley {
1455649ef022SMatthew Knepley   DM              cdm   = dm;
1456a712f3bbSMatthew G. Knepley   PetscFE         fe[3], fediv;
1457649ef022SMatthew Knepley   Parameter      *param;
1458649ef022SMatthew Knepley   DMPolytopeType  ct;
1459649ef022SMatthew Knepley   PetscInt        dim, cStart;
1460649ef022SMatthew Knepley   PetscBool       simplex;
1461649ef022SMatthew Knepley   PetscErrorCode  ierr;
1462649ef022SMatthew Knepley 
1463649ef022SMatthew Knepley   PetscFunctionBeginUser;
1464649ef022SMatthew Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1465649ef022SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
1466649ef022SMatthew Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
1467649ef022SMatthew Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
1468649ef022SMatthew Knepley   /* Create finite element */
1469a712f3bbSMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
1470649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
1471649ef022SMatthew Knepley 
1472a712f3bbSMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
1473649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
1474649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
1475649ef022SMatthew Knepley 
1476a712f3bbSMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
1477649ef022SMatthew Knepley   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
1478649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
1479649ef022SMatthew Knepley 
1480a712f3bbSMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "div_", PETSC_DEFAULT, &fediv);CHKERRQ(ierr);
1481a712f3bbSMatthew G. Knepley   ierr = PetscFECopyQuadrature(fe[0], fediv);CHKERRQ(ierr);
1482a712f3bbSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) fediv, "divergence");CHKERRQ(ierr);
1483a712f3bbSMatthew G. Knepley 
1484649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
1485444129b9SMatthew G. Knepley   ierr = DMSetField(dm, VEL,  NULL, (PetscObject) fe[VEL]);CHKERRQ(ierr);
1486444129b9SMatthew G. Knepley   ierr = DMSetField(dm, PRES, NULL, (PetscObject) fe[PRES]);CHKERRQ(ierr);
1487444129b9SMatthew G. Knepley   ierr = DMSetField(dm, TEMP, NULL, (PetscObject) fe[TEMP]);CHKERRQ(ierr);
1488649ef022SMatthew Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
1489649ef022SMatthew Knepley   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
1490649ef022SMatthew Knepley   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
1491649ef022SMatthew Knepley   while (cdm) {
1492649ef022SMatthew Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
1493649ef022SMatthew Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
1494649ef022SMatthew Knepley   }
1495444129b9SMatthew G. Knepley   ierr = PetscFEDestroy(&fe[VEL]);CHKERRQ(ierr);
1496444129b9SMatthew G. Knepley   ierr = PetscFEDestroy(&fe[PRES]);CHKERRQ(ierr);
1497444129b9SMatthew G. Knepley   ierr = PetscFEDestroy(&fe[TEMP]);CHKERRQ(ierr);
1498649ef022SMatthew Knepley 
1499a712f3bbSMatthew G. Knepley   ierr = DMClone(dm, &user->dmCell);CHKERRQ(ierr);
1500a712f3bbSMatthew G. Knepley   ierr = DMSetField(user->dmCell, 0, NULL, (PetscObject) fediv);CHKERRQ(ierr);
1501a712f3bbSMatthew G. Knepley   ierr = DMCreateDS(user->dmCell);CHKERRQ(ierr);
1502a712f3bbSMatthew G. Knepley   ierr = PetscFEDestroy(&fediv);CHKERRQ(ierr);
1503a712f3bbSMatthew G. Knepley 
1504444129b9SMatthew G. Knepley   if (user->hasNullSpace) {
1505649ef022SMatthew Knepley     PetscObject  pressure;
1506649ef022SMatthew Knepley     MatNullSpace nullspacePres;
1507649ef022SMatthew Knepley 
1508444129b9SMatthew G. Knepley     ierr = DMGetField(dm, PRES, NULL, &pressure);CHKERRQ(ierr);
1509649ef022SMatthew Knepley     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
1510649ef022SMatthew Knepley     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
1511649ef022SMatthew Knepley     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
1512649ef022SMatthew Knepley   }
1513649ef022SMatthew Knepley   PetscFunctionReturn(0);
1514649ef022SMatthew Knepley }
1515649ef022SMatthew Knepley 
1516649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
1517649ef022SMatthew Knepley {
1518649ef022SMatthew Knepley   Vec              vec;
1519649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
1520649ef022SMatthew Knepley   PetscErrorCode   ierr;
1521649ef022SMatthew Knepley 
1522649ef022SMatthew Knepley   PetscFunctionBeginUser;
1523444129b9SMatthew G. Knepley   if (ofield != PRES) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index %D, not %D", PRES, ofield);
1524649ef022SMatthew Knepley   funcs[nfield] = constant;
1525649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
1526649ef022SMatthew Knepley   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
1527649ef022SMatthew Knepley   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
1528649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
1529649ef022SMatthew Knepley   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
1530649ef022SMatthew Knepley   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
1531649ef022SMatthew Knepley   ierr = VecDestroy(&vec);CHKERRQ(ierr);
1532649ef022SMatthew Knepley   PetscFunctionReturn(0);
1533649ef022SMatthew Knepley }
1534649ef022SMatthew Knepley 
1535649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u)
1536649ef022SMatthew Knepley {
1537649ef022SMatthew Knepley   DM             dm;
1538444129b9SMatthew G. Knepley   AppCtx        *user;
1539649ef022SMatthew Knepley   MatNullSpace   nullsp;
1540649ef022SMatthew Knepley   PetscErrorCode ierr;
1541649ef022SMatthew Knepley 
1542649ef022SMatthew Knepley   PetscFunctionBegin;
1543649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
1544444129b9SMatthew G. Knepley   ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr);
1545444129b9SMatthew G. Knepley   if (!user->hasNullSpace) PetscFunctionReturn(0);
1546649ef022SMatthew Knepley   ierr = CreatePressureNullSpace(dm, 1, 1, &nullsp);CHKERRQ(ierr);
1547649ef022SMatthew Knepley   ierr = MatNullSpaceRemove(nullsp, u);CHKERRQ(ierr);
1548649ef022SMatthew Knepley   ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr);
1549649ef022SMatthew Knepley   PetscFunctionReturn(0);
1550649ef022SMatthew Knepley }
1551649ef022SMatthew Knepley 
1552649ef022SMatthew Knepley /* Make the discrete pressure discretely divergence free */
1553649ef022SMatthew Knepley static PetscErrorCode RemoveDiscretePressureNullspace(TS ts)
1554649ef022SMatthew Knepley {
1555649ef022SMatthew Knepley   Vec            u;
1556649ef022SMatthew Knepley   PetscErrorCode ierr;
1557649ef022SMatthew Knepley 
1558649ef022SMatthew Knepley   PetscFunctionBegin;
1559649ef022SMatthew Knepley   ierr = TSGetSolution(ts, &u);CHKERRQ(ierr);
1560649ef022SMatthew Knepley   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
1561649ef022SMatthew Knepley   PetscFunctionReturn(0);
1562649ef022SMatthew Knepley }
1563649ef022SMatthew Knepley 
1564649ef022SMatthew Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u)
1565649ef022SMatthew Knepley {
1566444129b9SMatthew G. Knepley   AppCtx        *user;
1567649ef022SMatthew Knepley   DM             dm;
1568649ef022SMatthew Knepley   PetscReal      t;
1569649ef022SMatthew Knepley   PetscErrorCode ierr;
1570649ef022SMatthew Knepley 
1571649ef022SMatthew Knepley   PetscFunctionBegin;
1572649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
1573649ef022SMatthew Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
1574649ef022SMatthew Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
1575444129b9SMatthew G. Knepley   ierr = DMGetApplicationContext(dm, (void **) &user);CHKERRQ(ierr);
1576649ef022SMatthew Knepley   ierr = RemoveDiscretePressureNullspace_Private(ts, u);CHKERRQ(ierr);
1577649ef022SMatthew Knepley   PetscFunctionReturn(0);
1578649ef022SMatthew Knepley }
1579649ef022SMatthew Knepley 
1580a712f3bbSMatthew G. Knepley static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux,
1581a712f3bbSMatthew G. Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
1582a712f3bbSMatthew G. Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
1583a712f3bbSMatthew G. Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar divu[])
1584a712f3bbSMatthew G. Knepley {
1585a712f3bbSMatthew G. Knepley   PetscInt d;
1586a712f3bbSMatthew G. Knepley 
1587a712f3bbSMatthew G. Knepley   divu[0] = 0.;
1588a712f3bbSMatthew G. Knepley   for (d = 0; d < dim; ++d) divu[0] += u_x[d*dim+d];
1589a712f3bbSMatthew G. Knepley }
1590a712f3bbSMatthew G. Knepley 
1591649ef022SMatthew Knepley static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
1592649ef022SMatthew Knepley {
1593649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1594649ef022SMatthew Knepley   void            *ctxs[3];
1595a712f3bbSMatthew G. Knepley   PetscPointFunc   diagnostics[1] = {divergence};
1596a712f3bbSMatthew G. Knepley   DM               dm, dmCell = ((AppCtx *) ctx)->dmCell;
1597649ef022SMatthew Knepley   PetscDS          ds;
1598a712f3bbSMatthew G. Knepley   Vec              v, divu;
1599a712f3bbSMatthew G. Knepley   PetscReal        ferrors[3], massFlux;
1600649ef022SMatthew Knepley   PetscInt         f;
1601649ef022SMatthew Knepley   PetscErrorCode   ierr;
1602649ef022SMatthew Knepley 
1603649ef022SMatthew Knepley   PetscFunctionBeginUser;
1604649ef022SMatthew Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
1605649ef022SMatthew Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
1606649ef022SMatthew Knepley 
1607649ef022SMatthew Knepley   for (f = 0; f < 3; ++f) {ierr = PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]);CHKERRQ(ierr);}
1608649ef022SMatthew Knepley   ierr = DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors);CHKERRQ(ierr);
1609a712f3bbSMatthew G. Knepley   ierr = DMGetGlobalVector(dmCell, &divu);CHKERRQ(ierr);
1610a712f3bbSMatthew G. Knepley   ierr = DMProjectField(dmCell, crtime, u, diagnostics, INSERT_VALUES, divu);CHKERRQ(ierr);
1611a712f3bbSMatthew G. Knepley   ierr = VecViewFromOptions(divu, NULL, "-divu_vec_view");CHKERRQ(ierr);
1612a712f3bbSMatthew G. Knepley   ierr = VecNorm(divu, NORM_2, &massFlux);CHKERRQ(ierr);
1613a712f3bbSMatthew G. Knepley   ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g, %2.3g] ||div u||: %2.3g\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1], (double) ferrors[2], (double) massFlux);CHKERRQ(ierr);
1614649ef022SMatthew Knepley 
1615649ef022SMatthew Knepley   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
1616649ef022SMatthew Knepley 
1617649ef022SMatthew Knepley   ierr = DMGetGlobalVector(dm, &v);CHKERRQ(ierr);
1618a712f3bbSMatthew G. Knepley   ierr = DMProjectFunction(dm, crtime, exactFuncs, ctxs, INSERT_ALL_VALUES, v);CHKERRQ(ierr);
1619649ef022SMatthew Knepley   ierr = PetscObjectSetName((PetscObject) v, "Exact Solution");CHKERRQ(ierr);
1620649ef022SMatthew Knepley   ierr = VecViewFromOptions(v, NULL, "-exact_vec_view");CHKERRQ(ierr);
1621649ef022SMatthew Knepley   ierr = DMRestoreGlobalVector(dm, &v);CHKERRQ(ierr);
1622649ef022SMatthew Knepley 
1623a712f3bbSMatthew G. Knepley   ierr = VecViewFromOptions(divu, NULL, "-div_vec_view");CHKERRQ(ierr);
1624a712f3bbSMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmCell, &divu);CHKERRQ(ierr);
1625a712f3bbSMatthew G. Knepley 
1626649ef022SMatthew Knepley   PetscFunctionReturn(0);
1627649ef022SMatthew Knepley }
1628649ef022SMatthew Knepley 
1629649ef022SMatthew Knepley int main(int argc, char **argv)
1630649ef022SMatthew Knepley {
1631649ef022SMatthew Knepley   DM              dm;   /* problem definition */
1632649ef022SMatthew Knepley   TS              ts;   /* timestepper */
1633649ef022SMatthew Knepley   Vec             u;    /* solution */
1634649ef022SMatthew Knepley   AppCtx          user; /* user-defined work context */
1635649ef022SMatthew Knepley   PetscReal       t;
1636649ef022SMatthew Knepley   PetscErrorCode  ierr;
1637649ef022SMatthew Knepley 
1638649ef022SMatthew Knepley   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
1639649ef022SMatthew Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
1640649ef022SMatthew Knepley   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
1641649ef022SMatthew Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
1642444129b9SMatthew G. Knepley   ierr = SetupParameters(dm, &user);CHKERRQ(ierr);
1643444129b9SMatthew G. Knepley   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
1644649ef022SMatthew Knepley   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
1645649ef022SMatthew Knepley   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
1646649ef022SMatthew Knepley   /* Setup problem */
1647649ef022SMatthew Knepley   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
1648649ef022SMatthew Knepley   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
1649649ef022SMatthew Knepley 
1650649ef022SMatthew Knepley   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
1651a712f3bbSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
1652444129b9SMatthew G. Knepley   if (user.hasNullSpace) {ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);}
1653649ef022SMatthew Knepley 
1654649ef022SMatthew Knepley   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user);CHKERRQ(ierr);
1655649ef022SMatthew Knepley   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user);CHKERRQ(ierr);
1656649ef022SMatthew Knepley   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user);CHKERRQ(ierr);
1657649ef022SMatthew Knepley   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
1658649ef022SMatthew Knepley   ierr = TSSetPreStep(ts, RemoveDiscretePressureNullspace);CHKERRQ(ierr);
1659649ef022SMatthew Knepley   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
1660649ef022SMatthew Knepley 
1661649ef022SMatthew Knepley   ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr); /* Must come after SetFromOptions() */
1662649ef022SMatthew Knepley   ierr = SetInitialConditions(ts, u);CHKERRQ(ierr);
1663649ef022SMatthew Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
1664649ef022SMatthew Knepley   ierr = DMSetOutputSequenceNumber(dm, 0, t);CHKERRQ(ierr);
1665649ef022SMatthew Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
1666649ef022SMatthew Knepley   ierr = TSMonitorSet(ts, MonitorError, &user, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1667649ef022SMatthew Knepley 
1668649ef022SMatthew Knepley   ierr = TSSolve(ts, u);CHKERRQ(ierr);
1669649ef022SMatthew Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
1670649ef022SMatthew Knepley 
1671649ef022SMatthew Knepley   ierr = VecDestroy(&u);CHKERRQ(ierr);
1672a712f3bbSMatthew G. Knepley   ierr = DMDestroy(&user.dmCell);CHKERRQ(ierr);
1673649ef022SMatthew Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1674649ef022SMatthew Knepley   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1675649ef022SMatthew Knepley   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
1676649ef022SMatthew Knepley   ierr = PetscFinalize();
1677649ef022SMatthew Knepley   return ierr;
1678649ef022SMatthew Knepley }
1679649ef022SMatthew Knepley 
1680649ef022SMatthew Knepley /*TEST
1681649ef022SMatthew Knepley 
1682444129b9SMatthew G. Knepley   testset:
1683649ef022SMatthew Knepley     requires: triangle !single
1684444129b9SMatthew G. Knepley     args: -dm_plex_separate_marker \
1685a712f3bbSMatthew G. Knepley           -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \
1686444129b9SMatthew G. Knepley           -snes_error_if_not_converged -snes_convergence_test correct_pressure \
1687649ef022SMatthew Knepley           -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
1688444129b9SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \
1689444129b9SMatthew G. Knepley           -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
1690649ef022SMatthew Knepley             -fieldsplit_0_pc_type lu \
1691649ef022SMatthew Knepley             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1692649ef022SMatthew Knepley 
1693444129b9SMatthew G. Knepley     test:
1694444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1
1695444129b9SMatthew G. Knepley       args: -sol_type quadratic \
1696444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1697444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1698444129b9SMatthew G. Knepley 
1699649ef022SMatthew Knepley     test:
1700649ef022SMatthew Knepley       # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0]
1701649ef022SMatthew Knepley       suffix: 2d_tri_p2_p1_p1_tconv
1702444129b9SMatthew G. Knepley       args: -sol_type cubic_trig \
1703649ef022SMatthew Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1704444129b9SMatthew G. Knepley             -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1
1705649ef022SMatthew Knepley 
1706649ef022SMatthew Knepley     test:
1707649ef022SMatthew Knepley       # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9]
1708649ef022SMatthew Knepley       suffix: 2d_tri_p2_p1_p1_sconv
1709444129b9SMatthew G. Knepley       args: -sol_type cubic \
1710649ef022SMatthew Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1711444129b9SMatthew G. Knepley             -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1
1712649ef022SMatthew Knepley 
1713649ef022SMatthew Knepley     test:
1714649ef022SMatthew Knepley       suffix: 2d_tri_p3_p2_p2
1715444129b9SMatthew G. Knepley       args: -sol_type cubic \
1716649ef022SMatthew Knepley             -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
1717444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1718649ef022SMatthew Knepley 
1719606d57d4SMatthew G. Knepley     test:
1720606d57d4SMatthew G. Knepley       # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1]
1721606d57d4SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_tg_sconv
1722444129b9SMatthew G. Knepley       args: -sol_type taylor_green \
1723606d57d4SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1724444129b9SMatthew G. Knepley             -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1
1725606d57d4SMatthew G. Knepley 
1726606d57d4SMatthew G. Knepley     test:
1727606d57d4SMatthew G. Knepley       # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2]
1728606d57d4SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_tg_tconv
1729444129b9SMatthew G. Knepley       args: -sol_type taylor_green \
1730606d57d4SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1731444129b9SMatthew G. Knepley             -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1
1732444129b9SMatthew G. Knepley 
1733444129b9SMatthew G. Knepley   testset:
1734444129b9SMatthew G. Knepley     requires: triangle !single
1735444129b9SMatthew G. Knepley     args: -dm_plex_separate_marker -mod_type conducting \
1736a712f3bbSMatthew G. Knepley           -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \
1737444129b9SMatthew G. Knepley           -snes_error_if_not_converged -snes_max_linear_solve_fail 5 \
1738444129b9SMatthew G. Knepley           -ksp_type fgmres -ksp_max_it 2 -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
1739444129b9SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \
1740444129b9SMatthew G. Knepley           -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
1741606d57d4SMatthew G. Knepley             -fieldsplit_0_pc_type lu \
1742606d57d4SMatthew G. Knepley             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1743606d57d4SMatthew G. Knepley 
1744444129b9SMatthew G. Knepley     test:
1745444129b9SMatthew G. Knepley       # At this resolution, the rhs is inconsistent on some Newton steps
1746444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_conduct
1747444129b9SMatthew G. Knepley       args: -sol_type quadratic \
1748444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1749444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
1750444129b9SMatthew G. Knepley             -pc_fieldsplit_schur_precondition full \
1751444129b9SMatthew G. Knepley               -fieldsplit_pressure_ksp_max_it 2 -fieldsplit_pressure_pc_type svd
1752444129b9SMatthew G. Knepley 
1753444129b9SMatthew G. Knepley     test:
1754444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p2_conduct_pipe
1755444129b9SMatthew G. Knepley       args: -sol_type pipe \
1756444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 2 \
1757444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1758444129b9SMatthew G. Knepley 
1759649ef022SMatthew Knepley TEST*/
1760