xref: /petsc/src/snes/tutorials/ex76.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1649ef022SMatthew Knepley static char help[] = "Low Mach Flow in 2d and 3d channels with finite elements.\n\
2649ef022SMatthew Knepley We solve the Low Mach flow problem in a rectangular\n\
3649ef022SMatthew Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4649ef022SMatthew Knepley 
5649ef022SMatthew Knepley /*F
6649ef022SMatthew Knepley This Low Mach flow is a steady-state isoviscous Navier-Stokes flow. We discretize using the
7649ef022SMatthew Knepley finite element method on an unstructured mesh. The weak form equations are
8649ef022SMatthew Knepley 
9649ef022SMatthew Knepley \begin{align*}
10649ef022SMatthew Knepley     < q, \nabla\cdot u >                                                                                     = 0
11649ef022SMatthew Knepley     <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
12649ef022SMatthew Knepley     < w, u \cdot \nabla T > - < \nabla w, \alpha \nabla T > - < w, Q >                                       = 0
13649ef022SMatthew Knepley \end{align*}
14649ef022SMatthew Knepley 
15649ef022SMatthew Knepley where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
16649ef022SMatthew Knepley 
17649ef022SMatthew Knepley For visualization, use
18649ef022SMatthew Knepley 
19649ef022SMatthew Knepley   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
20649ef022SMatthew Knepley F*/
21649ef022SMatthew Knepley 
22649ef022SMatthew Knepley #include <petscdmplex.h>
23649ef022SMatthew Knepley #include <petscsnes.h>
24649ef022SMatthew Knepley #include <petscds.h>
25649ef022SMatthew Knepley #include <petscbag.h>
26649ef022SMatthew Knepley 
27*9371c9d4SSatish Balay typedef enum {
28*9371c9d4SSatish Balay   SOL_QUADRATIC,
29*9371c9d4SSatish Balay   SOL_CUBIC,
30*9371c9d4SSatish Balay   NUM_SOL_TYPES
31*9371c9d4SSatish Balay } SolType;
32649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES + 1] = {"quadratic", "cubic", "unknown"};
33649ef022SMatthew Knepley 
34649ef022SMatthew Knepley typedef struct {
35649ef022SMatthew Knepley   PetscReal nu;    /* Kinematic viscosity */
36649ef022SMatthew Knepley   PetscReal theta; /* Angle of pipe wall to x-axis */
37649ef022SMatthew Knepley   PetscReal alpha; /* Thermal diffusivity */
38649ef022SMatthew Knepley   PetscReal T_in;  /* Inlet temperature*/
39649ef022SMatthew Knepley } Parameter;
40649ef022SMatthew Knepley 
41649ef022SMatthew Knepley typedef struct {
4230602db0SMatthew G. Knepley   PetscBool showError;
4330602db0SMatthew G. Knepley   PetscBag  bag;
44649ef022SMatthew Knepley   SolType   solType;
45649ef022SMatthew Knepley } AppCtx;
46649ef022SMatthew Knepley 
47*9371c9d4SSatish Balay static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
48649ef022SMatthew Knepley   PetscInt d;
49649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
50649ef022SMatthew Knepley   return 0;
51649ef022SMatthew Knepley }
52649ef022SMatthew Knepley 
53*9371c9d4SSatish Balay static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
54649ef022SMatthew Knepley   PetscInt d;
55649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
56649ef022SMatthew Knepley   return 0;
57649ef022SMatthew Knepley }
58649ef022SMatthew Knepley 
59649ef022SMatthew Knepley /*
60649ef022SMatthew Knepley   CASE: quadratic
61649ef022SMatthew Knepley   In 2D we use exact solution:
62649ef022SMatthew Knepley 
63649ef022SMatthew Knepley     u = x^2 + y^2
64649ef022SMatthew Knepley     v = 2x^2 - 2xy
65649ef022SMatthew Knepley     p = x + y - 1
66649ef022SMatthew Knepley     T = x + y
67649ef022SMatthew Knepley     f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -4\nu + 1>
68649ef022SMatthew Knepley     Q = 3x^2 + y^2 - 2xy
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley   so that
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley (1)  \nabla \cdot u  = 2x - 2x = 0
73649ef022SMatthew Knepley 
74649ef022SMatthew Knepley (2)  u \cdot \nabla u - \nu \Delta u + \nabla p - f
75649ef022SMatthew Knepley      = <2x^3 + 4x^2y -2xy^2, 4xy^2 + 2x^2y - 2y^3> -\nu <4, 4> + <1, 1> - <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -         4\nu + 1>  = 0
76649ef022SMatthew Knepley 
77649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0
78649ef022SMatthew Knepley */
79649ef022SMatthew Knepley 
80*9371c9d4SSatish Balay static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
81649ef022SMatthew Knepley   u[0] = X[0] * X[0] + X[1] * X[1];
82649ef022SMatthew Knepley   u[1] = 2.0 * X[0] * X[0] - 2.0 * X[0] * X[1];
83649ef022SMatthew Knepley   return 0;
84649ef022SMatthew Knepley }
85649ef022SMatthew Knepley 
86*9371c9d4SSatish Balay static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
87649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
88649ef022SMatthew Knepley   return 0;
89649ef022SMatthew Knepley }
90649ef022SMatthew Knepley 
91*9371c9d4SSatish Balay static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
92649ef022SMatthew Knepley   T[0] = X[0] + X[1];
93649ef022SMatthew Knepley   return 0;
94649ef022SMatthew Knepley }
95649ef022SMatthew Knepley 
96*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[]) {
97649ef022SMatthew Knepley   PetscInt        c, d;
98649ef022SMatthew Knepley   PetscInt        Nc = dim;
99649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
100649ef022SMatthew Knepley 
101649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
102649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
103649ef022SMatthew Knepley   }
104649ef022SMatthew Knepley   f0[0] -= (2 * X[0] * X[0] * X[0] + 4 * X[0] * X[0] * X[1] - 2 * X[0] * X[1] * X[1] - 4.0 * nu + 1);
105649ef022SMatthew Knepley   f0[1] -= (4 * X[0] * X[1] * X[1] + 2 * X[0] * X[0] * X[1] - 2 * X[1] * X[1] * X[1] - 4.0 * nu + 1);
106649ef022SMatthew Knepley }
107649ef022SMatthew Knepley 
108*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[]) {
109649ef022SMatthew Knepley   PetscInt d;
110649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d];
111649ef022SMatthew Knepley   f0[0] -= (3 * X[0] * X[0] + X[1] * X[1] - 2 * X[0] * X[1]);
112649ef022SMatthew Knepley }
113649ef022SMatthew Knepley 
114649ef022SMatthew Knepley /*
115649ef022SMatthew Knepley   CASE: cubic
116649ef022SMatthew Knepley   In 2D we use exact solution:
117649ef022SMatthew Knepley 
118649ef022SMatthew Knepley     u = x^3 + y^3
119649ef022SMatthew Knepley     v = 2x^3 - 3x^2y
120649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
121649ef022SMatthew Knepley     T = 1/2 x^2 + 1/2 y^2
122649ef022SMatthew Knepley     f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y>
123649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2
124649ef022SMatthew Knepley 
125649ef022SMatthew Knepley   so that
126649ef022SMatthew Knepley 
127649ef022SMatthew Knepley   \nabla \cdot u = 3x^2 - 3x^2 = 0
128649ef022SMatthew Knepley 
129649ef022SMatthew Knepley   u \cdot \nabla u - \nu \Delta u + \nabla p - f
130649ef022SMatthew Knepley   = <3x^5 + 6x^3y^2 - 6x^2y^3, 6x^2y^3 + 3x^4y - 6xy^4> - \nu<6x + 6y, 12x - 6y> + <3x, 3y> - <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y> = 0
131649ef022SMatthew Knepley 
132649ef022SMatthew Knepley   u \cdot \nabla T - \alpha\Delta T - Q = (x^3 + y^3) x + (2x^3 - 3x^2y) y - 2*\alpha - (x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2)   = 0
133649ef022SMatthew Knepley */
134649ef022SMatthew Knepley 
135*9371c9d4SSatish Balay static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx) {
136649ef022SMatthew Knepley   u[0] = X[0] * X[0] * X[0] + X[1] * X[1] * X[1];
137649ef022SMatthew Knepley   u[1] = 2.0 * X[0] * X[0] * X[0] - 3.0 * X[0] * X[0] * X[1];
138649ef022SMatthew Knepley   return 0;
139649ef022SMatthew Knepley }
140649ef022SMatthew Knepley 
141*9371c9d4SSatish Balay static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx) {
142649ef022SMatthew Knepley   p[0] = 3.0 * X[0] * X[0] / 2.0 + 3.0 * X[1] * X[1] / 2.0 - 1.0;
143649ef022SMatthew Knepley   return 0;
144649ef022SMatthew Knepley }
145649ef022SMatthew Knepley 
146*9371c9d4SSatish Balay static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx) {
147649ef022SMatthew Knepley   T[0] = X[0] * X[0] / 2.0 + X[1] * X[1] / 2.0;
148649ef022SMatthew Knepley   return 0;
149649ef022SMatthew Knepley }
150649ef022SMatthew Knepley 
151*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[]) {
152649ef022SMatthew Knepley   PetscInt        c, d;
153649ef022SMatthew Knepley   PetscInt        Nc = dim;
154649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
155649ef022SMatthew Knepley 
156649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
157649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
158649ef022SMatthew Knepley   }
159649ef022SMatthew Knepley   f0[0] -= (3 * X[0] * X[0] * X[0] * X[0] * X[0] + 6 * X[0] * X[0] * X[0] * X[1] * X[1] - 6 * X[0] * X[0] * X[1] * X[1] * X[1] - (6 * X[0] + 6 * X[1]) * nu + 3 * X[0]);
160649ef022SMatthew Knepley   f0[1] -= (6 * X[0] * X[0] * X[1] * X[1] * X[1] + 3 * X[0] * X[0] * X[0] * X[0] * X[1] - 6 * X[0] * X[1] * X[1] * X[1] * X[1] - (12 * X[0] - 6 * X[1]) * nu + 3 * X[1]);
161649ef022SMatthew Knepley }
162649ef022SMatthew Knepley 
163*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[]) {
164649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
165649ef022SMatthew Knepley   PetscInt        d;
166649ef022SMatthew Knepley 
167649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0] + d] * u_x[uOff_x[2] + d];
168649ef022SMatthew Knepley   f0[0] -= (X[0] * X[0] * X[0] * X[0] + X[0] * X[1] * X[1] * X[1] + 2.0 * X[0] * X[0] * X[0] * X[1] - 3.0 * X[0] * X[0] * X[1] * X[1] - 2.0 * alpha);
169649ef022SMatthew Knepley }
170649ef022SMatthew Knepley 
171*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[]) {
172649ef022SMatthew Knepley   PetscInt d;
173649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
174649ef022SMatthew Knepley }
175649ef022SMatthew Knepley 
176*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[]) {
177649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
178649ef022SMatthew Knepley   const PetscInt  Nc = dim;
179649ef022SMatthew Knepley   PetscInt        c, d;
180649ef022SMatthew Knepley 
181649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
182649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
183649ef022SMatthew Knepley       f1[c * dim + d] = nu * (u_x[c * dim + d] + u_x[d * dim + c]);
184649ef022SMatthew Knepley       //f1[c*dim+d] = nu*u_x[c*dim+d];
185649ef022SMatthew Knepley     }
186649ef022SMatthew Knepley     f1[c * dim + c] -= u[uOff[1]];
187649ef022SMatthew Knepley   }
188649ef022SMatthew Knepley }
189649ef022SMatthew Knepley 
190*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[]) {
191649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
192649ef022SMatthew Knepley   PetscInt        d;
193649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha * u_x[uOff_x[2] + d];
194649ef022SMatthew Knepley }
195649ef022SMatthew Knepley 
196649ef022SMatthew Knepley /* Jacobians */
197*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[]) {
198649ef022SMatthew Knepley   PetscInt d;
199649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0;
200649ef022SMatthew Knepley }
201649ef022SMatthew Knepley 
202*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[]) {
203649ef022SMatthew Knepley   const PetscInt Nc = dim;
204649ef022SMatthew Knepley   PetscInt       c, d;
205649ef022SMatthew Knepley 
206649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
207*9371c9d4SSatish Balay     for (d = 0; d < dim; ++d) { g0[c * Nc + d] = u_x[c * Nc + d]; }
208649ef022SMatthew Knepley   }
209649ef022SMatthew Knepley }
210649ef022SMatthew Knepley 
211*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[]) {
212649ef022SMatthew Knepley   PetscInt NcI = dim;
213649ef022SMatthew Knepley   PetscInt NcJ = dim;
214649ef022SMatthew Knepley   PetscInt c, d, e;
215649ef022SMatthew Knepley 
216649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
217649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
218649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
219*9371c9d4SSatish Balay         if (c == d) { g1[(c * NcJ + d) * dim + e] = u[e]; }
220649ef022SMatthew Knepley       }
221649ef022SMatthew Knepley     }
222649ef022SMatthew Knepley   }
223649ef022SMatthew Knepley }
224649ef022SMatthew Knepley 
225*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[]) {
226649ef022SMatthew Knepley   PetscInt d;
227649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0;
228649ef022SMatthew Knepley }
229649ef022SMatthew Knepley 
230*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[]) {
231649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
232649ef022SMatthew Knepley   const PetscInt  Nc = dim;
233649ef022SMatthew Knepley   PetscInt        c, d;
234649ef022SMatthew Knepley 
235649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
236649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
237649ef022SMatthew Knepley       g3[((c * Nc + c) * dim + d) * dim + d] += nu; // gradU
238649ef022SMatthew Knepley       g3[((c * Nc + d) * dim + d) * dim + c] += nu; // gradU transpose
239649ef022SMatthew Knepley     }
240649ef022SMatthew Knepley   }
241649ef022SMatthew Knepley }
242649ef022SMatthew Knepley 
243*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[]) {
244649ef022SMatthew Knepley   PetscInt d;
245649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2] + d];
246649ef022SMatthew Knepley }
247649ef022SMatthew Knepley 
248*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[]) {
249649ef022SMatthew Knepley   PetscInt d;
250649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0] + d];
251649ef022SMatthew Knepley }
252649ef022SMatthew Knepley 
253*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[]) {
254649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
255649ef022SMatthew Knepley   PetscInt        d;
256649ef022SMatthew Knepley 
257649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = alpha;
258649ef022SMatthew Knepley }
259649ef022SMatthew Knepley 
260*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
26130602db0SMatthew G. Knepley   PetscInt sol;
262649ef022SMatthew Knepley 
263649ef022SMatthew Knepley   PetscFunctionBeginUser;
264649ef022SMatthew Knepley   options->solType   = SOL_QUADRATIC;
265649ef022SMatthew Knepley   options->showError = PETSC_FALSE;
266d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
267649ef022SMatthew Knepley   sol = options->solType;
2689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL));
269649ef022SMatthew Knepley   options->solType = (SolType)sol;
2709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL));
271d0609cedSBarry Smith   PetscOptionsEnd();
272649ef022SMatthew Knepley   PetscFunctionReturn(0);
273649ef022SMatthew Knepley }
274649ef022SMatthew Knepley 
275*9371c9d4SSatish Balay static PetscErrorCode SetupParameters(AppCtx *user) {
276649ef022SMatthew Knepley   PetscBag   bag;
277649ef022SMatthew Knepley   Parameter *p;
278649ef022SMatthew Knepley 
279649ef022SMatthew Knepley   PetscFunctionBeginUser;
280649ef022SMatthew Knepley   /* setup PETSc parameter bag */
2819566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&p));
2829566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters"));
283649ef022SMatthew Knepley   bag = user->bag;
2849566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->nu, 1.0, "nu", "Kinematic viscosity"));
2859566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0, "alpha", "Thermal diffusivity"));
2869566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->theta, 0.0, "theta", "Angle of pipe wall to x-axis"));
287649ef022SMatthew Knepley   PetscFunctionReturn(0);
288649ef022SMatthew Knepley }
289649ef022SMatthew Knepley 
290*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
291649ef022SMatthew Knepley   PetscFunctionBeginUser;
2929566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2939566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
2949566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
295649ef022SMatthew Knepley   {
296649ef022SMatthew Knepley     Parameter   *param;
297649ef022SMatthew Knepley     Vec          coordinates;
298649ef022SMatthew Knepley     PetscScalar *coords;
299649ef022SMatthew Knepley     PetscReal    theta;
300649ef022SMatthew Knepley     PetscInt     cdim, N, bs, i;
301649ef022SMatthew Knepley 
3029566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(*dm, &cdim));
3039566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(*dm, &coordinates));
3049566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coordinates, &N));
3059566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(coordinates, &bs));
30663a3b9bcSJacob Faibussowitsch     PetscCheck(bs == cdim, comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim);
3079566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
3089566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **)&param));
309649ef022SMatthew Knepley     theta = param->theta;
310649ef022SMatthew Knepley     for (i = 0; i < N; i += cdim) {
311649ef022SMatthew Knepley       PetscScalar x = coords[i + 0];
312649ef022SMatthew Knepley       PetscScalar y = coords[i + 1];
313649ef022SMatthew Knepley 
314649ef022SMatthew Knepley       coords[i + 0] = PetscCosReal(theta) * x - PetscSinReal(theta) * y;
315649ef022SMatthew Knepley       coords[i + 1] = PetscSinReal(theta) * x + PetscCosReal(theta) * y;
316649ef022SMatthew Knepley     }
3179566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
3189566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinates(*dm, coordinates));
319649ef022SMatthew Knepley   }
3209566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
321649ef022SMatthew Knepley   PetscFunctionReturn(0);
322649ef022SMatthew Knepley }
323649ef022SMatthew Knepley 
324*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *user) {
325649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
326649ef022SMatthew Knepley   PetscDS    prob;
32745480ffeSMatthew G. Knepley   DMLabel    label;
328649ef022SMatthew Knepley   Parameter *ctx;
329649ef022SMatthew Knepley   PetscInt   id;
330649ef022SMatthew Knepley 
331649ef022SMatthew Knepley   PetscFunctionBeginUser;
3329566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
3339566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
334649ef022SMatthew Knepley   switch (user->solType) {
335649ef022SMatthew Knepley   case SOL_QUADRATIC:
3369566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v));
3379566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w));
338649ef022SMatthew Knepley 
339649ef022SMatthew Knepley     exactFuncs[0] = quadratic_u;
340649ef022SMatthew Knepley     exactFuncs[1] = linear_p;
341649ef022SMatthew Knepley     exactFuncs[2] = linear_T;
342649ef022SMatthew Knepley     break;
343649ef022SMatthew Knepley   case SOL_CUBIC:
3449566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v));
3459566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w));
346649ef022SMatthew Knepley 
347649ef022SMatthew Knepley     exactFuncs[0] = cubic_u;
348649ef022SMatthew Knepley     exactFuncs[1] = quadratic_p;
349649ef022SMatthew Knepley     exactFuncs[2] = quadratic_T;
350649ef022SMatthew Knepley     break;
35163a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
352649ef022SMatthew Knepley   }
353649ef022SMatthew Knepley 
3549566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL));
355649ef022SMatthew Knepley 
3569566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu, NULL, g3_vu));
3579566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_vp, NULL));
3589566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL, NULL));
3599566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL, NULL));
3609566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL, g3_wT));
361649ef022SMatthew Knepley   /* Setup constants */
362649ef022SMatthew Knepley   {
363649ef022SMatthew Knepley     Parameter  *param;
364649ef022SMatthew Knepley     PetscScalar constants[3];
365649ef022SMatthew Knepley 
3669566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **)&param));
367649ef022SMatthew Knepley 
368649ef022SMatthew Knepley     constants[0] = param->nu;
369649ef022SMatthew Knepley     constants[1] = param->alpha;
370649ef022SMatthew Knepley     constants[2] = param->theta;
3719566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 3, constants));
372649ef022SMatthew Knepley   }
373649ef022SMatthew Knepley   /* Setup Boundary Conditions */
3749566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&ctx));
375649ef022SMatthew Knepley   id = 3;
3769566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL));
377649ef022SMatthew Knepley   id = 1;
3789566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL));
379649ef022SMatthew Knepley   id = 2;
3809566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL));
381649ef022SMatthew Knepley   id = 4;
3829566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void))exactFuncs[0], NULL, ctx, NULL));
383649ef022SMatthew Knepley   id = 3;
3849566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL));
385649ef022SMatthew Knepley   id = 1;
3869566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL));
387649ef022SMatthew Knepley   id = 2;
3889566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL));
389649ef022SMatthew Knepley   id = 4;
3909566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void))exactFuncs[2], NULL, ctx, NULL));
391649ef022SMatthew Knepley 
392649ef022SMatthew Knepley   /*setup exact solution.*/
3939566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx));
3949566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx));
3959566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx));
396649ef022SMatthew Knepley   PetscFunctionReturn(0);
397649ef022SMatthew Knepley }
398649ef022SMatthew Knepley 
399*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user) {
400649ef022SMatthew Knepley   DM         cdm = dm;
401649ef022SMatthew Knepley   PetscFE    fe[3];
402649ef022SMatthew Knepley   Parameter *param;
403649ef022SMatthew Knepley   MPI_Comm   comm;
40430602db0SMatthew G. Knepley   PetscInt   dim;
40530602db0SMatthew G. Knepley   PetscBool  simplex;
406649ef022SMatthew Knepley 
407649ef022SMatthew Knepley   PetscFunctionBeginUser;
4089566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4099566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
410649ef022SMatthew Knepley   /* Create finite element */
4119566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
4129566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
4139566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
414649ef022SMatthew Knepley 
4159566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
4169566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
4179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
418649ef022SMatthew Knepley 
4199566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]));
4209566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
4219566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe[2], "temperature"));
422649ef022SMatthew Knepley 
423649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
4249566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
4259566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
4269566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe[2]));
4279566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
4289566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
4299566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **)&param));
430649ef022SMatthew Knepley   while (cdm) {
4319566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
4329566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0));
4339566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
434649ef022SMatthew Knepley   }
4359566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[0]));
4369566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[1]));
4379566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[2]));
438649ef022SMatthew Knepley   PetscFunctionReturn(0);
439649ef022SMatthew Knepley }
440649ef022SMatthew Knepley 
441*9371c9d4SSatish Balay static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace) {
442649ef022SMatthew Knepley   Vec vec;
443649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
444649ef022SMatthew Knepley 
445649ef022SMatthew Knepley   PetscFunctionBeginUser;
44663a3b9bcSJacob Faibussowitsch   PetscCheck(ofield == 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield);
447649ef022SMatthew Knepley   funcs[nfield] = constant;
4489566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &vec));
4499566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
4509566063dSJacob Faibussowitsch   PetscCall(VecNormalize(vec, NULL));
4519566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)vec, "Pressure Null Space"));
4529566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view"));
4539566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullSpace));
4549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec));
455649ef022SMatthew Knepley   PetscFunctionReturn(0);
456649ef022SMatthew Knepley }
457649ef022SMatthew Knepley 
458*9371c9d4SSatish Balay int main(int argc, char **argv) {
459649ef022SMatthew Knepley   SNES   snes; /* nonlinear solver */
460649ef022SMatthew Knepley   DM     dm;   /* problem definition */
461649ef022SMatthew Knepley   Vec    u, r; /* solution, residual vectors */
462649ef022SMatthew Knepley   AppCtx user; /* user-defined work context */
463649ef022SMatthew Knepley 
464327415f7SBarry Smith   PetscFunctionBeginUser;
4659566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4669566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
4679566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
4689566063dSJacob Faibussowitsch   PetscCall(SetupParameters(&user));
4699566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
4709566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
4719566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
4729566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &user));
473649ef022SMatthew Knepley   /* Setup problem */
4749566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &user));
4759566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
476649ef022SMatthew Knepley 
4779566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
4789566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Solution"));
4799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
480649ef022SMatthew Knepley 
4819566063dSJacob Faibussowitsch   PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace));
4829566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user));
483649ef022SMatthew Knepley 
4849566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
485649ef022SMatthew Knepley   {
486649ef022SMatthew Knepley     PetscDS ds;
487649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
488649ef022SMatthew Knepley     void *ctxs[3];
489649ef022SMatthew Knepley 
4909566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
4919566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
4929566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
4939566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]));
4949566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u));
4959566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)u, "Exact Solution"));
4969566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view"));
497649ef022SMatthew Knepley   }
4989566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
4999566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
5009566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
501649ef022SMatthew Knepley 
502649ef022SMatthew Knepley   if (user.showError) {
503649ef022SMatthew Knepley     PetscDS ds;
504649ef022SMatthew Knepley     Vec     r;
505649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
506649ef022SMatthew Knepley     void *ctxs[3];
507649ef022SMatthew Knepley 
5089566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
5099566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
5109566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
5119566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]));
5129566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &r));
5139566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r));
5149566063dSJacob Faibussowitsch     PetscCall(VecAXPY(r, -1.0, u));
5159566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)r, "Solution Error"));
5169566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-error_vec_view"));
5179566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &r));
518649ef022SMatthew Knepley   }
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Numerical Solution"));
5209566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
521649ef022SMatthew Knepley 
5229566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
5239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
5249566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
5259566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
5269566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
5279566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
528b122ec5aSJacob Faibussowitsch   return 0;
529649ef022SMatthew Knepley }
530649ef022SMatthew Knepley 
531649ef022SMatthew Knepley /*TEST
532649ef022SMatthew Knepley 
533649ef022SMatthew Knepley   test:
534649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1
535649ef022SMatthew Knepley     requires: triangle !single
536649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
537649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
538649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
539649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
540649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
541649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
542649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
543649ef022SMatthew Knepley 
544649ef022SMatthew Knepley   test:
545649ef022SMatthew Knepley     # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9]
546649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1_conv
547649ef022SMatthew Knepley     requires: triangle !single
548649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
549649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
550649ef022SMatthew Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \
551649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
552649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
553649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
554649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
555649ef022SMatthew Knepley 
556649ef022SMatthew Knepley   test:
557649ef022SMatthew Knepley     suffix: 2d_tri_p3_p2_p2
558649ef022SMatthew Knepley     requires: triangle !single
559649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
560649ef022SMatthew Knepley       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
561649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
562649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
563649ef022SMatthew Knepley       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
564649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
565649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
566649ef022SMatthew Knepley 
567649ef022SMatthew Knepley TEST*/
568