xref: /petsc/src/ts/tutorials/ex76.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
234e6a9dc0SMatthew Knepley To look at nonlinear solver convergence, use
244e6a9dc0SMatthew Knepley 
254e6a9dc0SMatthew Knepley   -dm_refine <k> -ts_max_steps 1 \
264e6a9dc0SMatthew Knepley   -ts_view -ts_monitor -snes_monitor -snes_converged_reason -ksp_converged_reason -fieldsplit_pressure_ksp_converged_reason
274e6a9dc0SMatthew Knepley 
28444129b9SMatthew G. Knepley [1] https://ubchrest.github.io/ablate/content/formulations/lowMachFlow/
29444129b9SMatthew G. Knepley [2] https://github.com/UBCHREST/ablate/blob/main/ablateCore/flow/lowMachFlow.c
30444129b9SMatthew 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.
31649ef022SMatthew Knepley F*/
32649ef022SMatthew Knepley 
33649ef022SMatthew Knepley #include <petscdmplex.h>
34649ef022SMatthew Knepley #include <petscsnes.h>
35649ef022SMatthew Knepley #include <petscts.h>
36649ef022SMatthew Knepley #include <petscds.h>
37649ef022SMatthew Knepley #include <petscbag.h>
38649ef022SMatthew Knepley 
39*9371c9d4SSatish Balay typedef enum {
40*9371c9d4SSatish Balay   MOD_INCOMPRESSIBLE,
41*9371c9d4SSatish Balay   MOD_CONDUCTING,
42*9371c9d4SSatish Balay   NUM_MOD_TYPES
43*9371c9d4SSatish Balay } ModType;
44444129b9SMatthew G. Knepley const char *modTypes[NUM_MOD_TYPES + 1] = {"incompressible", "conducting", "unknown"};
45444129b9SMatthew G. Knepley 
46*9371c9d4SSatish Balay typedef enum {
47*9371c9d4SSatish Balay   SOL_QUADRATIC,
48*9371c9d4SSatish Balay   SOL_CUBIC,
49*9371c9d4SSatish Balay   SOL_CUBIC_TRIG,
50*9371c9d4SSatish Balay   SOL_TAYLOR_GREEN,
51*9371c9d4SSatish Balay   SOL_PIPE,
52*9371c9d4SSatish Balay   SOL_PIPE_WIGGLY,
53*9371c9d4SSatish Balay   NUM_SOL_TYPES
54*9371c9d4SSatish Balay } SolType;
55367970cfSMatthew G. Knepley const char *solTypes[NUM_SOL_TYPES + 1] = {"quadratic", "cubic", "cubic_trig", "taylor_green", "pipe", "pipe_wiggly", "unknown"};
56444129b9SMatthew G. Knepley 
57444129b9SMatthew G. Knepley /* Fields */
58444129b9SMatthew G. Knepley const PetscInt VEL      = 0;
59444129b9SMatthew G. Knepley const PetscInt PRES     = 1;
60444129b9SMatthew G. Knepley const PetscInt TEMP     = 2;
61444129b9SMatthew G. Knepley /* Sources */
62444129b9SMatthew G. Knepley const PetscInt MOMENTUM = 0;
63444129b9SMatthew G. Knepley const PetscInt MASS     = 1;
64444129b9SMatthew G. Knepley const PetscInt ENERGY   = 2;
65444129b9SMatthew G. Knepley /* Constants */
66444129b9SMatthew G. Knepley const PetscInt STROUHAL = 0;
67444129b9SMatthew G. Knepley const PetscInt FROUDE   = 1;
68444129b9SMatthew G. Knepley const PetscInt REYNOLDS = 2;
69444129b9SMatthew G. Knepley const PetscInt PECLET   = 3;
70444129b9SMatthew G. Knepley const PetscInt P_TH     = 4;
71444129b9SMatthew G. Knepley const PetscInt MU       = 5;
72444129b9SMatthew G. Knepley const PetscInt NU       = 6;
73444129b9SMatthew G. Knepley const PetscInt C_P      = 7;
74444129b9SMatthew G. Knepley const PetscInt K        = 8;
75444129b9SMatthew G. Knepley const PetscInt ALPHA    = 9;
76444129b9SMatthew G. Knepley const PetscInt T_IN     = 10;
77444129b9SMatthew G. Knepley const PetscInt G_DIR    = 11;
78367970cfSMatthew G. Knepley const PetscInt EPSILON  = 12;
79649ef022SMatthew Knepley 
80649ef022SMatthew Knepley typedef struct {
81444129b9SMatthew G. Knepley   PetscReal Strouhal; /* Strouhal number */
82444129b9SMatthew G. Knepley   PetscReal Froude;   /* Froude number */
83444129b9SMatthew G. Knepley   PetscReal Reynolds; /* Reynolds number */
84444129b9SMatthew G. Knepley   PetscReal Peclet;   /* Peclet number */
85444129b9SMatthew G. Knepley   PetscReal p_th;     /* Thermodynamic pressure */
86444129b9SMatthew G. Knepley   PetscReal mu;       /* Dynamic viscosity */
87649ef022SMatthew Knepley   PetscReal nu;       /* Kinematic viscosity */
88444129b9SMatthew G. Knepley   PetscReal c_p;      /* Specific heat at constant pressure */
89444129b9SMatthew G. Knepley   PetscReal k;        /* Thermal conductivity */
90649ef022SMatthew Knepley   PetscReal alpha;    /* Thermal diffusivity */
91649ef022SMatthew Knepley   PetscReal T_in;     /* Inlet temperature */
92444129b9SMatthew G. Knepley   PetscReal g_dir;    /* Gravity direction */
93367970cfSMatthew G. Knepley   PetscReal epsilon;  /* Strength of perturbation */
94649ef022SMatthew Knepley } Parameter;
95649ef022SMatthew Knepley 
96649ef022SMatthew Knepley typedef struct {
97649ef022SMatthew Knepley   /* Problem definition */
98649ef022SMatthew Knepley   PetscBag  bag;          /* Holds problem parameters */
99444129b9SMatthew G. Knepley   ModType   modType;      /* Model type */
100649ef022SMatthew Knepley   SolType   solType;      /* MMS solution type */
101444129b9SMatthew G. Knepley   PetscBool hasNullSpace; /* Problem has the constant null space for pressure */
102a712f3bbSMatthew G. Knepley   /* Flow diagnostics */
103a712f3bbSMatthew G. Knepley   DM        dmCell; /* A DM with piecewise constant discretization */
104649ef022SMatthew Knepley } AppCtx;
105649ef022SMatthew Knepley 
106*9371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
107649ef022SMatthew Knepley   PetscInt d;
108649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
109649ef022SMatthew Knepley   return 0;
110649ef022SMatthew Knepley }
111649ef022SMatthew Knepley 
112*9371c9d4SSatish Balay static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
113649ef022SMatthew Knepley   PetscInt d;
114649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
115649ef022SMatthew Knepley   return 0;
116649ef022SMatthew Knepley }
117649ef022SMatthew Knepley 
118649ef022SMatthew Knepley /*
119649ef022SMatthew Knepley   CASE: quadratic
120649ef022SMatthew Knepley   In 2D we use exact solution:
121649ef022SMatthew Knepley 
122649ef022SMatthew Knepley     u = t + x^2 + y^2
123649ef022SMatthew Knepley     v = t + 2x^2 - 2xy
124649ef022SMatthew Knepley     p = x + y - 1
125444129b9SMatthew G. Knepley     T = t + x + y + 1
126649ef022SMatthew 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>
127649ef022SMatthew Knepley     Q = 1 + 2t + 3x^2 - 2xy + y^2
128649ef022SMatthew Knepley 
129649ef022SMatthew Knepley   so that
130649ef022SMatthew Knepley 
131649ef022SMatthew Knepley     \nabla \cdot u = 2x - 2x = 0
132649ef022SMatthew Knepley 
133649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
134649ef022SMatthew Knepley     = <1, 1> + <t + x^2 + y^2, t + 2x^2 - 2xy> . <<2x, 4x - 2y>, <2y, -2x>> - \nu <4, 4> + <1, 1>
135649ef022SMatthew 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>
136649ef022SMatthew 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>
137649ef022SMatthew Knepley 
138649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
139649ef022SMatthew Knepley     = 1 + <t + x^2 + y^2, t + 2x^2 - 2xy> . <1, 1> - \alpha 0
140649ef022SMatthew Knepley     = 1 + 2t + 3x^2 - 2xy + y^2
141649ef022SMatthew Knepley */
142649ef022SMatthew Knepley 
143*9371c9d4SSatish Balay static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
144649ef022SMatthew Knepley   u[0] = time + X[0] * X[0] + X[1] * X[1];
145649ef022SMatthew Knepley   u[1] = time + 2.0 * X[0] * X[0] - 2.0 * X[0] * X[1];
146649ef022SMatthew Knepley   return 0;
147649ef022SMatthew Knepley }
148*9371c9d4SSatish Balay static PetscErrorCode quadratic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
149649ef022SMatthew Knepley   u[0] = 1.0;
150649ef022SMatthew Knepley   u[1] = 1.0;
151649ef022SMatthew Knepley   return 0;
152649ef022SMatthew Knepley }
153649ef022SMatthew Knepley 
154*9371c9d4SSatish Balay static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
155649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
156649ef022SMatthew Knepley   return 0;
157649ef022SMatthew Knepley }
158649ef022SMatthew Knepley 
159*9371c9d4SSatish Balay static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
160444129b9SMatthew G. Knepley   T[0] = time + X[0] + X[1] + 1.0;
161649ef022SMatthew Knepley   return 0;
162649ef022SMatthew Knepley }
163*9371c9d4SSatish Balay static PetscErrorCode quadratic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
164649ef022SMatthew Knepley   T[0] = 1.0;
165649ef022SMatthew Knepley   return 0;
166649ef022SMatthew Knepley }
167649ef022SMatthew Knepley 
168*9371c9d4SSatish Balay static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
169444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
170649ef022SMatthew Knepley 
171444129b9SMatthew 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;
172444129b9SMatthew 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;
173649ef022SMatthew Knepley }
174649ef022SMatthew Knepley 
175*9371c9d4SSatish Balay static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
176444129b9SMatthew G. Knepley   f0[0] -= 2 * t + 1 + 3 * X[0] * X[0] - 2 * X[0] * X[1] + X[1] * X[1];
177444129b9SMatthew G. Knepley }
178444129b9SMatthew G. Knepley 
179444129b9SMatthew G. Knepley /*
180444129b9SMatthew G. Knepley   CASE: quadratic
181444129b9SMatthew G. Knepley   In 2D we use exact solution:
182444129b9SMatthew G. Knepley 
183444129b9SMatthew G. Knepley     u = t + x^2 + y^2
184444129b9SMatthew G. Knepley     v = t + 2x^2 - 2xy
185444129b9SMatthew G. Knepley     p = x + y - 1
186444129b9SMatthew G. Knepley     T = t + x + y + 1
187444129b9SMatthew G. Knepley   rho = p^{th} / T
188444129b9SMatthew G. Knepley 
189444129b9SMatthew G. Knepley   so that
190444129b9SMatthew G. Knepley 
191444129b9SMatthew G. Knepley     \nabla \cdot u = 2x - 2x = 0
192444129b9SMatthew G. Knepley     grad u = <<2 x, 4x - 2y>, <2 y, -2x>>
193444129b9SMatthew G. Knepley     epsilon(u) = 1/2 (grad u + grad u^T) = <<2x, 2x>, <2x, -2x>>
194444129b9SMatthew G. Knepley     epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u)
195444129b9SMatthew G. Knepley     div epsilon'(u) = <2, 2>
196444129b9SMatthew G. Knepley 
197444129b9SMatthew 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
198444129b9SMatthew 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>
199444129b9SMatthew 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>
200444129b9SMatthew G. Knepley 
201444129b9SMatthew G. Knepley   g = S rho_t + div (rho u)
202444129b9SMatthew G. Knepley     = -S pth T_t/T^2 + rho div (u) + u . grad rho
203444129b9SMatthew G. Knepley     = -S pth 1/T^2 - pth u . grad T / T^2
204444129b9SMatthew G. Knepley     = -pth / T^2 (S + 2t + 3 x^2 - 2xy + y^2)
205444129b9SMatthew G. Knepley 
206444129b9SMatthew G. Knepley   Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T
207444129b9SMatthew G. Knepley     = c_p S pth / T + c_p pth (2t + 3 x^2 - 2xy + y^2) / T - k/Pe 0
208444129b9SMatthew G. Knepley     = c_p pth / T (S + 2t + 3 x^2 - 2xy + y^2)
209444129b9SMatthew G. Knepley */
210*9371c9d4SSatish Balay static void f0_conduct_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
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 
224*9371c9d4SSatish Balay static void f0_conduct_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
225444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
226444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
227444129b9SMatthew G. Knepley 
228444129b9SMatthew 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.);
229444129b9SMatthew G. Knepley }
230444129b9SMatthew G. Knepley 
231*9371c9d4SSatish Balay static void f0_conduct_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
232444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
233444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
234444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
235444129b9SMatthew G. Knepley 
236444129b9SMatthew 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.);
237649ef022SMatthew Knepley }
238649ef022SMatthew Knepley 
239649ef022SMatthew Knepley /*
240649ef022SMatthew Knepley   CASE: cubic
241649ef022SMatthew Knepley   In 2D we use exact solution:
242649ef022SMatthew Knepley 
243649ef022SMatthew Knepley     u = t + x^3 + y^3
244649ef022SMatthew Knepley     v = t + 2x^3 - 3x^2y
245649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
246649ef022SMatthew Knepley     T = t + 1/2 x^2 + 1/2 y^2
247649ef022SMatthew Knepley     f = < t(3x^2 + 3y^2) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y) + 3x + 1,
248649ef022SMatthew Knepley           t(3x^2 - 6xy) + 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y + 1>
249649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 + xt + yt - 2\alpha + 1
250649ef022SMatthew Knepley 
251649ef022SMatthew Knepley   so that
252649ef022SMatthew Knepley 
253649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
254649ef022SMatthew Knepley 
255649ef022SMatthew Knepley   du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p - f
256649ef022SMatthew 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
257649ef022SMatthew Knepley 
258649ef022SMatthew 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
259649ef022SMatthew Knepley */
260*9371c9d4SSatish Balay static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
261649ef022SMatthew Knepley   u[0] = time + X[0] * X[0] * X[0] + X[1] * X[1] * X[1];
262649ef022SMatthew Knepley   u[1] = time + 2.0 * X[0] * X[0] * X[0] - 3.0 * X[0] * X[0] * X[1];
263649ef022SMatthew Knepley   return 0;
264649ef022SMatthew Knepley }
265*9371c9d4SSatish Balay static PetscErrorCode cubic_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
266649ef022SMatthew Knepley   u[0] = 1.0;
267649ef022SMatthew Knepley   u[1] = 1.0;
268649ef022SMatthew Knepley   return 0;
269649ef022SMatthew Knepley }
270649ef022SMatthew Knepley 
271*9371c9d4SSatish Balay static PetscErrorCode cubic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
272649ef022SMatthew Knepley   p[0] = 3.0 * X[0] * X[0] / 2.0 + 3.0 * X[1] * X[1] / 2.0 - 1.0;
273649ef022SMatthew Knepley   return 0;
274649ef022SMatthew Knepley }
275649ef022SMatthew Knepley 
276*9371c9d4SSatish Balay static PetscErrorCode cubic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
277649ef022SMatthew Knepley   T[0] = time + X[0] * X[0] / 2.0 + X[1] * X[1] / 2.0;
278649ef022SMatthew Knepley   return 0;
279649ef022SMatthew Knepley }
280*9371c9d4SSatish Balay static PetscErrorCode cubic_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
281649ef022SMatthew Knepley   T[0] = 1.0;
282649ef022SMatthew Knepley   return 0;
283649ef022SMatthew Knepley }
284649ef022SMatthew Knepley 
285*9371c9d4SSatish Balay static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
286444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
287649ef022SMatthew Knepley 
288649ef022SMatthew 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);
289649ef022SMatthew 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);
290649ef022SMatthew Knepley }
291649ef022SMatthew Knepley 
292*9371c9d4SSatish Balay static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
293444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
294649ef022SMatthew Knepley 
295444129b9SMatthew 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;
296649ef022SMatthew Knepley }
297649ef022SMatthew Knepley 
298649ef022SMatthew Knepley /*
299649ef022SMatthew Knepley   CASE: cubic-trigonometric
300649ef022SMatthew Knepley   In 2D we use exact solution:
301649ef022SMatthew Knepley 
302649ef022SMatthew Knepley     u = beta cos t + x^3 + y^3
303649ef022SMatthew Knepley     v = beta sin t + 2x^3 - 3x^2y
304649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
305649ef022SMatthew Knepley     T = 20 cos t + 1/2 x^2 + 1/2 y^2
306649ef022SMatthew 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,
307649ef022SMatthew Knepley           beta cos t (6x^2 - 6xy) - beta sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu(12x - 6y) + 3y>
308649ef022SMatthew Knepley     Q = beta cos t x + beta sin t (y - 1) + x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2\alpha
309649ef022SMatthew Knepley 
310649ef022SMatthew Knepley   so that
311649ef022SMatthew Knepley 
312649ef022SMatthew Knepley     \nabla \cdot u = 3x^2 - 3x^2 = 0
313649ef022SMatthew Knepley 
314649ef022SMatthew Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
315649ef022SMatthew 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>
316649ef022SMatthew 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>
317649ef022SMatthew Knepley     = <cos t (3x^2)       + sin t (3y^2 - 1) + 3x^5 + 6x^3y^2 - 6x^2y^3 - \nu (6x + 6y)  + 3x,
318649ef022SMatthew Knepley        cos t (6x^2 - 6xy) - sin t (3x^2)     + 3x^4y + 6x^2y^3 - 6xy^4  - \nu (12x - 6y) + 3y>
319649ef022SMatthew Knepley 
320649ef022SMatthew Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
321649ef022SMatthew Knepley     = -sin t + <cos t + x^3 + y^3, sin t + 2x^3 - 3x^2y> . <x, y> - 2 \alpha
322649ef022SMatthew Knepley     = -sin t + cos t (x) + x^4 + xy^3 + sin t (y) + 2x^3y - 3x^2y^2 - 2 \alpha
323649ef022SMatthew Knepley     = cos t x + sin t (y - 1) + (x^4 + 2x^3y - 3x^2y^2 + xy^3 - 2 \alpha)
324649ef022SMatthew Knepley */
325*9371c9d4SSatish Balay static PetscErrorCode cubic_trig_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
326649ef022SMatthew Knepley   u[0] = 100. * PetscCosReal(time) + X[0] * X[0] * X[0] + X[1] * X[1] * X[1];
327649ef022SMatthew Knepley   u[1] = 100. * PetscSinReal(time) + 2.0 * X[0] * X[0] * X[0] - 3.0 * X[0] * X[0] * X[1];
328649ef022SMatthew Knepley   return 0;
329649ef022SMatthew Knepley }
330*9371c9d4SSatish Balay static PetscErrorCode cubic_trig_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
331649ef022SMatthew Knepley   u[0] = -100. * PetscSinReal(time);
332649ef022SMatthew Knepley   u[1] = 100. * PetscCosReal(time);
333649ef022SMatthew Knepley   return 0;
334649ef022SMatthew Knepley }
335649ef022SMatthew Knepley 
336*9371c9d4SSatish Balay static PetscErrorCode cubic_trig_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
337649ef022SMatthew Knepley   p[0] = 3.0 * X[0] * X[0] / 2.0 + 3.0 * X[1] * X[1] / 2.0 - 1.0;
338649ef022SMatthew Knepley   return 0;
339649ef022SMatthew Knepley }
340649ef022SMatthew Knepley 
341*9371c9d4SSatish Balay static PetscErrorCode cubic_trig_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
342649ef022SMatthew Knepley   T[0] = 100. * PetscCosReal(time) + X[0] * X[0] / 2.0 + X[1] * X[1] / 2.0;
343649ef022SMatthew Knepley   return 0;
344649ef022SMatthew Knepley }
345*9371c9d4SSatish Balay static PetscErrorCode cubic_trig_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
346649ef022SMatthew Knepley   T[0] = -100. * PetscSinReal(time);
347649ef022SMatthew Knepley   return 0;
348649ef022SMatthew Knepley }
349649ef022SMatthew Knepley 
350*9371c9d4SSatish Balay static void f0_cubic_trig_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
351444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
352649ef022SMatthew Knepley 
353649ef022SMatthew 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];
354649ef022SMatthew 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];
355649ef022SMatthew Knepley }
356649ef022SMatthew Knepley 
357*9371c9d4SSatish Balay static void f0_cubic_trig_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
358444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
359649ef022SMatthew Knepley 
360444129b9SMatthew 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;
361649ef022SMatthew Knepley }
362649ef022SMatthew Knepley 
363606d57d4SMatthew G. Knepley /*
364444129b9SMatthew G. Knepley   CASE: Taylor-Green vortex
365606d57d4SMatthew G. Knepley   In 2D we use exact solution:
366606d57d4SMatthew G. Knepley 
367606d57d4SMatthew G. Knepley     u = 1 - cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)
368606d57d4SMatthew G. Knepley     v = 1 + sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)
369606d57d4SMatthew G. Knepley     p = -1/4 [cos(2 \pi(x - t)) + cos(2 \pi(y - t))] exp(-4 \pi^2 \nu t)
370606d57d4SMatthew G. Knepley     T = t + x + y
371606d57d4SMatthew 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))  >
372606d57d4SMatthew G. Knepley     Q = 3 + sin(\pi(x-y)) exp(-2\nu \pi^2 t)
373606d57d4SMatthew G. Knepley 
374606d57d4SMatthew G. Knepley   so that
375606d57d4SMatthew G. Knepley 
376606d57d4SMatthew 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
377606d57d4SMatthew G. Knepley 
378606d57d4SMatthew G. Knepley   f = du/dt + u \cdot \nabla u - \nu \Delta u + \nabla p
379606d57d4SMatthew 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),
380606d57d4SMatthew 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)>
381606d57d4SMatthew 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),
382606d57d4SMatthew 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)>
383606d57d4SMatthew 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),
384606d57d4SMatthew 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)>
385606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
386606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
387606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
388606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
389606d57d4SMatthew 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),
390606d57d4SMatthew 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)>
391606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
392606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
393606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
394606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
395606d57d4SMatthew G. Knepley     + <-\pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
396606d57d4SMatthew G. Knepley        -\pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
397606d57d4SMatthew G. Knepley     + <-2\pi^2 cos(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
398606d57d4SMatthew G. Knepley         2\pi^2 sin(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
399606d57d4SMatthew G. Knepley     + < \pi/2 sin(2\pi(x - t)) exp(-4 \pi^2 \nu t),
400606d57d4SMatthew G. Knepley         \pi/2 sin(2\pi(y - t)) exp(-4 \pi^2 \nu t)>
401606d57d4SMatthew 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),
402606d57d4SMatthew 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)>
403606d57d4SMatthew G. Knepley     + < \pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t),
404606d57d4SMatthew G. Knepley         \pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t)>
405606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)) exp(-2 \pi^2 \nu t),
406606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t)) exp(-2 \pi^2 \nu t)>
407606d57d4SMatthew G. Knepley     = < \pi cos(\pi(x - t)) cos(\pi(y - t)),
408606d57d4SMatthew G. Knepley         \pi sin(\pi(x - t)) sin(\pi(y - t))>
409606d57d4SMatthew G. Knepley     + <-\pi cos(\pi(x - t)) cos(\pi(y - t)),
410606d57d4SMatthew G. Knepley        -\pi sin(\pi(x - t)) sin(\pi(y - t))> = 0
411606d57d4SMatthew G. Knepley   Q = dT/dt + u \cdot \nabla T - \alpha \Delta T
412606d57d4SMatthew G. Knepley     = 1 + u \cdot <1, 1> - 0
413606d57d4SMatthew G. Knepley     = 1 + u + v
414606d57d4SMatthew G. Knepley */
415606d57d4SMatthew G. Knepley 
416*9371c9d4SSatish Balay static PetscErrorCode taylor_green_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
417606d57d4SMatthew G. Knepley   u[0] = 1 - PetscCosReal(PETSC_PI * (X[0] - time)) * PetscSinReal(PETSC_PI * (X[1] - time)) * PetscExpReal(-2 * PETSC_PI * PETSC_PI * time);
418606d57d4SMatthew G. Knepley   u[1] = 1 + PetscSinReal(PETSC_PI * (X[0] - time)) * PetscCosReal(PETSC_PI * (X[1] - time)) * PetscExpReal(-2 * PETSC_PI * PETSC_PI * time);
419606d57d4SMatthew G. Knepley   return 0;
420606d57d4SMatthew G. Knepley }
421*9371c9d4SSatish Balay static PetscErrorCode taylor_green_u_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
422*9371c9d4SSatish Balay   u[0] = -PETSC_PI * (PetscSinReal(PETSC_PI * (X[0] - time)) * PetscSinReal(PETSC_PI * (X[1] - time)) - PetscCosReal(PETSC_PI * (X[0] - time)) * PetscCosReal(PETSC_PI * (X[1] - time)) - 2 * PETSC_PI * PetscCosReal(PETSC_PI * (X[0] - time)) * PetscSinReal(PETSC_PI * (X[1] - time))) * PetscExpReal(-2 * PETSC_PI * PETSC_PI * time);
423*9371c9d4SSatish Balay   u[1] = PETSC_PI * (PetscSinReal(PETSC_PI * (X[0] - time)) * PetscSinReal(PETSC_PI * (X[1] - time)) - PetscCosReal(PETSC_PI * (X[0] - time)) * PetscCosReal(PETSC_PI * (X[1] - time)) - 2 * PETSC_PI * PetscSinReal(PETSC_PI * (X[0] - time)) * PetscCosReal(PETSC_PI * (X[1] - time))) * PetscExpReal(-2 * PETSC_PI * PETSC_PI * time);
424606d57d4SMatthew G. Knepley   return 0;
425606d57d4SMatthew G. Knepley }
426606d57d4SMatthew G. Knepley 
427*9371c9d4SSatish Balay static PetscErrorCode taylor_green_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
428606d57d4SMatthew 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);
429606d57d4SMatthew G. Knepley   return 0;
430606d57d4SMatthew G. Knepley }
431606d57d4SMatthew G. Knepley 
432*9371c9d4SSatish Balay static PetscErrorCode taylor_green_p_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
433*9371c9d4SSatish Balay   p[0] = PETSC_PI * (0.5 * (PetscSinReal(2 * PETSC_PI * (X[0] - time)) + PetscSinReal(2 * PETSC_PI * (X[1] - time))) + PETSC_PI * (PetscCosReal(2 * PETSC_PI * (X[0] - time)) + PetscCosReal(2 * PETSC_PI * (X[1] - time)))) * PetscExpReal(-4 * PETSC_PI * PETSC_PI * time);
434606d57d4SMatthew G. Knepley   return 0;
435606d57d4SMatthew G. Knepley }
436606d57d4SMatthew G. Knepley 
437*9371c9d4SSatish Balay static PetscErrorCode taylor_green_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
438606d57d4SMatthew G. Knepley   T[0] = time + X[0] + X[1];
439606d57d4SMatthew G. Knepley   return 0;
440606d57d4SMatthew G. Knepley }
441*9371c9d4SSatish Balay static PetscErrorCode taylor_green_T_t(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
442606d57d4SMatthew G. Knepley   T[0] = 1.0;
443606d57d4SMatthew G. Knepley   return 0;
444606d57d4SMatthew G. Knepley }
445606d57d4SMatthew G. Knepley 
446*9371c9d4SSatish Balay static void f0_taylor_green_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
447606d57d4SMatthew G. Knepley   PetscScalar vel[2];
448606d57d4SMatthew G. Knepley 
449606d57d4SMatthew G. Knepley   taylor_green_u(dim, t, X, Nf, vel, NULL);
450444129b9SMatthew G. Knepley   f0[0] -= 1.0 + vel[0] + vel[1];
451606d57d4SMatthew G. Knepley }
452606d57d4SMatthew G. Knepley 
453444129b9SMatthew G. Knepley /*
454444129b9SMatthew G. Knepley   CASE: Pipe flow
455444129b9SMatthew G. Knepley   Poiseuille flow, with the incoming fluid having a parabolic temperature profile and the side walls being held at T_in
456444129b9SMatthew G. Knepley 
457444129b9SMatthew G. Knepley     u = \Delta Re/(2 mu) y (1 - y)
458444129b9SMatthew G. Knepley     v = 0
459444129b9SMatthew G. Knepley     p = -\Delta x
460444129b9SMatthew G. Knepley     T = y (1 - y) + T_in
461444129b9SMatthew G. Knepley   rho = p^{th} / T
462444129b9SMatthew G. Knepley 
463444129b9SMatthew G. Knepley   so that
464444129b9SMatthew G. Knepley 
465444129b9SMatthew G. Knepley     \nabla \cdot u = 0 - 0 = 0
466444129b9SMatthew G. Knepley     grad u = \Delta Re/(2 mu) <<0, 0>, <1 - 2y, 0>>
467444129b9SMatthew G. Knepley     epsilon(u) = 1/2 (grad u + grad u^T) = \Delta Re/(4 mu) <<0, 1 - 2y>, <<1 - 2y, 0>>
468444129b9SMatthew G. Knepley     epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u)
469444129b9SMatthew G. Knepley     div epsilon'(u) = -\Delta Re/(2 mu) <1, 0>
470444129b9SMatthew G. Knepley 
471444129b9SMatthew 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
472444129b9SMatthew G. Knepley     = 0 + 0 - div (2\mu/Re \epsilon'(u) - pI) + rho / F^2 \hat y
473444129b9SMatthew G. Knepley     = -\Delta div <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> + rho / F^2 \hat y
474444129b9SMatthew G. Knepley     = \Delta <1, 0> - \Delta <1, 0> + rho/F^2 <0, 1>
475444129b9SMatthew G. Knepley     = rho/F^2 <0, 1>
476444129b9SMatthew G. Knepley 
477444129b9SMatthew G. Knepley   g = S rho_t + div (rho u)
478444129b9SMatthew G. Knepley     = 0 + rho div (u) + u . grad rho
479444129b9SMatthew G. Knepley     = 0 + 0 - pth u . grad T / T^2
480444129b9SMatthew G. Knepley     = 0
481444129b9SMatthew G. Knepley 
482444129b9SMatthew G. Knepley   Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T
483444129b9SMatthew G. Knepley     = 0 + c_p pth / T 0 + 2 k/Pe
484444129b9SMatthew G. Knepley     = 2 k/Pe
485444129b9SMatthew G. Knepley 
486444129b9SMatthew G. Knepley   The boundary conditions on the top and bottom are zero velocity and T_in temperature. The boundary term is
487444129b9SMatthew G. Knepley 
488444129b9SMatthew G. Knepley     (2\mu/Re \epsilon'(u) - p I) . n = \Delta <<x, (1 - 2y)/2>, <<(1 - 2y)/2, x>> . n
489444129b9SMatthew G. Knepley 
490444129b9SMatthew G. Knepley   so that
491444129b9SMatthew G. Knepley 
492444129b9SMatthew G. Knepley     x = 0: \Delta <<0, (1 - 2y)/2>, <<(1 - 2y)/2, 0>> . <-1, 0> = <0, (2y - 1)/2>
493444129b9SMatthew G. Knepley     x = 1: \Delta <<1, (1 - 2y)/2>, <<(1 - 2y)/2, 1>> . <1, 0> = <1, (1 - 2y)/2>
494444129b9SMatthew G. Knepley */
495*9371c9d4SSatish Balay static PetscErrorCode pipe_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
496444129b9SMatthew G. Knepley   Parameter *param = (Parameter *)ctx;
497444129b9SMatthew G. Knepley 
498444129b9SMatthew G. Knepley   u[0] = (0.5 * param->Reynolds / param->mu) * X[1] * (1.0 - X[1]);
499444129b9SMatthew G. Knepley   u[1] = 0.0;
500444129b9SMatthew G. Knepley   return 0;
501444129b9SMatthew G. Knepley }
502*9371c9d4SSatish Balay static PetscErrorCode pipe_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
503444129b9SMatthew G. Knepley   u[0] = 0.0;
504444129b9SMatthew G. Knepley   u[1] = 0.0;
505444129b9SMatthew G. Knepley   return 0;
506444129b9SMatthew G. Knepley }
507444129b9SMatthew G. Knepley 
508*9371c9d4SSatish Balay static PetscErrorCode pipe_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
509444129b9SMatthew G. Knepley   p[0] = -X[0];
510444129b9SMatthew G. Knepley   return 0;
511444129b9SMatthew G. Knepley }
512*9371c9d4SSatish Balay static PetscErrorCode pipe_p_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
513444129b9SMatthew G. Knepley   p[0] = 0.0;
514444129b9SMatthew G. Knepley   return 0;
515444129b9SMatthew G. Knepley }
516444129b9SMatthew G. Knepley 
517*9371c9d4SSatish Balay static PetscErrorCode pipe_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
518444129b9SMatthew G. Knepley   Parameter *param = (Parameter *)ctx;
519444129b9SMatthew G. Knepley 
520444129b9SMatthew G. Knepley   T[0] = X[1] * (1.0 - X[1]) + param->T_in;
521444129b9SMatthew G. Knepley   return 0;
522444129b9SMatthew G. Knepley }
523*9371c9d4SSatish Balay static PetscErrorCode pipe_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
524444129b9SMatthew G. Knepley   T[0] = 0.0;
525444129b9SMatthew G. Knepley   return 0;
526444129b9SMatthew G. Knepley }
527444129b9SMatthew G. Knepley 
528*9371c9d4SSatish Balay static void f0_conduct_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
529444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
530444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
531444129b9SMatthew G. Knepley   const PetscReal T_in = PetscRealPart(constants[T_IN]);
532444129b9SMatthew G. Knepley   const PetscReal rho  = p_th / (X[1] * (1. - X[1]) + T_in);
533444129b9SMatthew G. Knepley   const PetscInt  gd   = (PetscInt)PetscRealPart(constants[G_DIR]);
534444129b9SMatthew G. Knepley 
535444129b9SMatthew G. Knepley   f0[gd] -= rho / PetscSqr(F);
536444129b9SMatthew G. Knepley }
537444129b9SMatthew G. Knepley 
538*9371c9d4SSatish Balay static void f0_conduct_bd_pipe_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
539444129b9SMatthew G. Knepley   PetscReal sigma[4] = {X[0], 0.5 * (1. - 2. * X[1]), 0.5 * (1. - 2. * X[1]), X[0]};
540444129b9SMatthew G. Knepley   PetscInt  d, e;
541444129b9SMatthew G. Knepley 
542444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) {
543*9371c9d4SSatish Balay     for (e = 0; e < dim; ++e) { f0[d] -= sigma[d * dim + e] * n[e]; }
544444129b9SMatthew G. Knepley   }
545444129b9SMatthew G. Knepley }
546444129b9SMatthew G. Knepley 
547*9371c9d4SSatish Balay static void f0_conduct_pipe_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
548444129b9SMatthew G. Knepley   f0[0] += 0.0;
549444129b9SMatthew G. Knepley }
550444129b9SMatthew G. Knepley 
551*9371c9d4SSatish Balay static void f0_conduct_pipe_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
552444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
553444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
554444129b9SMatthew G. Knepley 
555444129b9SMatthew G. Knepley   f0[0] -= 2 * k / Pe;
556444129b9SMatthew G. Knepley }
557444129b9SMatthew G. Knepley 
558367970cfSMatthew G. Knepley /*
559367970cfSMatthew G. Knepley   CASE: Wiggly pipe flow
560367970cfSMatthew G. Knepley   Perturbed Poiseuille flow, with the incoming fluid having a perturbed parabolic temperature profile and the side walls being held at T_in
561367970cfSMatthew G. Knepley 
562367970cfSMatthew G. Knepley     u = \Delta Re/(2 mu) [y (1 - y) + a sin(pi y)]
563367970cfSMatthew G. Knepley     v = 0
564367970cfSMatthew G. Knepley     p = -\Delta x
565367970cfSMatthew G. Knepley     T = y (1 - y) + a sin(pi y) + T_in
566367970cfSMatthew G. Knepley   rho = p^{th} / T
567367970cfSMatthew G. Knepley 
568367970cfSMatthew G. Knepley   so that
569367970cfSMatthew G. Knepley 
570367970cfSMatthew G. Knepley     \nabla \cdot u = 0 - 0 = 0
571367970cfSMatthew G. Knepley     grad u = \Delta Re/(2 mu) <<0, 0>, <1 - 2y + a pi cos(pi y), 0>>
572367970cfSMatthew G. Knepley     epsilon(u) = 1/2 (grad u + grad u^T) = \Delta Re/(4 mu) <<0, 1 - 2y + a pi cos(pi y)>, <<1 - 2y + a pi cos(pi y), 0>>
573367970cfSMatthew G. Knepley     epsilon'(u) = epsilon(u) - 1/3 (div u) I = epsilon(u)
574367970cfSMatthew G. Knepley     div epsilon'(u) = -\Delta Re/(2 mu) <1 + a pi^2/2 sin(pi y), 0>
575367970cfSMatthew G. Knepley 
576367970cfSMatthew 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
577367970cfSMatthew G. Knepley     = 0 + 0 - div (2\mu/Re \epsilon'(u) - pI) + rho / F^2 \hat y
578367970cfSMatthew G. Knepley     = -\Delta div <<x, (1 - 2y)/2 + a pi/2 cos(pi y)>, <<(1 - 2y)/2 + a pi/2 cos(pi y), x>> + rho / F^2 \hat y
579367970cfSMatthew G. Knepley     = -\Delta <1 - 1 - a pi^2/2 sin(pi y), 0> + rho/F^2 <0, 1>
580367970cfSMatthew G. Knepley     = a \Delta pi^2/2 sin(pi y) <1, 0> + rho/F^2 <0, 1>
581367970cfSMatthew G. Knepley 
582367970cfSMatthew G. Knepley   g = S rho_t + div (rho u)
583367970cfSMatthew G. Knepley     = 0 + rho div (u) + u . grad rho
584367970cfSMatthew G. Knepley     = 0 + 0 - pth u . grad T / T^2
585367970cfSMatthew G. Knepley     = 0
586367970cfSMatthew G. Knepley 
587367970cfSMatthew G. Knepley   Q = rho c_p S dT/dt + rho c_p u . grad T - 1/Pe div k grad T
588367970cfSMatthew G. Knepley     = 0 + c_p pth / T 0 - k/Pe div <0, 1 - 2y + a pi cos(pi y)>
589367970cfSMatthew G. Knepley     = - k/Pe (-2 - a pi^2 sin(pi y))
590367970cfSMatthew G. Knepley     = 2 k/Pe (1 + a pi^2/2 sin(pi y))
591367970cfSMatthew G. Knepley 
592367970cfSMatthew G. Knepley   The boundary conditions on the top and bottom are zero velocity and T_in temperature. The boundary term is
593367970cfSMatthew G. Knepley 
594367970cfSMatthew G. Knepley     (2\mu/Re \epsilon'(u) - p I) . n = \Delta <<x, (1 - 2y)/2 + a pi/2 cos(pi y)>, <<(1 - 2y)/2 + a pi/2 cos(pi y), x>> . n
595367970cfSMatthew G. Knepley 
596367970cfSMatthew G. Knepley   so that
597367970cfSMatthew G. Knepley 
598367970cfSMatthew G. Knepley     x = 0: \Delta <<0, (1 - 2y)/2>, <<(1 - 2y)/2, 0>> . <-1, 0> = <0, (2y - 1)/2 - a pi/2 cos(pi y)>
599367970cfSMatthew G. Knepley     x = 1: \Delta <<1, (1 - 2y)/2>, <<(1 - 2y)/2, 1>> . < 1, 0> = <1, (1 - 2y)/2 + a pi/2 cos(pi y)>
600367970cfSMatthew G. Knepley */
601*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
602367970cfSMatthew G. Knepley   Parameter *param = (Parameter *)ctx;
603367970cfSMatthew G. Knepley 
604367970cfSMatthew G. Knepley   u[0] = (0.5 * param->Reynolds / param->mu) * (X[1] * (1.0 - X[1]) + param->epsilon * PetscSinReal(PETSC_PI * X[1]));
605367970cfSMatthew G. Knepley   u[1] = 0.0;
606367970cfSMatthew G. Knepley   return 0;
607367970cfSMatthew G. Knepley }
608*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_u_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
609367970cfSMatthew G. Knepley   u[0] = 0.0;
610367970cfSMatthew G. Knepley   u[1] = 0.0;
611367970cfSMatthew G. Knepley   return 0;
612367970cfSMatthew G. Knepley }
613367970cfSMatthew G. Knepley 
614*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
615367970cfSMatthew G. Knepley   p[0] = -X[0];
616367970cfSMatthew G. Knepley   return 0;
617367970cfSMatthew G. Knepley }
618*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_p_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
619367970cfSMatthew G. Knepley   p[0] = 0.0;
620367970cfSMatthew G. Knepley   return 0;
621367970cfSMatthew G. Knepley }
622367970cfSMatthew G. Knepley 
623*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_T(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
624367970cfSMatthew G. Knepley   Parameter *param = (Parameter *)ctx;
625367970cfSMatthew G. Knepley 
626367970cfSMatthew G. Knepley   T[0] = X[1] * (1.0 - X[1]) + param->epsilon * PetscSinReal(PETSC_PI * X[1]) + param->T_in;
627367970cfSMatthew G. Knepley   return 0;
628367970cfSMatthew G. Knepley }
629*9371c9d4SSatish Balay static PetscErrorCode pipe_wiggly_T_t(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
630367970cfSMatthew G. Knepley   T[0] = 0.0;
631367970cfSMatthew G. Knepley   return 0;
632367970cfSMatthew G. Knepley }
633367970cfSMatthew G. Knepley 
634*9371c9d4SSatish Balay static void f0_conduct_pipe_wiggly_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
635367970cfSMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
636367970cfSMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
637367970cfSMatthew G. Knepley   const PetscReal T_in = PetscRealPart(constants[T_IN]);
638367970cfSMatthew G. Knepley   const PetscReal eps  = PetscRealPart(constants[EPSILON]);
639367970cfSMatthew G. Knepley   const PetscReal rho  = p_th / (X[1] * (1. - X[1]) + T_in);
640367970cfSMatthew G. Knepley   const PetscInt  gd   = (PetscInt)PetscRealPart(constants[G_DIR]);
641367970cfSMatthew G. Knepley 
642367970cfSMatthew G. Knepley   f0[0] -= eps * 0.5 * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * X[1]);
643367970cfSMatthew G. Knepley   f0[gd] -= rho / PetscSqr(F);
644367970cfSMatthew G. Knepley }
645367970cfSMatthew G. Knepley 
646*9371c9d4SSatish Balay static void f0_conduct_bd_pipe_wiggly_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
647367970cfSMatthew G. Knepley   const PetscReal eps      = PetscRealPart(constants[EPSILON]);
648367970cfSMatthew G. Knepley   PetscReal       sigma[4] = {X[0], 0.5 * (1. - 2. * X[1]) + eps * 0.5 * PETSC_PI * PetscCosReal(PETSC_PI * X[1]), 0.5 * (1. - 2. * X[1]) + eps * 0.5 * PETSC_PI * PetscCosReal(PETSC_PI * X[1]), X[0]};
649367970cfSMatthew G. Knepley   PetscInt        d, e;
650367970cfSMatthew G. Knepley 
651367970cfSMatthew G. Knepley   for (d = 0; d < dim; ++d) {
652*9371c9d4SSatish Balay     for (e = 0; e < dim; ++e) { f0[d] -= sigma[d * dim + e] * n[e]; }
653367970cfSMatthew G. Knepley   }
654367970cfSMatthew G. Knepley }
655367970cfSMatthew G. Knepley 
656*9371c9d4SSatish Balay static void f0_conduct_pipe_wiggly_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
657367970cfSMatthew G. Knepley   f0[0] += 0.0;
658367970cfSMatthew G. Knepley }
659367970cfSMatthew G. Knepley 
660*9371c9d4SSatish Balay static void f0_conduct_pipe_wiggly_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
661367970cfSMatthew G. Knepley   const PetscReal k   = PetscRealPart(constants[K]);
662367970cfSMatthew G. Knepley   const PetscReal Pe  = PetscRealPart(constants[PECLET]);
663367970cfSMatthew G. Knepley   const PetscReal eps = PetscRealPart(constants[EPSILON]);
664367970cfSMatthew G. Knepley 
665367970cfSMatthew G. Knepley   f0[0] -= 2 * k / Pe * (1.0 + eps * 0.5 * PetscSqr(PETSC_PI) * PetscSinReal(PETSC_PI * X[1]));
666367970cfSMatthew G. Knepley }
667367970cfSMatthew G. Knepley 
668444129b9SMatthew G. Knepley /*      Physics Kernels      */
669444129b9SMatthew G. Knepley 
670*9371c9d4SSatish Balay static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
671649ef022SMatthew Knepley   PetscInt d;
672649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
673649ef022SMatthew Knepley }
674649ef022SMatthew Knepley 
675444129b9SMatthew 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 */
676*9371c9d4SSatish Balay static void f0_conduct_q(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
677444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
678444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
679444129b9SMatthew G. Knepley   PetscInt        d;
680444129b9SMatthew G. Knepley 
681444129b9SMatthew G. Knepley   // -\frac{S p^{th}}{T^2} \frac{\partial T}{\partial t}
682444129b9SMatthew G. Knepley   f0[0] += -u_t[uOff[TEMP]] * S * p_th / PetscSqr(u[uOff[TEMP]]);
683444129b9SMatthew G. Knepley 
684444129b9SMatthew G. Knepley   // \frac{p^{th}}{T} \nabla \cdot \vb{u}
685*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { f0[0] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + d * dim + d]; }
686444129b9SMatthew G. Knepley 
687444129b9SMatthew G. Knepley   // - \frac{p^{th}}{T^2} \vb{u} \cdot \nabla T
688*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { f0[0] -= p_th / (u[uOff[TEMP]] * u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; }
689444129b9SMatthew G. Knepley 
690444129b9SMatthew G. Knepley   // Add in any fixed source term
691*9371c9d4SSatish Balay   if (NfAux > 0) { f0[0] += a[aOff[MASS]]; }
692444129b9SMatthew G. Knepley }
693444129b9SMatthew G. Knepley 
694444129b9SMatthew G. Knepley /* \vb{u}_t + \vb{u} \cdot \nabla\vb{u} */
695*9371c9d4SSatish Balay static void f0_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
696444129b9SMatthew G. Knepley   const PetscInt Nc = dim;
697444129b9SMatthew G. Knepley   PetscInt       c, d;
698444129b9SMatthew G. Knepley 
699444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
700444129b9SMatthew G. Knepley     /* \vb{u}_t */
701444129b9SMatthew G. Knepley     f0[c] += u_t[uOff[VEL] + c];
702444129b9SMatthew G. Knepley     /* \vb{u} \cdot \nabla\vb{u} */
703444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) f0[c] += u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c * dim + d];
704444129b9SMatthew G. Knepley   }
705444129b9SMatthew G. Knepley }
706444129b9SMatthew G. Knepley 
707444129b9SMatthew G. Knepley /* \rho S \frac{\partial \vb{u}}{\partial t} + \rho \vb{u} \cdot \nabla \vb{u} + \rho \frac{\hat{\vb{z}}}{F^2} */
708*9371c9d4SSatish Balay static void f0_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
709444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
710444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
711444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
712444129b9SMatthew G. Knepley   const PetscReal rho  = p_th / PetscRealPart(u[uOff[TEMP]]);
713444129b9SMatthew G. Knepley   const PetscInt  gdir = (PetscInt)PetscRealPart(constants[G_DIR]);
714444129b9SMatthew G. Knepley   PetscInt        Nc   = dim;
715444129b9SMatthew G. Knepley   PetscInt        c, d;
716444129b9SMatthew G. Knepley 
717444129b9SMatthew G. Knepley   // \rho S \frac{\partial \vb{u}}{\partial t}
718*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { f0[d] = rho * S * u_t[uOff[VEL] + d]; }
719444129b9SMatthew G. Knepley 
720444129b9SMatthew G. Knepley   // \rho \vb{u} \cdot \nabla \vb{u}
721444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
722*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { f0[c] += rho * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c * dim + d]; }
723444129b9SMatthew G. Knepley   }
724444129b9SMatthew G. Knepley 
725444129b9SMatthew G. Knepley   // rho \hat{z}/F^2
726444129b9SMatthew G. Knepley   f0[gdir] += rho / (F * F);
727444129b9SMatthew G. Knepley 
728444129b9SMatthew G. Knepley   // Add in any fixed source term
729444129b9SMatthew G. Knepley   if (NfAux > 0) {
730*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { f0[d] += a[aOff[MOMENTUM] + d]; }
731444129b9SMatthew G. Knepley   }
732444129b9SMatthew G. Knepley }
733444129b9SMatthew G. Knepley 
734649ef022SMatthew Knepley /*f1_v = \nu[grad(u) + grad(u)^T] - pI */
735*9371c9d4SSatish Balay static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
736444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
737649ef022SMatthew Knepley   const PetscInt  Nc = dim;
738649ef022SMatthew Knepley   PetscInt        c, d;
739649ef022SMatthew Knepley 
740649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
741*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { f1[c * dim + d] = nu * (u_x[c * dim + d] + u_x[d * dim + c]); }
742649ef022SMatthew Knepley     f1[c * dim + c] -= u[uOff[1]];
743649ef022SMatthew Knepley   }
744649ef022SMatthew Knepley }
745649ef022SMatthew Knepley 
746444129b9SMatthew G. Knepley /* 2 \mu/Re (1/2 (\nabla \vb{u} + \nabla \vb{u}^T) - 1/3 (\nabla \cdot \vb{u}) I) - p I */
747*9371c9d4SSatish Balay static void f1_conduct_v(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
748444129b9SMatthew G. Knepley   const PetscReal Re    = PetscRealPart(constants[REYNOLDS]);
749444129b9SMatthew G. Knepley   const PetscReal mu    = PetscRealPart(constants[MU]);
750444129b9SMatthew G. Knepley   const PetscReal coef  = mu / Re;
751444129b9SMatthew G. Knepley   PetscReal       u_div = 0.0;
752444129b9SMatthew G. Knepley   const PetscInt  Nc    = dim;
753444129b9SMatthew G. Knepley   PetscInt        c, d;
754444129b9SMatthew G. Knepley 
755*9371c9d4SSatish Balay   for (c = 0; c < Nc; ++c) { u_div += PetscRealPart(u_x[uOff_x[VEL] + c * dim + c]); }
756444129b9SMatthew G. Knepley 
757444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
758444129b9SMatthew G. Knepley     // 2 \mu/Re 1/2 (\nabla \vb{u} + \nabla \vb{u}^T
759*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { f1[c * dim + d] += coef * (u_x[uOff_x[VEL] + c * dim + d] + u_x[uOff_x[VEL] + d * dim + c]); }
760444129b9SMatthew G. Knepley     // -2/3 \mu/Re (\nabla \cdot \vb{u}) I
761444129b9SMatthew G. Knepley     f1[c * dim + c] -= 2.0 * coef / 3.0 * u_div;
762444129b9SMatthew G. Knepley   }
763444129b9SMatthew G. Knepley 
764444129b9SMatthew G. Knepley   // -p I
765*9371c9d4SSatish Balay   for (c = 0; c < Nc; ++c) { f1[c * dim + c] -= u[uOff[PRES]]; }
766444129b9SMatthew G. Knepley }
767444129b9SMatthew G. Knepley 
768444129b9SMatthew G. Knepley /* T_t + \vb{u} \cdot \nabla T */
769*9371c9d4SSatish Balay static void f0_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
770444129b9SMatthew G. Knepley   PetscInt d;
771444129b9SMatthew G. Knepley 
772444129b9SMatthew G. Knepley   /* T_t */
773444129b9SMatthew G. Knepley   f0[0] += u_t[uOff[TEMP]];
774444129b9SMatthew G. Knepley   /* \vb{u} \cdot \nabla T */
775444129b9SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d];
776444129b9SMatthew G. Knepley }
777444129b9SMatthew G. Knepley 
778444129b9SMatthew 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 */
779*9371c9d4SSatish Balay static void f0_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
780444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
781444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
782444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
783444129b9SMatthew G. Knepley   PetscInt        d;
784444129b9SMatthew G. Knepley 
785444129b9SMatthew G. Knepley   // \frac{C_p S p^{th}}{T} \frac{\partial T}{\partial t}
786444129b9SMatthew G. Knepley   f0[0] = c_p * S * p_th / u[uOff[TEMP]] * u_t[uOff[TEMP]];
787444129b9SMatthew G. Knepley 
788444129b9SMatthew G. Knepley   // \frac{C_p p^{th}}{T} \vb{u} \cdot \nabla T
789*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { f0[0] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; }
790444129b9SMatthew G. Knepley 
791444129b9SMatthew G. Knepley   // Add in any fixed source term
792*9371c9d4SSatish Balay   if (NfAux > 0) { f0[0] += a[aOff[ENERGY]]; }
793444129b9SMatthew G. Knepley }
794444129b9SMatthew G. Knepley 
795*9371c9d4SSatish Balay static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
796444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
797649ef022SMatthew Knepley   PetscInt        d;
798444129b9SMatthew G. Knepley 
799649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha * u_x[uOff_x[2] + d];
800649ef022SMatthew Knepley }
801649ef022SMatthew Knepley 
802444129b9SMatthew G. Knepley /* \frac{k}{Pe} \nabla T */
803*9371c9d4SSatish Balay static void f1_conduct_w(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
804444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
805444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
806444129b9SMatthew G. Knepley   PetscInt        d;
807444129b9SMatthew G. Knepley 
808444129b9SMatthew G. Knepley   // \frac{k}{Pe} \nabla T
809*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { f1[d] = k / Pe * u_x[uOff_x[TEMP] + d]; }
810444129b9SMatthew G. Knepley }
811444129b9SMatthew G. Knepley 
812*9371c9d4SSatish Balay static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
813649ef022SMatthew Knepley   PetscInt d;
814649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
815649ef022SMatthew Knepley }
816649ef022SMatthew Knepley 
817*9371c9d4SSatish Balay static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
818649ef022SMatthew Knepley   PetscInt       c, d;
819649ef022SMatthew Knepley   const PetscInt Nc = dim;
820649ef022SMatthew Knepley 
821649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift;
822649ef022SMatthew Knepley 
823649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
824*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { g0[c * Nc + d] += u_x[c * Nc + d]; }
825649ef022SMatthew Knepley   }
826649ef022SMatthew Knepley }
827649ef022SMatthew Knepley 
828*9371c9d4SSatish Balay static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
829649ef022SMatthew Knepley   PetscInt NcI = dim;
830649ef022SMatthew Knepley   PetscInt NcJ = dim;
831649ef022SMatthew Knepley   PetscInt c, d, e;
832649ef022SMatthew Knepley 
833649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
834649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
835649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
836*9371c9d4SSatish Balay         if (c == d) { g1[(c * NcJ + d) * dim + e] += u[e]; }
837649ef022SMatthew Knepley       }
838649ef022SMatthew Knepley     }
839649ef022SMatthew Knepley   }
840649ef022SMatthew Knepley }
841649ef022SMatthew Knepley 
842*9371c9d4SSatish Balay static void g0_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
843444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
844444129b9SMatthew G. Knepley   PetscInt        d;
845444129b9SMatthew G. Knepley 
846444129b9SMatthew G. Knepley   // - \phi_i \frac{p^{th}}{T^2} \frac{\partial T}{\partial x_c} \psi_{j, u_c}
847*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g0[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u_x[uOff_x[TEMP] + d]; }
848444129b9SMatthew G. Knepley }
849444129b9SMatthew G. Knepley 
850*9371c9d4SSatish Balay static void g1_conduct_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
851444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
852444129b9SMatthew G. Knepley   PetscInt        d;
853444129b9SMatthew G. Knepley 
854444129b9SMatthew G. Knepley   // \phi_i \frac{p^{th}}{T} \frac{\partial \psi_{u_c,j}}{\partial x_c}
855*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g1[d * dim + d] = p_th / u[uOff[TEMP]]; }
856444129b9SMatthew G. Knepley }
857444129b9SMatthew G. Knepley 
858*9371c9d4SSatish Balay static void g0_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
859444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
860444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
861444129b9SMatthew G. Knepley   PetscInt        d;
862444129b9SMatthew G. Knepley 
863444129b9SMatthew G. Knepley   // - \phi_i \frac{S p^{th}}{T^2} \psi_j
864444129b9SMatthew G. Knepley   g0[0] -= S * p_th / PetscSqr(u[uOff[TEMP]]) * u_tShift;
865444129b9SMatthew G. Knepley   // \phi_i 2 \frac{S p^{th}}{T^3} T_t \psi_j
866444129b9SMatthew G. Knepley   g0[0] += 2.0 * S * p_th / PetscPowScalarInt(u[uOff[TEMP]], 3) * u_t[uOff[TEMP]];
867444129b9SMatthew 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)
868*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { 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]); }
869444129b9SMatthew G. Knepley }
870444129b9SMatthew G. Knepley 
871*9371c9d4SSatish Balay static void g1_conduct_qT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
872444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
873444129b9SMatthew G. Knepley   PetscInt        d;
874444129b9SMatthew G. Knepley 
875444129b9SMatthew G. Knepley   // - \phi_i \frac{p^{th}}{T^2} \vb{u} \cdot \nabla \psi_j
876*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g1[d] = -p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d]; }
877444129b9SMatthew G. Knepley }
878444129b9SMatthew G. Knepley 
879*9371c9d4SSatish Balay static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) {
880649ef022SMatthew Knepley   PetscInt d;
881649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0;
882649ef022SMatthew Knepley }
883649ef022SMatthew Knepley 
884*9371c9d4SSatish Balay static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
885444129b9SMatthew G. Knepley   const PetscReal nu = PetscRealPart(constants[NU]);
886649ef022SMatthew Knepley   const PetscInt  Nc = dim;
887649ef022SMatthew Knepley   PetscInt        c, d;
888649ef022SMatthew Knepley 
889649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
890649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
891606d57d4SMatthew G. Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += nu;
892606d57d4SMatthew G. Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += nu;
893649ef022SMatthew Knepley     }
894649ef022SMatthew Knepley   }
895649ef022SMatthew Knepley }
896649ef022SMatthew Knepley 
897*9371c9d4SSatish Balay static void g0_conduct_vT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
898444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
899444129b9SMatthew G. Knepley   const PetscReal F    = PetscRealPart(constants[FROUDE]);
900444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
901444129b9SMatthew G. Knepley   const PetscInt  gdir = (PetscInt)PetscRealPart(constants[G_DIR]);
902444129b9SMatthew G. Knepley   const PetscInt  Nc   = dim;
903444129b9SMatthew G. Knepley   PetscInt        c, d;
904444129b9SMatthew G. Knepley 
905444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vb{u}_t \frac{p^{th} S}{T^2} \psi_j
906*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g0[d] -= p_th * S / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[VEL] + d]; }
907444129b9SMatthew G. Knepley 
908444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vb{u} \cdot \nabla \vb{u} \frac{p^{th}}{T^2} \psi_j
909444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
910*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { g0[c] -= p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[VEL] + c * dim + d]; }
911444129b9SMatthew G. Knepley   }
912444129b9SMatthew G. Knepley 
913444129b9SMatthew G. Knepley   // - \vb{\phi}_i \cdot \vu{z} \frac{p^{th}}{T^2 F^2} \psi_j
914444129b9SMatthew G. Knepley   g0[gdir] -= p_th / PetscSqr(u[uOff[TEMP]] * F);
915444129b9SMatthew G. Knepley }
916444129b9SMatthew G. Knepley 
917*9371c9d4SSatish Balay static void g0_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
918444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
919444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
920444129b9SMatthew G. Knepley   const PetscInt  Nc   = dim;
921444129b9SMatthew G. Knepley   PetscInt        c, d;
922444129b9SMatthew G. Knepley 
923444129b9SMatthew G. Knepley   // \vb{\phi}_i \cdot S \rho \psi_j
924*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g0[d * dim + d] = S * p_th / u[uOff[TEMP]] * u_tShift; }
925444129b9SMatthew G. Knepley 
926444129b9SMatthew G. Knepley   // \phi^c_i \cdot \rho \frac{\partial u^c}{\partial x^d} \psi^d_j
927444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
928*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { g0[c * Nc + d] += p_th / u[uOff[TEMP]] * u_x[uOff_x[VEL] + c * Nc + d]; }
929444129b9SMatthew G. Knepley   }
930444129b9SMatthew G. Knepley }
931444129b9SMatthew G. Knepley 
932*9371c9d4SSatish Balay static void g1_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
933444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
934444129b9SMatthew G. Knepley   const PetscInt  NcI  = dim;
935444129b9SMatthew G. Knepley   const PetscInt  NcJ  = dim;
936444129b9SMatthew G. Knepley   PetscInt        c, d, e;
937444129b9SMatthew G. Knepley 
938444129b9SMatthew G. Knepley   // \phi^c_i \rho u^e \frac{\partial \psi^d_j}{\partial x^e}
939444129b9SMatthew G. Knepley   for (c = 0; c < NcI; ++c) {
940444129b9SMatthew G. Knepley     for (d = 0; d < NcJ; ++d) {
941444129b9SMatthew G. Knepley       for (e = 0; e < dim; ++e) {
942*9371c9d4SSatish Balay         if (c == d) { g1[(c * NcJ + d) * dim + e] += p_th / u[uOff[TEMP]] * u[uOff[VEL] + e]; }
943444129b9SMatthew G. Knepley       }
944444129b9SMatthew G. Knepley     }
945444129b9SMatthew G. Knepley   }
946444129b9SMatthew G. Knepley }
947444129b9SMatthew G. Knepley 
948*9371c9d4SSatish Balay static void g3_conduct_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
949444129b9SMatthew G. Knepley   const PetscReal Re = PetscRealPart(constants[REYNOLDS]);
950444129b9SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[MU]);
951444129b9SMatthew G. Knepley   const PetscInt  Nc = dim;
952444129b9SMatthew G. Knepley   PetscInt        c, d;
953444129b9SMatthew G. Knepley 
954444129b9SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
955444129b9SMatthew G. Knepley     for (d = 0; d < dim; ++d) {
956444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^c_i}{\partial x^d}
957444129b9SMatthew G. Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += mu / Re; // gradU
958444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} \mu/Re \frac{\partial \psi^d_i}{\partial x^c}
959444129b9SMatthew G. Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += mu / Re; // gradU transpose
960444129b9SMatthew G. Knepley       // \frac{\partial \phi^c_i}{\partial x^d} -2/3 \mu/Re \frac{\partial \psi^d_i}{\partial x^c}
961444129b9SMatthew G. Knepley       g3[((c * Nc + d) * dim + c) * dim + d] -= 2.0 / 3.0 * mu / Re;
962444129b9SMatthew G. Knepley     }
963444129b9SMatthew G. Knepley   }
964444129b9SMatthew G. Knepley }
965444129b9SMatthew G. Knepley 
966*9371c9d4SSatish Balay static void g2_conduct_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[]) {
967444129b9SMatthew G. Knepley   PetscInt d;
968*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g2[d * dim + d] = -1.0; }
969444129b9SMatthew G. Knepley }
970444129b9SMatthew G. Knepley 
971*9371c9d4SSatish Balay static void g0_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
972a712f3bbSMatthew G. Knepley   g0[0] = u_tShift;
973649ef022SMatthew Knepley }
974649ef022SMatthew Knepley 
975*9371c9d4SSatish Balay static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
976649ef022SMatthew Knepley   PetscInt d;
977649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2] + d];
978649ef022SMatthew Knepley }
979649ef022SMatthew Knepley 
980*9371c9d4SSatish Balay static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
981649ef022SMatthew Knepley   PetscInt d;
982649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0] + d];
983649ef022SMatthew Knepley }
984649ef022SMatthew Knepley 
985*9371c9d4SSatish Balay static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
986444129b9SMatthew G. Knepley   const PetscReal alpha = PetscRealPart(constants[ALPHA]);
987649ef022SMatthew Knepley   PetscInt        d;
988649ef022SMatthew Knepley 
989649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = alpha;
990649ef022SMatthew Knepley }
991649ef022SMatthew Knepley 
992*9371c9d4SSatish Balay static void g0_conduct_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
993444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
994444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
995444129b9SMatthew G. Knepley   PetscInt        d;
996444129b9SMatthew G. Knepley 
997444129b9SMatthew G. Knepley   // \phi_i \frac{C_p p^{th}}{T} \nabla T \cdot \psi_j
998*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g0[d] = c_p * p_th / u[uOff[TEMP]] * u_x[uOff_x[TEMP] + d]; }
999444129b9SMatthew G. Knepley }
1000444129b9SMatthew G. Knepley 
1001*9371c9d4SSatish Balay static void g0_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
1002444129b9SMatthew G. Knepley   const PetscReal S    = PetscRealPart(constants[STROUHAL]);
1003444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1004444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
1005444129b9SMatthew G. Knepley   PetscInt        d;
1006444129b9SMatthew G. Knepley 
1007444129b9SMatthew G. Knepley   // \psi_i C_p S p^{th}\T \psi_{j}
1008444129b9SMatthew G. Knepley   g0[0] += c_p * S * p_th / u[uOff[TEMP]] * u_tShift;
1009444129b9SMatthew G. Knepley   // - \phi_i C_p S p^{th}/T^2 T_t \psi_j
1010444129b9SMatthew G. Knepley   g0[0] -= c_p * S * p_th / PetscSqr(u[uOff[TEMP]]) * u_t[uOff[TEMP]];
1011444129b9SMatthew G. Knepley   // - \phi_i C_p p^{th}/T^2 \vb{u} \cdot \nabla T \psi_j
1012*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g0[0] -= c_p * p_th / PetscSqr(u[uOff[TEMP]]) * u[uOff[VEL] + d] * u_x[uOff_x[TEMP] + d]; }
1013444129b9SMatthew G. Knepley }
1014444129b9SMatthew G. Knepley 
1015*9371c9d4SSatish Balay static void g1_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[]) {
1016444129b9SMatthew G. Knepley   const PetscReal p_th = PetscRealPart(constants[P_TH]);
1017444129b9SMatthew G. Knepley   const PetscReal c_p  = PetscRealPart(constants[C_P]);
1018444129b9SMatthew G. Knepley   PetscInt        d;
1019444129b9SMatthew G. Knepley 
1020444129b9SMatthew G. Knepley   // \phi_i C_p p^{th}/T \vb{u} \cdot \nabla \psi_j
1021*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g1[d] += c_p * p_th / u[uOff[TEMP]] * u[uOff[VEL] + d]; }
1022444129b9SMatthew G. Knepley }
1023444129b9SMatthew G. Knepley 
1024*9371c9d4SSatish Balay static void g3_conduct_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
1025444129b9SMatthew G. Knepley   const PetscReal Pe = PetscRealPart(constants[PECLET]);
1026444129b9SMatthew G. Knepley   const PetscReal k  = PetscRealPart(constants[K]);
1027444129b9SMatthew G. Knepley   PetscInt        d;
1028444129b9SMatthew G. Knepley 
1029444129b9SMatthew G. Knepley   // \nabla \phi_i \frac{k}{Pe} \nabla \phi_j
1030*9371c9d4SSatish Balay   for (d = 0; d < dim; ++d) { g3[d * dim + d] = k / Pe; }
1031444129b9SMatthew G. Knepley }
1032444129b9SMatthew G. Knepley 
1033*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
1034444129b9SMatthew G. Knepley   PetscInt mod, sol;
1035649ef022SMatthew Knepley 
1036649ef022SMatthew Knepley   PetscFunctionBeginUser;
1037444129b9SMatthew G. Knepley   options->modType      = MOD_INCOMPRESSIBLE;
1038649ef022SMatthew Knepley   options->solType      = SOL_QUADRATIC;
1039444129b9SMatthew G. Knepley   options->hasNullSpace = PETSC_TRUE;
1040367970cfSMatthew G. Knepley   options->dmCell       = NULL;
1041649ef022SMatthew Knepley 
1042d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Low Mach flow Problem Options", "DMPLEX");
1043444129b9SMatthew G. Knepley   mod = options->modType;
10449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-mod_type", "The model type", "ex76.c", modTypes, NUM_MOD_TYPES, modTypes[options->modType], &mod, NULL));
1045444129b9SMatthew G. Knepley   options->modType = (ModType)mod;
1046649ef022SMatthew Knepley   sol              = options->solType;
10479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex76.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL));
1048649ef022SMatthew Knepley   options->solType = (SolType)sol;
1049d0609cedSBarry Smith   PetscOptionsEnd();
1050649ef022SMatthew Knepley   PetscFunctionReturn(0);
1051649ef022SMatthew Knepley }
1052649ef022SMatthew Knepley 
1053*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(DM dm, AppCtx *user) {
1054649ef022SMatthew Knepley   PetscBag   bag;
1055649ef022SMatthew Knepley   Parameter *p;
1056444129b9SMatthew G. Knepley   PetscReal  dir;
1057444129b9SMatthew G. Knepley   PetscInt   dim;
1058649ef022SMatthew Knepley 
1059649ef022SMatthew Knepley   PetscFunctionBeginUser;
10609566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1061444129b9SMatthew G. Knepley   dir = (PetscReal)(dim - 1);
1062649ef022SMatthew Knepley   /* setup PETSc parameter bag */
10639566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&p));
10649566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(user->bag, "par", "Low Mach flow parameters"));
1065649ef022SMatthew Knepley   bag = user->bag;
10669566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->Strouhal, 1.0, "S", "Strouhal number"));
10679566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->Froude, 1.0, "Fr", "Froude number"));
10689566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->Reynolds, 1.0, "Re", "Reynolds number"));
10699566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->Peclet, 1.0, "Pe", "Peclet number"));
10709566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->p_th, 1.0, "p_th", "Thermodynamic pressure"));
10719566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->mu, 1.0, "mu", "Dynamic viscosity"));
10729566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity"));
10739566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->c_p, 1.0, "c_p", "Specific heat at constant pressure"));
10749566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->k, 1.0, "k", "Thermal conductivity"));
10759566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity"));
10769566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->T_in, 1.0, "T_in", "Inlet temperature"));
10779566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->g_dir, dir, "g_dir", "Gravity direction"));
10789566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->epsilon, 1.0, "epsilon", "Perturbation strength"));
1079649ef022SMatthew Knepley   PetscFunctionReturn(0);
1080649ef022SMatthew Knepley }
1081649ef022SMatthew Knepley 
1082*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
1083649ef022SMatthew Knepley   PetscFunctionBeginUser;
10849566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
10859566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
10869566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
10879566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1088649ef022SMatthew Knepley   PetscFunctionReturn(0);
1089649ef022SMatthew Knepley }
1090649ef022SMatthew Knepley 
1091*9371c9d4SSatish Balay static PetscErrorCode UniformBoundaryConditions(DM dm, DMLabel label, PetscSimplePointFunc exactFuncs[], PetscSimplePointFunc exactFuncs_t[], AppCtx *user) {
1092444129b9SMatthew G. Knepley   PetscDS  ds;
1093444129b9SMatthew G. Knepley   PetscInt id;
1094444129b9SMatthew G. Knepley   void    *ctx;
1095444129b9SMatthew G. Knepley 
1096444129b9SMatthew G. Knepley   PetscFunctionBeginUser;
10979566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
10989566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, &ctx));
1099444129b9SMatthew G. Knepley   id = 3;
11009566063dSJacob Faibussowitsch   PetscCall(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));
1101444129b9SMatthew G. Knepley   id = 1;
11029566063dSJacob Faibussowitsch   PetscCall(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));
1103444129b9SMatthew G. Knepley   id = 2;
11049566063dSJacob Faibussowitsch   PetscCall(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));
1105444129b9SMatthew G. Knepley   id = 4;
11069566063dSJacob Faibussowitsch   PetscCall(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));
1107444129b9SMatthew G. Knepley   id = 3;
11089566063dSJacob Faibussowitsch   PetscCall(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));
1109444129b9SMatthew G. Knepley   id = 1;
11109566063dSJacob Faibussowitsch   PetscCall(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));
1111444129b9SMatthew G. Knepley   id = 2;
11129566063dSJacob Faibussowitsch   PetscCall(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));
1113444129b9SMatthew G. Knepley   id = 4;
11149566063dSJacob Faibussowitsch   PetscCall(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));
1115444129b9SMatthew G. Knepley   PetscFunctionReturn(0);
1116444129b9SMatthew G. Knepley }
1117444129b9SMatthew G. Knepley 
1118*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
111945480ffeSMatthew G. Knepley   PetscSimplePointFunc exactFuncs[3];
112045480ffeSMatthew G. Knepley   PetscSimplePointFunc exactFuncs_t[3];
1121444129b9SMatthew G. Knepley   PetscDS              ds;
1122444129b9SMatthew G. Knepley   PetscWeakForm        wf;
112345480ffeSMatthew G. Knepley   DMLabel              label;
1124649ef022SMatthew Knepley   Parameter           *ctx;
1125444129b9SMatthew G. Knepley   PetscInt             id, bd;
1126649ef022SMatthew Knepley 
1127649ef022SMatthew Knepley   PetscFunctionBeginUser;
11289566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
11299566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
11309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetWeakForm(ds, &wf));
1131444129b9SMatthew G. Knepley 
1132444129b9SMatthew G. Knepley   switch (user->modType) {
1133444129b9SMatthew G. Knepley   case MOD_INCOMPRESSIBLE:
11349566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, VEL, f0_v, f1_v));
11359566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, PRES, f0_q, NULL));
11369566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, TEMP, f0_w, f1_w));
1137444129b9SMatthew G. Knepley 
11389566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, VEL, VEL, g0_vu, g1_vu, NULL, g3_vu));
11399566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, VEL, PRES, NULL, NULL, g2_vp, NULL));
11409566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, PRES, VEL, NULL, g1_qu, NULL, NULL));
11419566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, TEMP, VEL, g0_wu, NULL, NULL, NULL));
11429566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, TEMP, TEMP, g0_wT, g1_wT, NULL, g3_wT));
1143444129b9SMatthew G. Knepley 
1144649ef022SMatthew Knepley     switch (user->solType) {
1145649ef022SMatthew Knepley     case SOL_QUADRATIC:
11469566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_quadratic_v, 0, NULL));
11479566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_quadratic_w, 0, NULL));
1148649ef022SMatthew Knepley 
1149444129b9SMatthew G. Knepley       exactFuncs[VEL]    = quadratic_u;
1150444129b9SMatthew G. Knepley       exactFuncs[PRES]   = quadratic_p;
1151444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = quadratic_T;
1152444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = quadratic_u_t;
1153444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = NULL;
1154444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = quadratic_T_t;
1155444129b9SMatthew G. Knepley 
11569566063dSJacob Faibussowitsch       PetscCall(UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user));
1157649ef022SMatthew Knepley       break;
1158649ef022SMatthew Knepley     case SOL_CUBIC:
11599566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_cubic_v, 0, NULL));
11609566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_w, 0, NULL));
1161649ef022SMatthew Knepley 
1162444129b9SMatthew G. Knepley       exactFuncs[VEL]    = cubic_u;
1163444129b9SMatthew G. Knepley       exactFuncs[PRES]   = cubic_p;
1164444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = cubic_T;
1165444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = cubic_u_t;
1166444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = NULL;
1167444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = cubic_T_t;
1168444129b9SMatthew G. Knepley 
11699566063dSJacob Faibussowitsch       PetscCall(UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user));
1170649ef022SMatthew Knepley       break;
1171649ef022SMatthew Knepley     case SOL_CUBIC_TRIG:
11729566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_cubic_trig_v, 0, NULL));
11739566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_cubic_trig_w, 0, NULL));
1174649ef022SMatthew Knepley 
1175444129b9SMatthew G. Knepley       exactFuncs[VEL]    = cubic_trig_u;
1176444129b9SMatthew G. Knepley       exactFuncs[PRES]   = cubic_trig_p;
1177444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = cubic_trig_T;
1178444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = cubic_trig_u_t;
1179444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = NULL;
1180444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = cubic_trig_T_t;
1181444129b9SMatthew G. Knepley 
11829566063dSJacob Faibussowitsch       PetscCall(UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user));
1183649ef022SMatthew Knepley       break;
1184606d57d4SMatthew G. Knepley     case SOL_TAYLOR_GREEN:
11859566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_taylor_green_w, 0, NULL));
1186606d57d4SMatthew G. Knepley 
1187444129b9SMatthew G. Knepley       exactFuncs[VEL]    = taylor_green_u;
1188444129b9SMatthew G. Knepley       exactFuncs[PRES]   = taylor_green_p;
1189444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = taylor_green_T;
1190444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = taylor_green_u_t;
1191444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = taylor_green_p_t;
1192444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = taylor_green_T_t;
1193444129b9SMatthew G. Knepley 
11949566063dSJacob Faibussowitsch       PetscCall(UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user));
1195606d57d4SMatthew G. Knepley       break;
119663a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
1197649ef022SMatthew Knepley     }
1198444129b9SMatthew G. Knepley     break;
1199444129b9SMatthew G. Knepley   case MOD_CONDUCTING:
12009566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, VEL, f0_conduct_v, f1_conduct_v));
12019566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, PRES, f0_conduct_q, NULL));
12029566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, TEMP, f0_conduct_w, f1_conduct_w));
1203649ef022SMatthew Knepley 
12049566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, VEL, VEL, g0_conduct_vu, g1_conduct_vu, NULL, g3_conduct_vu));
12059566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, VEL, PRES, NULL, NULL, g2_conduct_vp, NULL));
12069566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, VEL, TEMP, g0_conduct_vT, NULL, NULL, NULL));
12079566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, PRES, VEL, g0_conduct_qu, g1_conduct_qu, NULL, NULL));
12089566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, PRES, TEMP, g0_conduct_qT, g1_conduct_qT, NULL, NULL));
12099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, TEMP, VEL, g0_conduct_wu, NULL, NULL, NULL));
12109566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, TEMP, TEMP, g0_conduct_wT, g1_conduct_wT, NULL, g3_conduct_wT));
1211649ef022SMatthew Knepley 
1212444129b9SMatthew G. Knepley     switch (user->solType) {
1213444129b9SMatthew G. Knepley     case SOL_QUADRATIC:
12149566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_quadratic_v, 0, NULL));
12159566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_quadratic_q, 0, NULL));
12169566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_quadratic_w, 0, NULL));
1217444129b9SMatthew G. Knepley 
1218444129b9SMatthew G. Knepley       exactFuncs[VEL]    = quadratic_u;
1219444129b9SMatthew G. Knepley       exactFuncs[PRES]   = quadratic_p;
1220444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = quadratic_T;
1221444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = quadratic_u_t;
1222444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = NULL;
1223444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = quadratic_T_t;
1224444129b9SMatthew G. Knepley 
12259566063dSJacob Faibussowitsch       PetscCall(UniformBoundaryConditions(dm, label, exactFuncs, exactFuncs_t, user));
1226444129b9SMatthew G. Knepley       break;
1227444129b9SMatthew G. Knepley     case SOL_PIPE:
1228444129b9SMatthew G. Knepley       user->hasNullSpace = PETSC_FALSE;
12299566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_pipe_v, 0, NULL));
12309566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_pipe_q, 0, NULL));
12319566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_pipe_w, 0, NULL));
1232444129b9SMatthew G. Knepley 
1233444129b9SMatthew G. Knepley       exactFuncs[VEL]    = pipe_u;
1234444129b9SMatthew G. Knepley       exactFuncs[PRES]   = pipe_p;
1235444129b9SMatthew G. Knepley       exactFuncs[TEMP]   = pipe_T;
1236444129b9SMatthew G. Knepley       exactFuncs_t[VEL]  = pipe_u_t;
1237444129b9SMatthew G. Knepley       exactFuncs_t[PRES] = pipe_p_t;
1238444129b9SMatthew G. Knepley       exactFuncs_t[TEMP] = pipe_T_t;
1239444129b9SMatthew G. Knepley 
12409566063dSJacob Faibussowitsch       PetscCall(PetscBagGetData(user->bag, (void **)&ctx));
1241444129b9SMatthew G. Knepley       id = 2;
12429566063dSJacob Faibussowitsch       PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
12439566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
12449566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL));
1245444129b9SMatthew G. Knepley       id = 4;
12469566063dSJacob Faibussowitsch       PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
12479566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
12489566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_v, 0, NULL));
1249444129b9SMatthew G. Knepley       id = 4;
12509566063dSJacob Faibussowitsch       PetscCall(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));
1251444129b9SMatthew G. Knepley       id = 3;
12529566063dSJacob Faibussowitsch       PetscCall(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));
12539566063dSJacob Faibussowitsch       PetscCall(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));
1254444129b9SMatthew G. Knepley       id = 1;
12559566063dSJacob Faibussowitsch       PetscCall(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));
12569566063dSJacob Faibussowitsch       PetscCall(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));
1257444129b9SMatthew G. Knepley       break;
1258367970cfSMatthew G. Knepley     case SOL_PIPE_WIGGLY:
1259367970cfSMatthew G. Knepley       user->hasNullSpace = PETSC_FALSE;
12609566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, VEL, 0, 1, f0_conduct_pipe_wiggly_v, 0, NULL));
12619566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, PRES, 0, 1, f0_conduct_pipe_wiggly_q, 0, NULL));
12629566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, TEMP, 0, 1, f0_conduct_pipe_wiggly_w, 0, NULL));
1263367970cfSMatthew G. Knepley 
1264367970cfSMatthew G. Knepley       exactFuncs[VEL]    = pipe_wiggly_u;
1265367970cfSMatthew G. Knepley       exactFuncs[PRES]   = pipe_wiggly_p;
1266367970cfSMatthew G. Knepley       exactFuncs[TEMP]   = pipe_wiggly_T;
1267367970cfSMatthew G. Knepley       exactFuncs_t[VEL]  = pipe_wiggly_u_t;
1268367970cfSMatthew G. Knepley       exactFuncs_t[PRES] = pipe_wiggly_p_t;
1269367970cfSMatthew G. Knepley       exactFuncs_t[TEMP] = pipe_wiggly_T_t;
1270367970cfSMatthew G. Knepley 
12719566063dSJacob Faibussowitsch       PetscCall(PetscBagGetData(user->bag, (void **)&ctx));
1272367970cfSMatthew G. Knepley       id = 2;
12739566063dSJacob Faibussowitsch       PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
12749566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
12759566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_wiggly_v, 0, NULL));
1276367970cfSMatthew G. Knepley       id = 4;
12779566063dSJacob Faibussowitsch       PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall", label, 1, &id, 0, 0, NULL, NULL, NULL, ctx, &bd));
12789566063dSJacob Faibussowitsch       PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
12799566063dSJacob Faibussowitsch       PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, VEL, 0, 0, f0_conduct_bd_pipe_wiggly_v, 0, NULL));
1280367970cfSMatthew G. Knepley       id = 4;
12819566063dSJacob Faibussowitsch       PetscCall(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));
1282367970cfSMatthew G. Knepley       id = 3;
12839566063dSJacob Faibussowitsch       PetscCall(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));
12849566063dSJacob Faibussowitsch       PetscCall(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));
1285367970cfSMatthew G. Knepley       id = 1;
12869566063dSJacob Faibussowitsch       PetscCall(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));
12879566063dSJacob Faibussowitsch       PetscCall(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));
1288367970cfSMatthew G. Knepley       break;
128963a3b9bcSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
1290444129b9SMatthew G. Knepley     }
1291444129b9SMatthew G. Knepley     break;
129263a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)ds), PETSC_ERR_ARG_WRONG, "Unsupported model type: %s (%d)", modTypes[PetscMin(user->modType, NUM_MOD_TYPES)], user->modType);
1293444129b9SMatthew G. Knepley   }
1294649ef022SMatthew Knepley   /* Setup constants */
1295649ef022SMatthew Knepley   {
1296649ef022SMatthew Knepley     Parameter  *param;
1297367970cfSMatthew G. Knepley     PetscScalar constants[13];
1298649ef022SMatthew Knepley 
12999566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **)&param));
1300649ef022SMatthew Knepley 
1301444129b9SMatthew G. Knepley     constants[STROUHAL] = param->Strouhal;
1302444129b9SMatthew G. Knepley     constants[FROUDE]   = param->Froude;
1303444129b9SMatthew G. Knepley     constants[REYNOLDS] = param->Reynolds;
1304444129b9SMatthew G. Knepley     constants[PECLET]   = param->Peclet;
1305444129b9SMatthew G. Knepley     constants[P_TH]     = param->p_th;
1306444129b9SMatthew G. Knepley     constants[MU]       = param->mu;
1307444129b9SMatthew G. Knepley     constants[NU]       = param->nu;
1308444129b9SMatthew G. Knepley     constants[C_P]      = param->c_p;
1309444129b9SMatthew G. Knepley     constants[K]        = param->k;
1310444129b9SMatthew G. Knepley     constants[ALPHA]    = param->alpha;
1311444129b9SMatthew G. Knepley     constants[T_IN]     = param->T_in;
1312444129b9SMatthew G. Knepley     constants[G_DIR]    = param->g_dir;
1313367970cfSMatthew G. Knepley     constants[EPSILON]  = param->epsilon;
13149566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(ds, 13, constants));
1315649ef022SMatthew Knepley   }
1316649ef022SMatthew Knepley 
13179566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&ctx));
13189566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, VEL, exactFuncs[VEL], ctx));
13199566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, PRES, exactFuncs[PRES], ctx));
13209566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, TEMP, exactFuncs[TEMP], ctx));
13219566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, VEL, exactFuncs_t[VEL], ctx));
13229566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, PRES, exactFuncs_t[PRES], ctx));
13239566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, TEMP, exactFuncs_t[TEMP], ctx));
1324649ef022SMatthew Knepley   PetscFunctionReturn(0);
1325649ef022SMatthew Knepley }
1326649ef022SMatthew Knepley 
1327*9371c9d4SSatish Balay static PetscErrorCode CreateCellDM(DM dm, AppCtx *user) {
1328367970cfSMatthew G. Knepley   PetscFE        fe, fediv;
1329367970cfSMatthew G. Knepley   DMPolytopeType ct;
1330367970cfSMatthew G. Knepley   PetscInt       dim, cStart;
1331367970cfSMatthew G. Knepley   PetscBool      simplex;
1332367970cfSMatthew G. Knepley 
1333367970cfSMatthew G. Knepley   PetscFunctionBeginUser;
13349566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
13369566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1337367970cfSMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
1338367970cfSMatthew G. Knepley 
13399566063dSJacob Faibussowitsch   PetscCall(DMGetField(dm, VEL, NULL, (PetscObject *)&fe));
13409566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "div_", PETSC_DEFAULT, &fediv));
13419566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe, fediv));
13429566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fediv, "divergence"));
1343367970cfSMatthew G. Knepley 
13449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user->dmCell));
13459566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &user->dmCell));
13469566063dSJacob Faibussowitsch   PetscCall(DMSetField(user->dmCell, 0, NULL, (PetscObject)fediv));
13479566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(user->dmCell));
13489566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fediv));
1349367970cfSMatthew G. Knepley   PetscFunctionReturn(0);
1350367970cfSMatthew G. Knepley }
1351367970cfSMatthew G. Knepley 
1352*9371c9d4SSatish Balay static PetscErrorCode GetCellDM(DM dm, AppCtx *user, DM *dmCell) {
1353367970cfSMatthew G. Knepley   PetscInt cStart, cEnd, cellStart = -1, cellEnd = -1;
1354367970cfSMatthew G. Knepley 
1355367970cfSMatthew G. Knepley   PetscFunctionBeginUser;
13569566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSimplexOrBoxCells(dm, 0, &cStart, &cEnd));
13579566063dSJacob Faibussowitsch   if (user->dmCell) PetscCall(DMPlexGetSimplexOrBoxCells(user->dmCell, 0, &cellStart, &cellEnd));
13589566063dSJacob Faibussowitsch   if (cStart != cellStart || cEnd != cellEnd) PetscCall(CreateCellDM(dm, user));
1359367970cfSMatthew G. Knepley   *dmCell = user->dmCell;
1360367970cfSMatthew G. Knepley   PetscFunctionReturn(0);
1361367970cfSMatthew G. Knepley }
1362367970cfSMatthew G. Knepley 
1363*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) {
1364649ef022SMatthew Knepley   DM             cdm = dm;
1365367970cfSMatthew G. Knepley   PetscFE        fe[3];
1366649ef022SMatthew Knepley   Parameter     *param;
1367649ef022SMatthew Knepley   DMPolytopeType ct;
1368649ef022SMatthew Knepley   PetscInt       dim, cStart;
1369649ef022SMatthew Knepley   PetscBool      simplex;
1370649ef022SMatthew Knepley 
1371649ef022SMatthew Knepley   PetscFunctionBeginUser;
13729566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13739566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
13749566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
1375649ef022SMatthew Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
1376649ef022SMatthew Knepley   /* Create finite element */
13779566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
13789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
1379649ef022SMatthew Knepley 
13809566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
13819566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
13829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
1383649ef022SMatthew Knepley 
13849566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]));
13859566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
13869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[2], "temperature"));
1387649ef022SMatthew Knepley 
1388649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
13899566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, VEL, NULL, (PetscObject)fe[VEL]));
13909566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, PRES, NULL, (PetscObject)fe[PRES]));
13919566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, TEMP, NULL, (PetscObject)fe[TEMP]));
13929566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
13939566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
13949566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&param));
1395649ef022SMatthew Knepley   while (cdm) {
13969566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
13979566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
1398649ef022SMatthew Knepley   }
13999566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[VEL]));
14009566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[PRES]));
14019566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[TEMP]));
1402649ef022SMatthew Knepley 
1403444129b9SMatthew G. Knepley   if (user->hasNullSpace) {
1404649ef022SMatthew Knepley     PetscObject  pressure;
1405649ef022SMatthew Knepley     MatNullSpace nullspacePres;
1406649ef022SMatthew Knepley 
14079566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm, PRES, NULL, &pressure));
14089566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres));
14099566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nullspacePres));
14109566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullspacePres));
1411649ef022SMatthew Knepley   }
1412649ef022SMatthew Knepley   PetscFunctionReturn(0);
1413649ef022SMatthew Knepley }
1414649ef022SMatthew Knepley 
1415*9371c9d4SSatish Balay static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) {
1416649ef022SMatthew Knepley   Vec vec;
1417649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
1418649ef022SMatthew Knepley 
1419649ef022SMatthew Knepley   PetscFunctionBeginUser;
142063a3b9bcSJacob Faibussowitsch   PetscCheck(ofield == PRES, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index %" PetscInt_FMT ", not %" PetscInt_FMT, PRES, ofield);
1421649ef022SMatthew Knepley   funcs[nfield] = constant;
14229566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &vec));
14239566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
14249566063dSJacob Faibussowitsch   PetscCall(VecNormalize(vec, NULL));
14259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)vec, "Pressure Null Space"));
14269566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view"));
14279566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace));
14289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec));
1429649ef022SMatthew Knepley   PetscFunctionReturn(0);
1430649ef022SMatthew Knepley }
1431649ef022SMatthew Knepley 
1432*9371c9d4SSatish Balay static PetscErrorCode RemoveDiscretePressureNullspace_Private(TS ts, Vec u) {
1433649ef022SMatthew Knepley   DM           dm;
1434444129b9SMatthew G. Knepley   AppCtx      *user;
1435649ef022SMatthew Knepley   MatNullSpace nullsp;
1436649ef022SMatthew Knepley 
14377510d9b0SBarry Smith   PetscFunctionBeginUser;
14389566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14399566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &user));
1440444129b9SMatthew G. Knepley   if (!user->hasNullSpace) PetscFunctionReturn(0);
14419566063dSJacob Faibussowitsch   PetscCall(CreatePressureNullSpace(dm, 1, 1, &nullsp));
14429566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceRemove(nullsp, u));
14439566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceDestroy(&nullsp));
1444649ef022SMatthew Knepley   PetscFunctionReturn(0);
1445649ef022SMatthew Knepley }
1446649ef022SMatthew Knepley 
1447649ef022SMatthew Knepley /* Make the discrete pressure discretely divergence free */
1448*9371c9d4SSatish Balay static PetscErrorCode RemoveDiscretePressureNullspace(TS ts) {
1449649ef022SMatthew Knepley   Vec u;
1450649ef022SMatthew Knepley 
14517510d9b0SBarry Smith   PetscFunctionBeginUser;
14529566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &u));
14539566063dSJacob Faibussowitsch   PetscCall(RemoveDiscretePressureNullspace_Private(ts, u));
1454649ef022SMatthew Knepley   PetscFunctionReturn(0);
1455649ef022SMatthew Knepley }
1456649ef022SMatthew Knepley 
1457*9371c9d4SSatish Balay static void divergence(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar divu[]) {
1458367970cfSMatthew G. Knepley   PetscInt d;
1459367970cfSMatthew G. Knepley 
1460367970cfSMatthew G. Knepley   divu[0] = 0.;
1461367970cfSMatthew G. Knepley   for (d = 0; d < dim; ++d) divu[0] += u_x[d * dim + d];
1462367970cfSMatthew G. Knepley }
1463367970cfSMatthew G. Knepley 
1464*9371c9d4SSatish Balay static PetscErrorCode SetInitialConditions(TS ts, Vec u) {
1465444129b9SMatthew G. Knepley   AppCtx   *user;
1466649ef022SMatthew Knepley   DM        dm;
1467649ef022SMatthew Knepley   PetscReal t;
1468649ef022SMatthew Knepley 
14697510d9b0SBarry Smith   PetscFunctionBeginUser;
14709566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14719566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
14729566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
14739566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &user));
14749566063dSJacob Faibussowitsch   PetscCall(RemoveDiscretePressureNullspace_Private(ts, u));
1475649ef022SMatthew Knepley   PetscFunctionReturn(0);
1476649ef022SMatthew Knepley }
1477649ef022SMatthew Knepley 
1478*9371c9d4SSatish Balay static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) {
1479649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
1480649ef022SMatthew Knepley   void          *ctxs[3];
1481a712f3bbSMatthew G. Knepley   PetscPointFunc diagnostics[1] = {divergence};
1482367970cfSMatthew G. Knepley   DM             dm, dmCell = NULL;
1483649ef022SMatthew Knepley   PetscDS        ds;
1484a712f3bbSMatthew G. Knepley   Vec            v, divu;
1485a712f3bbSMatthew G. Knepley   PetscReal      ferrors[3], massFlux;
1486649ef022SMatthew Knepley   PetscInt       f;
1487649ef022SMatthew Knepley 
1488649ef022SMatthew Knepley   PetscFunctionBeginUser;
14899566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
14909566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
1491649ef022SMatthew Knepley 
14929566063dSJacob Faibussowitsch   for (f = 0; f < 3; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &exactFuncs[f], &ctxs[f]));
14939566063dSJacob Faibussowitsch   PetscCall(DMComputeL2FieldDiff(dm, crtime, exactFuncs, ctxs, u, ferrors));
14949566063dSJacob Faibussowitsch   PetscCall(GetCellDM(dm, (AppCtx *)ctx, &dmCell));
14959566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmCell, &divu));
14969566063dSJacob Faibussowitsch   PetscCall(DMProjectField(dmCell, crtime, u, diagnostics, INSERT_VALUES, divu));
14979566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(divu, NULL, "-divu_vec_view"));
14989566063dSJacob Faibussowitsch   PetscCall(VecNorm(divu, NORM_2, &massFlux));
14999566063dSJacob Faibussowitsch   PetscCall(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));
1500649ef022SMatthew Knepley 
15019566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
1502649ef022SMatthew Knepley 
15039566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &v));
15049566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, crtime, exactFuncs, ctxs, INSERT_ALL_VALUES, v));
15059566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)v, "Exact Solution"));
15069566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(v, NULL, "-exact_vec_view"));
15079566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &v));
1508649ef022SMatthew Knepley 
15099566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(divu, NULL, "-div_vec_view"));
15109566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmCell, &divu));
1511a712f3bbSMatthew G. Knepley 
1512649ef022SMatthew Knepley   PetscFunctionReturn(0);
1513649ef022SMatthew Knepley }
1514649ef022SMatthew Knepley 
1515*9371c9d4SSatish Balay int main(int argc, char **argv) {
1516649ef022SMatthew Knepley   DM        dm;   /* problem definition */
1517649ef022SMatthew Knepley   TS        ts;   /* timestepper */
1518649ef022SMatthew Knepley   Vec       u;    /* solution */
1519649ef022SMatthew Knepley   AppCtx    user; /* user-defined work context */
1520649ef022SMatthew Knepley   PetscReal t;
1521649ef022SMatthew Knepley 
1522327415f7SBarry Smith   PetscFunctionBeginUser;
15239566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15249566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
15259566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
15269566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
15279566063dSJacob Faibussowitsch   PetscCall(SetupParameters(dm, &user));
15289566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
15299566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
15309566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &user));
1531649ef022SMatthew Knepley   /* Setup problem */
15329566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &user));
15339566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
1534649ef022SMatthew Knepley 
15359566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
15369566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution"));
15379566063dSJacob Faibussowitsch   if (user.hasNullSpace) PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace));
1538649ef022SMatthew Knepley 
15399566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &user));
15409566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &user));
15419566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &user));
15429566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
15439566063dSJacob Faibussowitsch   PetscCall(TSSetPreStep(ts, RemoveDiscretePressureNullspace));
15449566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1545649ef022SMatthew Knepley 
15469566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions)); /* Must come after SetFromOptions() */
15479566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(ts, u));
15489566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
15499566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, 0, t));
15509566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
15519566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MonitorError, &user, NULL));
1552649ef022SMatthew Knepley 
15539566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
15549566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
1555649ef022SMatthew Knepley 
15569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
15579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmCell));
15589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
15599566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
15609566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
15619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1562b122ec5aSJacob Faibussowitsch   return 0;
1563649ef022SMatthew Knepley }
1564649ef022SMatthew Knepley 
1565649ef022SMatthew Knepley /*TEST
1566649ef022SMatthew Knepley 
1567444129b9SMatthew G. Knepley   testset:
1568649ef022SMatthew Knepley     requires: triangle !single
1569444129b9SMatthew G. Knepley     args: -dm_plex_separate_marker \
1570a712f3bbSMatthew G. Knepley           -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \
1571444129b9SMatthew G. Knepley           -snes_error_if_not_converged -snes_convergence_test correct_pressure \
1572649ef022SMatthew Knepley           -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
1573444129b9SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \
1574444129b9SMatthew G. Knepley           -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
1575649ef022SMatthew Knepley             -fieldsplit_0_pc_type lu \
1576649ef022SMatthew Knepley             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1577649ef022SMatthew Knepley 
1578444129b9SMatthew G. Knepley     test:
1579444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1
1580444129b9SMatthew G. Knepley       args: -sol_type quadratic \
1581444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1582444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1583444129b9SMatthew G. Knepley 
1584649ef022SMatthew Knepley     test:
1585649ef022SMatthew Knepley       # Using -dm_refine 5 -convest_num_refine 2 gives L_2 convergence rate: [0.89, 0.011, 1.0]
1586649ef022SMatthew Knepley       suffix: 2d_tri_p2_p1_p1_tconv
1587444129b9SMatthew G. Knepley       args: -sol_type cubic_trig \
1588649ef022SMatthew Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1589444129b9SMatthew G. Knepley             -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1
1590649ef022SMatthew Knepley 
1591649ef022SMatthew Knepley     test:
1592649ef022SMatthew Knepley       # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.5, 1.9]
1593649ef022SMatthew Knepley       suffix: 2d_tri_p2_p1_p1_sconv
1594444129b9SMatthew G. Knepley       args: -sol_type cubic \
1595649ef022SMatthew Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1596444129b9SMatthew G. Knepley             -ts_max_steps 1 -ts_dt 1e-4 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1
1597649ef022SMatthew Knepley 
1598649ef022SMatthew Knepley     test:
1599649ef022SMatthew Knepley       suffix: 2d_tri_p3_p2_p2
1600444129b9SMatthew G. Knepley       args: -sol_type cubic \
1601649ef022SMatthew Knepley             -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
1602444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1603649ef022SMatthew Knepley 
1604606d57d4SMatthew G. Knepley     test:
1605606d57d4SMatthew G. Knepley       # Using -dm_refine 3 -convest_num_refine 3 gives L_2 convergence rate: [3.0, 2.1, 3.1]
1606606d57d4SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_tg_sconv
1607444129b9SMatthew G. Knepley       args: -sol_type taylor_green \
1608606d57d4SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1609444129b9SMatthew G. Knepley             -ts_max_steps 1 -ts_dt 1e-8 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1
1610606d57d4SMatthew G. Knepley 
1611606d57d4SMatthew G. Knepley     test:
1612606d57d4SMatthew G. Knepley       # Using -dm_refine 3 -convest_num_refine 2 gives L_2 convergence rate: [1.2, 1.5, 1.2]
1613606d57d4SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_tg_tconv
1614444129b9SMatthew G. Knepley       args: -sol_type taylor_green \
1615606d57d4SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1616444129b9SMatthew G. Knepley             -ts_max_steps 4 -ts_dt 0.1 -ts_convergence_estimate -convest_num_refine 1
1617444129b9SMatthew G. Knepley 
1618444129b9SMatthew G. Knepley   testset:
1619444129b9SMatthew G. Knepley     requires: triangle !single
1620444129b9SMatthew G. Knepley     args: -dm_plex_separate_marker -mod_type conducting \
1621a712f3bbSMatthew G. Knepley           -div_petscdualspace_lagrange_use_moments -div_petscdualspace_lagrange_moment_order 3 \
1622444129b9SMatthew G. Knepley           -snes_error_if_not_converged -snes_max_linear_solve_fail 5 \
162382894d03SBarry Smith           -ksp_type fgmres -ksp_max_it 2 -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 \
1624444129b9SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 \
1625444129b9SMatthew G. Knepley           -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
1626606d57d4SMatthew G. Knepley             -fieldsplit_0_pc_type lu \
1627606d57d4SMatthew G. Knepley             -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
1628606d57d4SMatthew G. Knepley 
1629444129b9SMatthew G. Knepley     test:
1630444129b9SMatthew G. Knepley       # At this resolution, the rhs is inconsistent on some Newton steps
1631444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p1_conduct
1632444129b9SMatthew G. Knepley       args: -sol_type quadratic \
1633444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
1634444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1 \
1635444129b9SMatthew G. Knepley             -pc_fieldsplit_schur_precondition full \
1636444129b9SMatthew G. Knepley               -fieldsplit_pressure_ksp_max_it 2 -fieldsplit_pressure_pc_type svd
1637444129b9SMatthew G. Knepley 
1638444129b9SMatthew G. Knepley     test:
1639444129b9SMatthew G. Knepley       suffix: 2d_tri_p2_p1_p2_conduct_pipe
1640444129b9SMatthew G. Knepley       args: -sol_type pipe \
1641444129b9SMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 2 \
1642444129b9SMatthew G. Knepley             -dmts_check .001 -ts_max_steps 4 -ts_dt 0.1
1643444129b9SMatthew G. Knepley 
1644367970cfSMatthew G. Knepley     test:
1645367970cfSMatthew G. Knepley       suffix: 2d_tri_p2_p1_p2_conduct_pipe_wiggly_sconv
1646367970cfSMatthew G. Knepley       args: -sol_type pipe_wiggly -Fr 1e10 \
1647367970cfSMatthew G. Knepley             -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 2 \
1648367970cfSMatthew G. Knepley             -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
1649367970cfSMatthew G. Knepley             -ts_max_steps 1 -ts_dt 1e10 \
1650367970cfSMatthew G. Knepley             -ksp_atol 1e-12 -ksp_max_it 300 \
1651367970cfSMatthew G. Knepley               -fieldsplit_pressure_ksp_atol 1e-14
1652367970cfSMatthew G. Knepley 
1653649ef022SMatthew Knepley TEST*/
1654