xref: /petsc/src/snes/tutorials/ex76.c (revision 60a1b3afd779eccf5fbace1d4ebce5aca36800a9)
1 static char help[] = "Low Mach Flow in 2d and 3d channels with finite elements.\n\
2 We solve the Low Mach flow problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4 
5 /*F
6 This Low Mach flow is a steady-state isoviscous Navier-Stokes flow. We discretize using the
7 finite element method on an unstructured mesh. The weak form equations are
8 
9 \begin{align*}
10     < q, \nabla\cdot u >                                                                                     = 0
11     <v, u \cdot \nabla u> + < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p >  - < v, f  >  = 0
12     < w, u \cdot \nabla T > - < \nabla w, \alpha \nabla T > - < w, Q >                                       = 0
13 \end{align*}
14 
15 where $\nu$ is the kinematic viscosity and $\alpha$ is thermal diffusivity.
16 
17 For visualization, use
18 
19   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
20 F*/
21 
22 #include <petscdmplex.h>
23 #include <petscsnes.h>
24 #include <petscds.h>
25 #include <petscbag.h>
26 
27 typedef enum {SOL_QUADRATIC, SOL_CUBIC, NUM_SOL_TYPES} SolType;
28 const char *solTypes[NUM_SOL_TYPES+1] = {"quadratic", "cubic",  "unknown"};
29 
30 typedef struct {
31   PetscReal nu;      /* Kinematic viscosity */
32   PetscReal theta;   /* Angle of pipe wall to x-axis */
33   PetscReal alpha;   /* Thermal diffusivity */
34   PetscReal T_in;    /* Inlet temperature*/
35 } Parameter;
36 
37 typedef struct {
38   PetscBool showError;
39   PetscBag  bag;
40   SolType   solType;
41 } AppCtx;
42 
43 static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
44 {
45   PetscInt d;
46   for (d = 0; d < Nc; ++d) u[d] = 0.0;
47   return 0;
48 }
49 
50 static PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51 {
52   PetscInt d;
53   for (d = 0; d < Nc; ++d) u[d] = 1.0;
54   return 0;
55 }
56 
57 /*
58   CASE: quadratic
59   In 2D we use exact solution:
60 
61     u = x^2 + y^2
62     v = 2x^2 - 2xy
63     p = x + y - 1
64     T = x + y
65     f = <2x^3 + 4x^2y - 2xy^2 -4\nu + 1,  4xy^2 + 2x^2y - 2y^3 -4\nu + 1>
66     Q = 3x^2 + y^2 - 2xy
67 
68   so that
69 
70 (1)  \nabla \cdot u  = 2x - 2x = 0
71 
72 (2)  u \cdot \nabla u - \nu \Delta u + \nabla p - f
73      = <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
74 
75 (3) u \cdot \nabla T - \alpha \Delta T - Q = 3x^2 + y^2 - 2xy - \alpha*0 - 3x^2 - y^2 + 2xy = 0
76 */
77 
78 static PetscErrorCode quadratic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
79 {
80   u[0] = X[0]*X[0] + X[1]*X[1];
81   u[1] = 2.0*X[0]*X[0] - 2.0*X[0]*X[1];
82   return 0;
83 }
84 
85 static PetscErrorCode linear_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
86 {
87   p[0] = X[0] + X[1] - 1.0;
88   return 0;
89 }
90 
91 static PetscErrorCode linear_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
92 {
93   T[0] = X[0] + X[1];
94   return 0;
95 }
96 
97 static void f0_quadratic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
98                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
99                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
101 {
102   PetscInt                   c, d;
103   PetscInt                   Nc = dim;
104   const PetscReal    nu = PetscRealPart(constants[0]);
105 
106   for (c=0; c<Nc; ++c) {
107     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
108   }
109   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);
110   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);
111 }
112 
113 static void f0_quadratic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116                            PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
117 {
118   PetscInt d;
119   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
120   f0[0] -= (3*X[0]*X[0] + X[1]*X[1] - 2*X[0]*X[1]);
121 }
122 
123 
124 /*
125   CASE: cubic
126   In 2D we use exact solution:
127 
128     u = x^3 + y^3
129     v = 2x^3 - 3x^2y
130     p = 3/2 x^2 + 3/2 y^2 - 1
131     T = 1/2 x^2 + 1/2 y^2
132     f = <3x^5 + 6x^3y^2 - 6x^2y^3 - \nu(6x + 6y), 6x^2y^3 + 3x^4y - 6xy^4 - \nu(12x - 6y) + 3y>
133     Q = x^4 + xy^3 + 2x^3y - 3x^2y^2 - 2
134 
135   so that
136 
137   \nabla \cdot u = 3x^2 - 3x^2 = 0
138 
139   u \cdot \nabla u - \nu \Delta u + \nabla p - f
140   = <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
141 
142   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
143 */
144 
145 static PetscErrorCode cubic_u(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
146 {
147   u[0] = X[0]*X[0]*X[0] + X[1]*X[1]*X[1];
148   u[1] = 2.0*X[0]*X[0]*X[0] - 3.0*X[0]*X[0]*X[1];
149   return 0;
150 }
151 
152 static PetscErrorCode quadratic_p(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
153 {
154   p[0] = 3.0*X[0]*X[0]/2.0 + 3.0*X[1]*X[1]/2.0 - 1.0;
155   return 0;
156 }
157 
158 static PetscErrorCode quadratic_T(PetscInt Dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *T, void *ctx)
159 {
160   T[0] = X[0]*X[0]/2.0 + X[1]*X[1]/2.0;
161   return 0;
162 }
163 
164 
165 static void f0_cubic_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
166                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
167                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
168                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
169 {
170   PetscInt                   c, d;
171   PetscInt                   Nc = dim;
172   const PetscReal    nu = PetscRealPart(constants[0]);
173 
174   for (c=0; c<Nc; ++c) {
175     for (d=0; d<dim; ++d) f0[c] += u[d]*u_x[c*dim+d];
176   }
177   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]);
178   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]);
179 }
180 
181 static void f0_cubic_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
182                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
183                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
184                        PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
185 {
186   const PetscReal alpha = PetscRealPart(constants[1]);
187   PetscInt        d;
188 
189   for (d = 0, f0[0] = 0; d < dim; ++d) f0[0] += u[uOff[0]+d]*u_x[uOff_x[2]+d];
190   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);
191 }
192 
193 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
194                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
195                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
196                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
197 {
198   PetscInt d;
199   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
200 }
201 
202 static void f1_v(PetscInt dim, PetscInt Nf, PetscInt NfAux,
203                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
204                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
205                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
206 {
207   const PetscReal nu = PetscRealPart(constants[0]);
208   const PetscInt  Nc = dim;
209   PetscInt        c, d;
210 
211   for (c = 0; c < Nc; ++c) {
212     for (d = 0; d < dim; ++d) {
213       f1[c*dim+d] = nu*(u_x[c*dim+d] + u_x[d*dim+c]);
214       //f1[c*dim+d] = nu*u_x[c*dim+d];
215     }
216     f1[c*dim+c] -= u[uOff[1]];
217   }
218 }
219 
220 static void f1_w(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223                  PetscReal t, const PetscReal X[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
224 {
225   const PetscReal alpha = PetscRealPart(constants[1]);
226   PetscInt d;
227   for (d = 0; d < dim; ++d) f1[d] = alpha*u_x[uOff_x[2]+d];
228 }
229 
230 /*Jacobians*/
231 static void g1_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
232                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
233                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
234                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
235 {
236   PetscInt d;
237   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0;
238 }
239 
240 static void g0_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
241                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
242                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
243                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
244 {
245   const PetscInt  Nc = dim;
246    PetscInt            c, d;
247 
248   for (c = 0; c < Nc; ++c) {
249     for (d = 0; d < dim; ++d) {
250       g0[c*Nc+d] = u_x[ c*Nc+d];
251     }
252   }
253 }
254 
255 static void g1_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
256                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
257                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
258                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
259 {
260   PetscInt NcI = dim;
261   PetscInt NcJ = dim;
262   PetscInt c, d, e;
263 
264   for (c = 0; c < NcI; ++c) {
265     for (d = 0; d < NcJ; ++d) {
266       for (e = 0; e < dim; ++e) {
267         if (c == d) {
268           g1[(c*NcJ+d)*dim+e] = u[e];
269         }
270       }
271     }
272   }
273 }
274 
275 static void g2_vp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
276                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
277                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
278                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
279 {
280   PetscInt d;
281   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0;
282 }
283 
284 static void g3_vu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
285                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
286                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
287                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
288 {
289    const PetscReal      nu = PetscRealPart(constants[0]);
290    const PetscInt         Nc = dim;
291    PetscInt                     c, d;
292 
293   for (c = 0; c < Nc; ++c) {
294     for (d = 0; d < dim; ++d) {
295       g3[((c*Nc+c)*dim+d)*dim+d] += nu; // gradU
296       g3[((c*Nc+d)*dim+d)*dim+c] += nu; // gradU transpose
297     }
298   }
299 }
300 
301 static void g0_wu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
302                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
303                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
304                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
305 {
306   PetscInt d;
307   for (d = 0; d < dim; ++d) g0[d] = u_x[uOff_x[2]+d];
308 }
309 
310 
311 static void g1_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
312                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
313                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
314                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
315 {
316   PetscInt d;
317   for (d = 0; d < dim; ++d) g1[d] = u[uOff[0]+d];
318 }
319 
320 static void g3_wT(PetscInt dim, PetscInt Nf, PetscInt NfAux,
321                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
322                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
323                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
324 {
325   const PetscReal alpha = PetscRealPart(constants[1]);
326   PetscInt        d;
327 
328   for (d = 0; d < dim; ++d) g3[d*dim+d] = alpha;
329 }
330 
331 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
332 {
333   PetscInt       sol;
334   PetscErrorCode ierr;
335 
336 
337   PetscFunctionBeginUser;
338   options->solType   = SOL_QUADRATIC;
339   options->showError = PETSC_FALSE;
340 
341   ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr);
342   sol = options->solType;
343   ierr = PetscOptionsEList("-sol_type", "The solution type", "ex62.c", solTypes, NUM_SOL_TYPES, solTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
344   options->solType = (SolType) sol;
345   ierr = PetscOptionsBool("-show_error", "Output the error for verification", "ex62.c", options->showError, &options->showError, NULL);CHKERRQ(ierr);
346   ierr = PetscOptionsEnd();
347   PetscFunctionReturn(0);
348 }
349 
350 static PetscErrorCode SetupParameters(AppCtx *user)
351 {
352   PetscBag       bag;
353   Parameter     *p;
354   PetscErrorCode ierr;
355 
356   PetscFunctionBeginUser;
357   /* setup PETSc parameter bag */
358   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
359   ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr);
360   bag  = user->bag;
361   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0,   "nu",      "Kinematic viscosity");CHKERRQ(ierr);
362   ierr = PetscBagRegisterReal(bag, &p->alpha, 1.0,   "alpha",   "Thermal diffusivity");CHKERRQ(ierr);
363   ierr = PetscBagRegisterReal(bag, &p->theta, 0.0,   "theta",   "Angle of pipe wall to x-axis");CHKERRQ(ierr);
364   PetscFunctionReturn(0);
365 }
366 
367 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
368 {
369   PetscErrorCode ierr;
370 
371   PetscFunctionBeginUser;
372   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
373   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
374   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
375   {
376     Parameter   *param;
377     Vec          coordinates;
378     PetscScalar *coords;
379     PetscReal    theta;
380     PetscInt     cdim, N, bs, i;
381 
382     ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr);
383     ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr);
384     ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
385     ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr);
386     if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim);
387     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
388     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
389     theta = param->theta;
390     for (i = 0; i < N; i += cdim) {
391       PetscScalar x = coords[i+0];
392       PetscScalar y = coords[i+1];
393 
394       coords[i+0] = PetscCosReal(theta)*x - PetscSinReal(theta)*y;
395       coords[i+1] = PetscSinReal(theta)*x + PetscCosReal(theta)*y;
396     }
397     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
398     ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr);
399   }
400   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
401   PetscFunctionReturn(0);
402 }
403 
404 static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
405 {
406   PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
407   PetscDS          prob;
408   DMLabel          label;
409   Parameter       *ctx;
410   PetscInt         id;
411   PetscErrorCode   ierr;
412 
413   PetscFunctionBeginUser;
414   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
415   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
416   switch(user->solType){
417   case SOL_QUADRATIC:
418     ierr = PetscDSSetResidual(prob, 0, f0_quadratic_v, f1_v);CHKERRQ(ierr);
419     ierr = PetscDSSetResidual(prob, 2, f0_quadratic_w, f1_w);CHKERRQ(ierr);
420 
421     exactFuncs[0] = quadratic_u;
422     exactFuncs[1] = linear_p;
423     exactFuncs[2] = linear_T;
424     break;
425   case SOL_CUBIC:
426     ierr = PetscDSSetResidual(prob, 0, f0_cubic_v, f1_v);CHKERRQ(ierr);
427     ierr = PetscDSSetResidual(prob, 2, f0_cubic_w, f1_w);CHKERRQ(ierr);
428 
429     exactFuncs[0] = cubic_u;
430     exactFuncs[1] = quadratic_p;
431     exactFuncs[2] = quadratic_T;
432     break;
433    default: SETERRQ2(PetscObjectComm((PetscObject) prob), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", solTypes[PetscMin(user->solType, NUM_SOL_TYPES)], user->solType);
434   }
435 
436   ierr = PetscDSSetResidual(prob, 1, f0_q, NULL);CHKERRQ(ierr);
437 
438   ierr = PetscDSSetJacobian(prob, 0, 0, g0_vu, g1_vu,  NULL,  g3_vu);CHKERRQ(ierr);
439   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_vp, NULL);CHKERRQ(ierr);
440   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_qu, NULL,  NULL);CHKERRQ(ierr);
441   ierr = PetscDSSetJacobian(prob, 2, 0, g0_wu, NULL, NULL,  NULL);CHKERRQ(ierr);
442   ierr = PetscDSSetJacobian(prob, 2, 2, NULL, g1_wT, NULL,  g3_wT);CHKERRQ(ierr);
443   /* Setup constants */
444   {
445     Parameter  *param;
446     PetscScalar constants[3];
447 
448     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
449 
450     constants[0] = param->nu;
451     constants[1] = param->alpha;
452     constants[2] = param->theta;
453     ierr = PetscDSSetConstants(prob, 3, constants);CHKERRQ(ierr);
454   }
455   /* Setup Boundary Conditions */
456   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
457   id   = 3;
458   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall velocity",    label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
459   id   = 1;
460   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall velocity", label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
461   id   = 2;
462   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall velocity",  label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
463   id   = 4;
464   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall velocity",   label, 1, &id, 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
465   id   = 3;
466   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall temp",    label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
467   id   = 1;
468   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall temp", label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
469   id   = 2;
470   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "right wall temp",  label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
471   id   = 4;
472   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "left wall temp",   label, 1, &id, 2, 0, NULL, (void (*)(void)) exactFuncs[2], NULL, ctx, NULL);CHKERRQ(ierr);
473 
474   /*setup exact solution.*/
475   ierr = PetscDSSetExactSolution(prob, 0, exactFuncs[0], ctx);CHKERRQ(ierr);
476   ierr = PetscDSSetExactSolution(prob, 1, exactFuncs[1], ctx);CHKERRQ(ierr);
477   ierr = PetscDSSetExactSolution(prob, 2, exactFuncs[2], ctx);CHKERRQ(ierr);
478   PetscFunctionReturn(0);
479 }
480 
481 static PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
482 {
483   DM              cdm   = dm;
484   PetscFE         fe[3];
485   Parameter      *param;
486   MPI_Comm        comm;
487   PetscInt        dim;
488   PetscBool       simplex;
489   PetscErrorCode  ierr;
490 
491   PetscFunctionBeginUser;
492   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
493   ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr);
494   /* Create finite element */
495   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
496   ierr = PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
497   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
498 
499   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
500   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
501   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
502 
503   ierr = PetscFECreateDefault(comm, dim, 1, simplex, "temp_", PETSC_DEFAULT, &fe[2]);CHKERRQ(ierr);
504   ierr = PetscFECopyQuadrature(fe[0], fe[2]);CHKERRQ(ierr);
505   ierr = PetscObjectSetName((PetscObject) fe[2], "temperature");CHKERRQ(ierr);
506 
507   /* Set discretization and boundary conditions for each mesh */
508   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
509   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
510   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
511   ierr = DMCreateDS(dm);CHKERRQ(ierr);
512   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
513   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
514   while (cdm) {
515     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
516     ierr = DMPlexCreateBasisRotation(cdm, param->theta, 0.0, 0.0);CHKERRQ(ierr);
517     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
518   }
519   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
520   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
521   ierr = PetscFEDestroy(&fe[2]);CHKERRQ(ierr);
522   PetscFunctionReturn(0);
523 }
524 
525 static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt ofield, PetscInt nfield, MatNullSpace *nullSpace)
526 {
527   Vec              vec;
528   PetscErrorCode (*funcs[3])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {zero, zero, zero};
529   PetscErrorCode   ierr;
530 
531   PetscFunctionBeginUser;
532   if (ofield != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Nullspace must be for pressure field at index 1, not %D", ofield);
533   funcs[nfield] = constant;
534   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
535   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
536   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
537   ierr = PetscObjectSetName((PetscObject) vec, "Pressure Null Space");CHKERRQ(ierr);
538   ierr = VecViewFromOptions(vec, NULL, "-pressure_nullspace_view");CHKERRQ(ierr);
539   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject) dm), PETSC_FALSE, 1, &vec, nullSpace);CHKERRQ(ierr);
540   ierr = VecDestroy(&vec);CHKERRQ(ierr);
541   PetscFunctionReturn(0);
542 }
543 
544 int main(int argc, char **argv)
545 {
546   SNES            snes;                 /* nonlinear solver */
547   DM              dm;                   /* problem definition */
548   Vec             u, r;                 /* solution, residual vectors */
549   AppCtx          user;                 /* user-defined work context */
550   PetscErrorCode  ierr;
551 
552   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
553   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
554   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
555   ierr = SetupParameters(&user);CHKERRQ(ierr);
556   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
557   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
558   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
559   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
560   /* Setup problem */
561   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
562   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
563 
564   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
565   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
566   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
567 
568   ierr = DMSetNullSpaceConstructor(dm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
569   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
570 
571   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
572   {
573     PetscDS          ds;
574     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
575     void            *ctxs[3];
576 
577     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
578     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr);
579     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr);
580     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr);
581     ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
582     ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr);
583     ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr);
584   }
585   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
586   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
587   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
588 
589   if (user.showError) {
590     PetscDS          ds;
591     Vec              r;
592     PetscErrorCode (*exactFuncs[3])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
593     void            *ctxs[3];
594 
595     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
596     ierr = PetscDSGetExactSolution(ds, 0, &exactFuncs[0], &ctxs[0]);CHKERRQ(ierr);
597     ierr = PetscDSGetExactSolution(ds, 1, &exactFuncs[1], &ctxs[1]);CHKERRQ(ierr);
598     ierr = PetscDSGetExactSolution(ds, 2, &exactFuncs[2], &ctxs[2]);CHKERRQ(ierr);
599     ierr = DMGetGlobalVector(dm, &r);CHKERRQ(ierr);
600     ierr = DMProjectFunction(dm, 0.0, exactFuncs, ctxs, INSERT_ALL_VALUES, r);CHKERRQ(ierr);
601     ierr = VecAXPY(r, -1.0, u);CHKERRQ(ierr);
602     ierr = PetscObjectSetName((PetscObject) r, "Solution Error");CHKERRQ(ierr);
603     ierr = VecViewFromOptions(r, NULL, "-error_vec_view");CHKERRQ(ierr);
604     ierr = DMRestoreGlobalVector(dm, &r);CHKERRQ(ierr);
605   }
606   ierr = PetscObjectSetName((PetscObject) u, "Numerical Solution");CHKERRQ(ierr);
607   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
608 
609   ierr = VecDestroy(&u);CHKERRQ(ierr);
610   ierr = VecDestroy(&r);CHKERRQ(ierr);
611   ierr = DMDestroy(&dm);CHKERRQ(ierr);
612   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
613   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
614   ierr = PetscFinalize();
615   return ierr;
616 }
617 
618 /*TEST
619 
620   test:
621     suffix: 2d_tri_p2_p1_p1
622     requires: triangle !single
623     args: -dm_plex_separate_marker -sol_type quadratic -dm_refine 0 \
624       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
625       -dmsnes_check .001 -snes_error_if_not_converged \
626       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
627       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
628         -fieldsplit_0_pc_type lu \
629         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
630 
631   test:
632     # Using -dm_refine 2 -convest_num_refine 3 gives L_2 convergence rate: [2.9, 2.3, 1.9]
633     suffix: 2d_tri_p2_p1_p1_conv
634     requires: triangle !single
635     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
636       -vel_petscspace_degree 2 -pres_petscspace_degree 1 -temp_petscspace_degree 1 \
637       -snes_error_if_not_converged -snes_convergence_test correct_pressure -snes_convergence_estimate -convest_num_refine 1 \
638       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
639       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
640         -fieldsplit_0_pc_type lu \
641         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
642 
643   test:
644     suffix: 2d_tri_p3_p2_p2
645     requires: triangle !single
646     args: -dm_plex_separate_marker -sol_type cubic -dm_refine 0 \
647       -vel_petscspace_degree 3 -pres_petscspace_degree 2 -temp_petscspace_degree 2 \
648       -dmsnes_check .001 -snes_error_if_not_converged \
649       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
650       -pc_type fieldsplit -pc_fieldsplit_0_fields 0,2 -pc_fieldsplit_1_fields 1 -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
651         -fieldsplit_0_pc_type lu \
652         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
653 
654 TEST*/
655