xref: /petsc/src/snes/tutorials/ex76.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
27649ef022SMatthew Knepley typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType;
28649ef022SMatthew Knepley const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic",  "unknown"};
29649ef022SMatthew Knepley 
30649ef022SMatthew Knepley typedef struct {
31649ef022SMatthew Knepley   PetscReal nu;      /* Kinematic viscosity */
32649ef022SMatthew Knepley   PetscReal theta;   /* Angle of pipe wall to x-axis */
33649ef022SMatthew Knepley   PetscReal alpha;   /* Thermal diffusivity */
34649ef022SMatthew Knepley   PetscReal T_in;    /* Inlet temperature*/
35649ef022SMatthew Knepley } Parameter;
36649ef022SMatthew Knepley 
37649ef022SMatthew Knepley typedef struct {
3830602db0SMatthew G. Knepley   PetscBool showError;
3930602db0SMatthew G. Knepley   PetscBag  bag;
40649ef022SMatthew Knepley   SolType   solType;
41649ef022SMatthew Knepley } AppCtx;
42649ef022SMatthew Knepley 
43649ef022SMatthew Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
44649ef022SMatthew Knepley {
45649ef022SMatthew Knepley   PetscInt d;
46649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 0.0;
47649ef022SMatthew Knepley   return 0;
48649ef022SMatthew Knepley }
49649ef022SMatthew Knepley 
50649ef022SMatthew Knepley static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51649ef022SMatthew Knepley {
52649ef022SMatthew Knepley   PetscInt d;
53649ef022SMatthew Knepley   for (d = 0; d < Nc; ++d) u[d] = 1.0;
54649ef022SMatthew Knepley   return 0;
55649ef022SMatthew Knepley }
56649ef022SMatthew Knepley 
57649ef022SMatthew Knepley /*
58649ef022SMatthew Knepley   CASE: quadratic
59649ef022SMatthew Knepley   In 2D we use exact solution:
60649ef022SMatthew Knepley 
61649ef022SMatthew Knepley     u = x^2 + y^2
62649ef022SMatthew Knepley     v = 2x^2 - 2xy
63649ef022SMatthew Knepley     p = x + y - 1
64649ef022SMatthew Knepley     T = x + y
65649ef022SMatthew Knepley     f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -4\nu + 1>
66649ef022SMatthew Knepley     Q = 3x^2 + y^2 - 2xy
67649ef022SMatthew Knepley 
68649ef022SMatthew Knepley   so that
69649ef022SMatthew Knepley 
70649ef022SMatthew Knepley (1)  \nabla \cdot u  = 2x - 2x = 0
71649ef022SMatthew Knepley 
72649ef022SMatthew Knepley (2)  u \cdot \nabla u - \nu \Delta u + \nabla p - f
73649ef022SMatthew 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
74649ef022SMatthew Knepley 
75649ef022SMatthew Knepley (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0
76649ef022SMatthew Knepley */
77649ef022SMatthew Knepley 
78649ef022SMatthew Knepley static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
79649ef022SMatthew Knepley {
80649ef022SMatthew Knepley   u[0] = X[0]*X[0] + X[1]*X[1];
81649ef022SMatthew Knepley   u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
82649ef022SMatthew Knepley   return 0;
83649ef022SMatthew Knepley }
84649ef022SMatthew Knepley 
85649ef022SMatthew Knepley static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
86649ef022SMatthew Knepley {
87649ef022SMatthew Knepley   p[0] = X[0] + X[1] - 1.0;
88649ef022SMatthew Knepley   return 0;
89649ef022SMatthew Knepley }
90649ef022SMatthew Knepley 
91649ef022SMatthew Knepley static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
92649ef022SMatthew Knepley {
93649ef022SMatthew Knepley   T[0] = X[0] + X[1];
94649ef022SMatthew Knepley   return 0;
95649ef022SMatthew Knepley }
96649ef022SMatthew Knepley 
97649ef022SMatthew Knepley static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
98649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
99649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
101649ef022SMatthew Knepley {
102649ef022SMatthew Knepley   PetscInt                   c, d;
103649ef022SMatthew Knepley   PetscInt                   Nc = dim;
104649ef022SMatthew Knepley   const PetscReal    nu = PetscRealPart(constants[0]);
105649ef022SMatthew Knepley 
106649ef022SMatthew Knepley   for (c=0; c<Nc; ++c) {
107649ef022SMatthew Knepley     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
108649ef022SMatthew Knepley   }
109649ef022SMatthew 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);
110649ef022SMatthew 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);
111649ef022SMatthew Knepley }
112649ef022SMatthew Knepley 
113649ef022SMatthew Knepley static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114649ef022SMatthew Knepley                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115649ef022SMatthew Knepley                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116649ef022SMatthew Knepley                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
117649ef022SMatthew Knepley {
118649ef022SMatthew Knepley   PetscInt d;
119649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
120649ef022SMatthew Knepley   f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]);
121649ef022SMatthew Knepley }
122649ef022SMatthew Knepley 
123649ef022SMatthew Knepley /*
124649ef022SMatthew Knepley   CASE: cubic
125649ef022SMatthew Knepley   In 2D we use exact solution:
126649ef022SMatthew Knepley 
127649ef022SMatthew Knepley     u = x^3 + y^3
128649ef022SMatthew Knepley     v = 2x^3 - 3x^2y
129649ef022SMatthew Knepley     p = 3/2 x^2 + 3/2 y^2 - 1
130649ef022SMatthew Knepley     T = 1/2 x^2 + 1/2 y^2
131649ef022SMatthew Knepley     f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y>
132649ef022SMatthew Knepley     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2
133649ef022SMatthew Knepley 
134649ef022SMatthew Knepley   so that
135649ef022SMatthew Knepley 
136649ef022SMatthew Knepley   \nabla \cdot u = 3x^2 - 3x^2 = 0
137649ef022SMatthew Knepley 
138649ef022SMatthew Knepley   u \cdot \nabla u - \nu \Delta u + \nabla p - f
139649ef022SMatthew 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
140649ef022SMatthew Knepley 
141649ef022SMatthew 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
142649ef022SMatthew Knepley */
143649ef022SMatthew Knepley 
144649ef022SMatthew Knepley static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
145649ef022SMatthew Knepley {
146649ef022SMatthew Knepley   u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
147649ef022SMatthew Knepley   u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
148649ef022SMatthew Knepley   return 0;
149649ef022SMatthew Knepley }
150649ef022SMatthew Knepley 
151649ef022SMatthew Knepley static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
152649ef022SMatthew Knepley {
153649ef022SMatthew Knepley   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
154649ef022SMatthew Knepley   return 0;
155649ef022SMatthew Knepley }
156649ef022SMatthew Knepley 
157649ef022SMatthew Knepley static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
158649ef022SMatthew Knepley {
159649ef022SMatthew Knepley   T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
160649ef022SMatthew Knepley   return 0;
161649ef022SMatthew Knepley }
162649ef022SMatthew Knepley 
163649ef022SMatthew Knepley static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
164649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
165649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
166649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
167649ef022SMatthew Knepley {
168649ef022SMatthew Knepley   PetscInt                   c, d;
169649ef022SMatthew Knepley   PetscInt                   Nc = dim;
170649ef022SMatthew Knepley   const PetscReal    nu = PetscRealPart(constants[0]);
171649ef022SMatthew Knepley 
172649ef022SMatthew Knepley   for (c=0; c<Nc; ++c) {
173649ef022SMatthew Knepley     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
174649ef022SMatthew Knepley   }
175649ef022SMatthew 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]);
176649ef022SMatthew 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]);
177649ef022SMatthew Knepley }
178649ef022SMatthew Knepley 
179649ef022SMatthew Knepley static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
180649ef022SMatthew Knepley                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
181649ef022SMatthew Knepley                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
182649ef022SMatthew Knepley                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
183649ef022SMatthew Knepley {
184649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
185649ef022SMatthew Knepley   PetscInt        d;
186649ef022SMatthew Knepley 
187649ef022SMatthew Knepley   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
188649ef022SMatthew 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);
189649ef022SMatthew Knepley }
190649ef022SMatthew Knepley 
191649ef022SMatthew Knepley static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
192649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
193649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
194649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
195649ef022SMatthew Knepley {
196649ef022SMatthew Knepley   PetscInt d;
197649ef022SMatthew Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
198649ef022SMatthew Knepley }
199649ef022SMatthew Knepley 
200649ef022SMatthew Knepley static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
201649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
202649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
203649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
204649ef022SMatthew Knepley {
205649ef022SMatthew Knepley   const PetscReal nu = PetscRealPart(constants[0]);
206649ef022SMatthew Knepley   const PetscInt  Nc = dim;
207649ef022SMatthew Knepley   PetscInt        c, d;
208649ef022SMatthew Knepley 
209649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
210649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
211649ef022SMatthew Knepley       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
212649ef022SMatthew Knepley       //f1[c*dim+d] = nu*u_x[c*dim+d];
213649ef022SMatthew Knepley     }
214649ef022SMatthew Knepley     f1[c*dim+c] -= u[uOff[1]];
215649ef022SMatthew Knepley   }
216649ef022SMatthew Knepley }
217649ef022SMatthew Knepley 
218649ef022SMatthew Knepley static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
219649ef022SMatthew Knepley                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
220649ef022SMatthew Knepley                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
221649ef022SMatthew Knepley                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
222649ef022SMatthew Knepley {
223649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
224649ef022SMatthew Knepley   PetscInt d;
225649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
226649ef022SMatthew Knepley }
227649ef022SMatthew Knepley 
228649ef022SMatthew Knepley /*Jacobians*/
229649ef022SMatthew Knepley static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
230649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
231649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
232649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
233649ef022SMatthew Knepley {
234649ef022SMatthew Knepley   PetscInt d;
235649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
236649ef022SMatthew Knepley }
237649ef022SMatthew Knepley 
238649ef022SMatthew Knepley static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
239649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
240649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
241649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
242649ef022SMatthew Knepley {
243649ef022SMatthew Knepley   const PetscInt  Nc = dim;
244649ef022SMatthew Knepley    PetscInt            c, d;
245649ef022SMatthew Knepley 
246649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
247649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
248649ef022SMatthew Knepley       g0[c*Nc+d] = u_x[ c*Nc+d];
249649ef022SMatthew Knepley     }
250649ef022SMatthew Knepley   }
251649ef022SMatthew Knepley }
252649ef022SMatthew Knepley 
253649ef022SMatthew Knepley static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
254649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
255649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
256649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
257649ef022SMatthew Knepley {
258649ef022SMatthew Knepley   PetscInt NcI = dim;
259649ef022SMatthew Knepley   PetscInt NcJ = dim;
260649ef022SMatthew Knepley   PetscInt c, d, e;
261649ef022SMatthew Knepley 
262649ef022SMatthew Knepley   for (c = 0; c < NcI; ++c) {
263649ef022SMatthew Knepley     for (d = 0; d < NcJ; ++d) {
264649ef022SMatthew Knepley       for (e = 0; e < dim; ++e) {
265649ef022SMatthew Knepley         if (c == d) {
266649ef022SMatthew Knepley           g1[(c*NcJ+d)*dim+e] = u[e];
267649ef022SMatthew Knepley         }
268649ef022SMatthew Knepley       }
269649ef022SMatthew Knepley     }
270649ef022SMatthew Knepley   }
271649ef022SMatthew Knepley }
272649ef022SMatthew Knepley 
273649ef022SMatthew Knepley static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
274649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
275649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
276649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
277649ef022SMatthew Knepley {
278649ef022SMatthew Knepley   PetscInt d;
279649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
280649ef022SMatthew Knepley }
281649ef022SMatthew Knepley 
282649ef022SMatthew Knepley static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
283649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
284649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
285649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
286649ef022SMatthew Knepley {
287649ef022SMatthew Knepley    const PetscReal      nu = PetscRealPart(constants[0]);
288649ef022SMatthew Knepley    const PetscInt         Nc = dim;
289649ef022SMatthew Knepley    PetscInt                     c, d;
290649ef022SMatthew Knepley 
291649ef022SMatthew Knepley   for (c = 0; c < Nc; ++c) {
292649ef022SMatthew Knepley     for (d = 0; d < dim; ++d) {
293649ef022SMatthew Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU
294649ef022SMatthew Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose
295649ef022SMatthew Knepley     }
296649ef022SMatthew Knepley   }
297649ef022SMatthew Knepley }
298649ef022SMatthew Knepley 
299649ef022SMatthew Knepley static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
300649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
301649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
302649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
303649ef022SMatthew Knepley {
304649ef022SMatthew Knepley   PetscInt d;
305649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
306649ef022SMatthew Knepley }
307649ef022SMatthew Knepley 
308649ef022SMatthew Knepley static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
309649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
310649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
311649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
312649ef022SMatthew Knepley {
313649ef022SMatthew Knepley   PetscInt d;
314649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
315649ef022SMatthew Knepley }
316649ef022SMatthew Knepley 
317649ef022SMatthew Knepley static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
318649ef022SMatthew Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
319649ef022SMatthew Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
320649ef022SMatthew Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
321649ef022SMatthew Knepley {
322649ef022SMatthew Knepley   const PetscReal alpha = PetscRealPart(constants[1]);
323649ef022SMatthew Knepley   PetscInt        d;
324649ef022SMatthew Knepley 
325649ef022SMatthew Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
326649ef022SMatthew Knepley }
327649ef022SMatthew Knepley 
328649ef022SMatthew Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
329649ef022SMatthew Knepley {
33030602db0SMatthew G. Knepley   PetscInt       sol;
331649ef022SMatthew Knepley 
332649ef022SMatthew Knepley   PetscFunctionBeginUser;
333649ef022SMatthew Knepley   options->solType   = SOL_QUADRATIC;
334649ef022SMatthew Knepley   options->showError = PETSC_FALSE;
335d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");
336649ef022SMatthew Knepley   sol = options->solType;
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL));
338649ef022SMatthew Knepley   options->solType = (SolType) sol;
3399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL));
340d0609cedSBarry Smith   PetscOptionsEnd();
341649ef022SMatthew Knepley   PetscFunctionReturn(0);
342649ef022SMatthew Knepley }
343649ef022SMatthew Knepley 
344649ef022SMatthew Knepley static PetscErrorCode SetupParameters(AppCtx *user)
345649ef022SMatthew Knepley {
346649ef022SMatthew Knepley   PetscBag       bag;
347649ef022SMatthew Knepley   Parameter     *p;
348649ef022SMatthew Knepley 
349649ef022SMatthew Knepley   PetscFunctionBeginUser;
350649ef022SMatthew Knepley   /* setup PETSc parameter bag */
3519566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **) &p));
3529566063dSJacob Faibussowitsch   PetscCall(PetscBagSetName(user->bag, "par", "Poiseuille flow parameters"));
353649ef022SMatthew Knepley   bag  = user->bag;
3549566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->nu,    1.0,   "nu",      "Kinematic viscosity"));
3559566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->alpha, 1.0,   "alpha",   "Thermal diffusivity"));
3569566063dSJacob Faibussowitsch   PetscCall(PetscBagRegisterReal(bag, &p->theta, 0.0,   "theta",   "Angle of pipe wall to x-axis"));
357649ef022SMatthew Knepley   PetscFunctionReturn(0);
358649ef022SMatthew Knepley }
359649ef022SMatthew Knepley 
360649ef022SMatthew Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
361649ef022SMatthew Knepley {
362649ef022SMatthew Knepley   PetscFunctionBeginUser;
3639566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
3649566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
3659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
366649ef022SMatthew Knepley   {
367649ef022SMatthew Knepley     Parameter   *param;
368649ef022SMatthew Knepley     Vec          coordinates;
369649ef022SMatthew Knepley     PetscScalar *coords;
370649ef022SMatthew Knepley     PetscReal    theta;
371649ef022SMatthew Knepley     PetscInt     cdim, N, bs, i;
372649ef022SMatthew Knepley 
3739566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(*dm, &cdim));
3749566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(*dm, &coordinates));
3759566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(coordinates, &N));
3769566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(coordinates, &bs));
37763a3b9bcSJacob Faibussowitsch     PetscCheck(bs == cdim,comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %" PetscInt_FMT " != embedding dimension %" PetscInt_FMT, bs, cdim);
3789566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &coords));
3799566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **) &param));
380649ef022SMatthew Knepley     theta = param->theta;
381649ef022SMatthew Knepley     for (i = 0; i < N; i += cdim) {
382649ef022SMatthew Knepley       PetscScalar x = coords[i+0];
383649ef022SMatthew Knepley       PetscScalar y = coords[i+1];
384649ef022SMatthew Knepley 
385649ef022SMatthew Knepley       coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y;
386649ef022SMatthew Knepley       coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y;
387649ef022SMatthew Knepley     }
3889566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &coords));
3899566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinates(*dm, coordinates));
390649ef022SMatthew Knepley   }
3919566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
392649ef022SMatthew Knepley   PetscFunctionReturn(0);
393649ef022SMatthew Knepley }
394649ef022SMatthew Knepley 
395649ef022SMatthew Knepley static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
396649ef022SMatthew Knepley {
397649ef022SMatthew Knepley   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
398649ef022SMatthew Knepley   PetscDS          prob;
39945480ffeSMatthew G. Knepley   DMLabel          label;
400649ef022SMatthew Knepley   Parameter       *ctx;
401649ef022SMatthew Knepley   PetscInt         id;
402649ef022SMatthew Knepley 
403649ef022SMatthew Knepley   PetscFunctionBeginUser;
4049566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
4059566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
406649ef022SMatthew Knepley   switch(user->solType) {
407649ef022SMatthew Knepley   case SOL_QUADRATIC:
4089566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v));
4099566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w));
410649ef022SMatthew Knepley 
411649ef022SMatthew Knepley     exactFuncs[0] = quadratic_u;
412649ef022SMatthew Knepley     exactFuncs[1] = linear_p;
413649ef022SMatthew Knepley     exactFuncs[2] = linear_T;
414649ef022SMatthew Knepley     break;
415649ef022SMatthew Knepley   case SOL_CUBIC:
4169566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v));
4179566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w));
418649ef022SMatthew Knepley 
419649ef022SMatthew Knepley     exactFuncs[0] = cubic_u;
420649ef022SMatthew Knepley     exactFuncs[1] = quadratic_p;
421649ef022SMatthew Knepley     exactFuncs[2] = quadratic_T;
422649ef022SMatthew Knepley     break;
42363a3b9bcSJacob Faibussowitsch    default: SETERRQ(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%d)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
424649ef022SMatthew Knepley   }
425649ef022SMatthew Knepley 
4269566063dSJacob Faibussowitsch   PetscCall(PetscDSSetResidual(prob, 1, f0_q, NULL));
427649ef022SMatthew Knepley 
4289566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu));
4299566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL));
4309566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL));
4319566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL));
4329566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL,  g3_wT));
433649ef022SMatthew Knepley   /* Setup constants */
434649ef022SMatthew Knepley   {
435649ef022SMatthew Knepley     Parameter  *param;
436649ef022SMatthew Knepley     PetscScalar constants[3];
437649ef022SMatthew Knepley 
4389566063dSJacob Faibussowitsch     PetscCall(PetscBagGetData(user->bag, (void **) &param));
439649ef022SMatthew Knepley 
440649ef022SMatthew Knepley     constants[0] = param->nu;
441649ef022SMatthew Knepley     constants[1] = param->alpha;
442649ef022SMatthew Knepley     constants[2] = param->theta;
4439566063dSJacob Faibussowitsch     PetscCall(PetscDSSetConstants(prob, 3, constants));
444649ef022SMatthew Knepley   }
445649ef022SMatthew Knepley   /* Setup Boundary Conditions */
4469566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **) &ctx));
447649ef022SMatthew Knepley   id   = 3;
4489566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity",    label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL));
449649ef022SMatthew Knepley   id   = 1;
4509566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL));
451649ef022SMatthew Knepley   id   = 2;
4529566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity",  label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL));
453649ef022SMatthew Knepley   id   = 4;
4549566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity",   label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL));
455649ef022SMatthew Knepley   id   = 3;
4569566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp",    label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL));
457649ef022SMatthew Knepley   id   = 1;
4589566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL));
459649ef022SMatthew Knepley   id   = 2;
4609566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp",  label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL));
461649ef022SMatthew Knepley   id   = 4;
4629566063dSJacob Faibussowitsch   PetscCall(PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp",   label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL));
463649ef022SMatthew Knepley 
464649ef022SMatthew Knepley   /*setup exact solution.*/
4659566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx));
4669566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx));
4679566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx));
468649ef022SMatthew Knepley   PetscFunctionReturn(0);
469649ef022SMatthew Knepley }
470649ef022SMatthew Knepley 
471649ef022SMatthew Knepley static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
472649ef022SMatthew Knepley {
473649ef022SMatthew Knepley   DM              cdm   = dm;
474649ef022SMatthew Knepley   PetscFE         fe[3];
475649ef022SMatthew Knepley   Parameter      *param;
476649ef022SMatthew Knepley   MPI_Comm        comm;
47730602db0SMatthew G. Knepley   PetscInt        dim;
47830602db0SMatthew G. Knepley   PetscBool       simplex;
479649ef022SMatthew Knepley 
480649ef022SMatthew Knepley   PetscFunctionBeginUser;
4819566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
4829566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
483649ef022SMatthew Knepley   /* Create finite element */
4849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
4859566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
4869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe[0], "velocity"));
487649ef022SMatthew Knepley 
4889566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
4899566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
4909566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe[1], "pressure"));
491649ef022SMatthew Knepley 
4929566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]));
4939566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe[0], fe[2]));
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe[2], "temperature"));
495649ef022SMatthew Knepley 
496649ef022SMatthew Knepley   /* Set discretization and boundary conditions for each mesh */
4979566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
4989566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
4999566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 2, NULL, (PetscObject) fe[2]));
5009566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
5019566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, user));
5029566063dSJacob Faibussowitsch   PetscCall(PetscBagGetData(user->bag, (void **) &param));
503649ef022SMatthew Knepley   while (cdm) {
5049566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
5059566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0));
5069566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
507649ef022SMatthew Knepley   }
5089566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[0]));
5099566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[1]));
5109566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe[2]));
511649ef022SMatthew Knepley   PetscFunctionReturn(0);
512649ef022SMatthew Knepley }
513649ef022SMatthew Knepley 
514649ef022SMatthew Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
515649ef022SMatthew Knepley {
516649ef022SMatthew Knepley   Vec              vec;
517649ef022SMatthew Knepley   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
518649ef022SMatthew Knepley 
519649ef022SMatthew Knepley   PetscFunctionBeginUser;
52063a3b9bcSJacob Faibussowitsch   PetscCheck(ofield == 1,PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %" PetscInt_FMT, ofield);
521649ef022SMatthew Knepley   funcs[nfield] = constant;
5229566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &vec));
5239566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec));
5249566063dSJacob Faibussowitsch   PetscCall(VecNormalize(vec, NULL));
5259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) vec, "Pressure Null Space"));
5269566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(vec, NULL, "-pressure_nullspace_view"));
5279566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace));
5289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vec));
529649ef022SMatthew Knepley   PetscFunctionReturn(0);
530649ef022SMatthew Knepley }
531649ef022SMatthew Knepley 
532649ef022SMatthew Knepley int main(int argc, char **argv)
533649ef022SMatthew Knepley {
534649ef022SMatthew Knepley   SNES            snes;                 /* nonlinear solver */
535649ef022SMatthew Knepley   DM              dm;                   /* problem definition */
536649ef022SMatthew Knepley   Vec             u, r;                 /* solution, residual vectors */
537649ef022SMatthew Knepley   AppCtx          user;                 /* user-defined work context */
538649ef022SMatthew Knepley 
539*327415f7SBarry Smith   PetscFunctionBeginUser;
5409566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
5419566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
5429566063dSJacob Faibussowitsch   PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag));
5439566063dSJacob Faibussowitsch   PetscCall(SetupParameters(&user));
5449566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
5459566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
5469566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, dm));
5479566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &user));
548649ef022SMatthew Knepley   /* Setup problem */
5499566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &user));
5509566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateClosureIndex(dm, NULL));
551649ef022SMatthew Knepley 
5529566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
5539566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "Solution"));
5549566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
555649ef022SMatthew Knepley 
5569566063dSJacob Faibussowitsch   PetscCall(DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace));
5579566063dSJacob Faibussowitsch   PetscCall(DMPlexSetSNESLocalFEM(dm,&user,&user,&user));
558649ef022SMatthew Knepley 
5599566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
560649ef022SMatthew Knepley   {
561649ef022SMatthew Knepley     PetscDS          ds;
562649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
563649ef022SMatthew Knepley     void            *ctxs[3];
564649ef022SMatthew Knepley 
5659566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
5669566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
5679566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
5689566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]));
5699566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u));
5709566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) u, "Exact Solution"));
5719566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(u, NULL, "-exact_vec_view"));
572649ef022SMatthew Knepley   }
5739566063dSJacob Faibussowitsch   PetscCall(DMSNESCheckFromOptions(snes, u));
5749566063dSJacob Faibussowitsch   PetscCall(VecSet(u, 0.0));
5759566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, u));
576649ef022SMatthew Knepley 
577649ef022SMatthew Knepley   if (user.showError) {
578649ef022SMatthew Knepley     PetscDS          ds;
579649ef022SMatthew Knepley     Vec              r;
580649ef022SMatthew Knepley     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
581649ef022SMatthew Knepley     void            *ctxs[3];
582649ef022SMatthew Knepley 
5839566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
5849566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]));
5859566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]));
5869566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]));
5879566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &r));
5889566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r));
5899566063dSJacob Faibussowitsch     PetscCall(VecAXPY(r, -1.0, u));
5909566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) r, "Solution Error"));
5919566063dSJacob Faibussowitsch     PetscCall(VecViewFromOptions(r, NULL, "-error_vec_view"));
5929566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &r));
593649ef022SMatthew Knepley   }
5949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "Numerical Solution"));
5959566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
596649ef022SMatthew Knepley 
5979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
5989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
5999566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
6009566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
6019566063dSJacob Faibussowitsch   PetscCall(PetscBagDestroy(&user.bag));
6029566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
603b122ec5aSJacob Faibussowitsch   return 0;
604649ef022SMatthew Knepley }
605649ef022SMatthew Knepley 
606649ef022SMatthew Knepley /*TEST
607649ef022SMatthew Knepley 
608649ef022SMatthew Knepley   test:
609649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1
610649ef022SMatthew Knepley     requires: triangle !single
611649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
612649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
613649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
614649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
615649ef022SMatthew 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 \
616649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
617649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
618649ef022SMatthew Knepley 
619649ef022SMatthew Knepley   test:
620649ef022SMatthew Knepley     # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9]
621649ef022SMatthew Knepley     suffix: 2d_tri_p2_p1_p1_conv
622649ef022SMatthew Knepley     requires: triangle !single
623649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
624649ef022SMatthew Knepley       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
625649ef022SMatthew Knepley       -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \
626649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
627649ef022SMatthew 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 \
628649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
629649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
630649ef022SMatthew Knepley 
631649ef022SMatthew Knepley   test:
632649ef022SMatthew Knepley     suffix: 2d_tri_p3_p2_p2
633649ef022SMatthew Knepley     requires: triangle !single
634649ef022SMatthew Knepley     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
635649ef022SMatthew Knepley       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
636649ef022SMatthew Knepley       -dmsnes_check .001 -snes_error_if_not_converged \
637649ef022SMatthew Knepley       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
638649ef022SMatthew 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 \
639649ef022SMatthew Knepley         -fieldsplit_0_pc_type lu \
640649ef022SMatthew Knepley         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
641649ef022SMatthew Knepley 
642649ef022SMatthew Knepley TEST*/
643