xref: /petsc/src/ts/tutorials/ex46.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Time dependent Navier-Stokes problem in 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the Navier-Stokes in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown This example supports discretized auxiliary fields (Re) as well as\n\
5c4762a1bSJed Brown multilevel nonlinear solvers.\n\
6c4762a1bSJed Brown Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscdmplex.h>
9c4762a1bSJed Brown #include <petscsnes.h>
10c4762a1bSJed Brown #include <petscts.h>
11c4762a1bSJed Brown #include <petscds.h>
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown   Navier-Stokes equation:
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   du/dt + u . grad u - \Delta u - grad p = f
17c4762a1bSJed Brown   div u  = 0
18c4762a1bSJed Brown */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscInt mms;
22c4762a1bSJed Brown } AppCtx;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #define REYN 400.0
25c4762a1bSJed Brown 
26c4762a1bSJed Brown /* MMS1
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   u = t + x^2 + y^2;
29c4762a1bSJed Brown   v = t + 2*x^2 - 2*x*y;
30c4762a1bSJed Brown   p = x + y - 1;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   f_x = -2*t*(x + y) + 2*x*y^2 - 4*x^2*y - 2*x^3 + 4.0/Re - 1.0
33c4762a1bSJed Brown   f_y = -2*t*x       + 2*y^3 - 4*x*y^2 - 2*x^2*y + 4.0/Re - 1.0
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   so that
36c4762a1bSJed Brown 
37c4762a1bSJed Brown     u_t + u \cdot \nabla u - 1/Re \Delta u + \nabla p + f = <1, 1> + <t (2x + 2y) + 2x^3 + 4x^2y - 2xy^2, t 2x + 2x^2y + 4xy^2 - 2y^3> - 1/Re <4, 4> + <1, 1>
38c4762a1bSJed Brown                                                     + <-t (2x + 2y) + 2xy^2 - 4x^2y - 2x^3 + 4/Re - 1, -2xt + 2y^3 - 4xy^2 - 2x^2y + 4/Re - 1> = 0
39c4762a1bSJed Brown     \nabla \cdot u                                  = 2x - 2x = 0
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   where
42c4762a1bSJed Brown 
43c4762a1bSJed Brown     <u, v> . <<u_x, v_x>, <u_y, v_y>> = <u u_x + v u_y, u v_x + v v_y>
44c4762a1bSJed Brown */
45c4762a1bSJed Brown PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
46c4762a1bSJed Brown {
47c4762a1bSJed Brown   u[0] = time + x[0]*x[0] + x[1]*x[1];
48c4762a1bSJed Brown   u[1] = time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
49c4762a1bSJed Brown   return 0;
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   *p = x[0] + x[1] - 1.0;
55c4762a1bSJed Brown   return 0;
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown /* MMS 2*/
59c4762a1bSJed Brown 
60c4762a1bSJed Brown static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   u[0] = PetscSinReal(time + x[0])*PetscSinReal(time + x[1]);
63c4762a1bSJed Brown   u[1] = PetscCosReal(time + x[0])*PetscCosReal(time + x[1]);
64c4762a1bSJed Brown   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   *p = PetscSinReal(time + x[0] - x[1]);
70c4762a1bSJed Brown   return 0;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
74c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
75c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
76c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   const PetscReal Re    = REYN;
79c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
80c4762a1bSJed Brown   PetscInt        c, d;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   for (c = 0; c < Ncomp; ++c) {
83c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
84c4762a1bSJed Brown       f0[c] += u[d] * u_x[c*dim+d];
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown   }
87c4762a1bSJed Brown   f0[0] += u_t[0];
88c4762a1bSJed Brown   f0[1] += u_t[1];
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   f0[0] += -2.0*t*(x[0] + x[1]) + 2.0*x[0]*x[1]*x[1] - 4.0*x[0]*x[0]*x[1] - 2.0*x[0]*x[0]*x[0] + 4.0/Re - 1.0;
91c4762a1bSJed Brown   f0[1] += -2.0*t*x[0]          + 2.0*x[1]*x[1]*x[1] - 4.0*x[0]*x[1]*x[1] - 2.0*x[0]*x[0]*x[1] + 4.0/Re - 1.0;
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
94c4762a1bSJed Brown static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   const PetscReal Re    = REYN;
100c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
101c4762a1bSJed Brown   PetscInt        c, d;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   for (c = 0; c < Ncomp; ++c) {
104c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
105c4762a1bSJed Brown       f0[c] += u[d] * u_x[c*dim+d];
106c4762a1bSJed Brown     }
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown   f0[0] += u_t[0];
109c4762a1bSJed Brown   f0[1] += u_t[1];
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   f0[0] -= ( Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[0]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscSinReal(t + x[0])*PetscSinReal(t + x[1]))/Re;
112c4762a1bSJed Brown   f0[1] -= (-Re*((1.0L/2.0L)*PetscSinReal(2*t + 2*x[1]) + PetscSinReal(2*t + x[0] + x[1]) + PetscCosReal(t + x[0] - x[1])) + 2.0*PetscCosReal(t + x[0])*PetscCosReal(t + x[1]))/Re;
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
116c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
117c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
118c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
119c4762a1bSJed Brown {
120c4762a1bSJed Brown   const PetscReal Re    = REYN;
121c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
122c4762a1bSJed Brown   PetscInt        comp, d;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
125c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
126c4762a1bSJed Brown       f1[comp*dim+d] = 1.0/Re * u_x[comp*dim+d];
127c4762a1bSJed Brown     }
128c4762a1bSJed Brown     f1[comp*dim+comp] -= u[Ncomp];
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
133c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
134c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
135c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   PetscInt d;
138c4762a1bSJed Brown   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
142c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
143c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
144c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
145c4762a1bSJed Brown {
146c4762a1bSJed Brown   PetscInt d;
147c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = 0.0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /*
151c4762a1bSJed Brown   (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i)
152c4762a1bSJed Brown */
153c4762a1bSJed Brown static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
154c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
155c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
156c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
157c4762a1bSJed Brown {
158c4762a1bSJed Brown   PetscInt NcI = dim, NcJ = dim;
159c4762a1bSJed Brown   PetscInt fc, gc;
160c4762a1bSJed Brown   PetscInt d;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
163c4762a1bSJed Brown     g0[d*dim+d] = u_tShift;
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   for (fc = 0; fc < NcI; ++fc) {
167c4762a1bSJed Brown     for (gc = 0; gc < NcJ; ++gc) {
168c4762a1bSJed Brown       g0[fc*NcJ+gc] += u_x[fc*NcJ+gc];
169c4762a1bSJed Brown     }
170c4762a1bSJed Brown   }
171c4762a1bSJed Brown }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown /*
174c4762a1bSJed Brown   (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
175c4762a1bSJed Brown */
176c4762a1bSJed Brown static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
177c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
178c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
179c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
180c4762a1bSJed Brown {
181c4762a1bSJed Brown   PetscInt NcI = dim;
182c4762a1bSJed Brown   PetscInt NcJ = dim;
183c4762a1bSJed Brown   PetscInt fc, gc, dg;
184c4762a1bSJed Brown   for (fc = 0; fc < NcI; ++fc) {
185c4762a1bSJed Brown     for (gc = 0; gc < NcJ; ++gc) {
186c4762a1bSJed Brown       for (dg = 0; dg < dim; ++dg) {
187c4762a1bSJed Brown         /* kronecker delta */
188c4762a1bSJed Brown         if (fc == gc) {
189c4762a1bSJed Brown           g1[(fc*NcJ+gc)*dim+dg] += u[dg];
190c4762a1bSJed Brown         }
191c4762a1bSJed Brown       }
192c4762a1bSJed Brown     }
193c4762a1bSJed Brown   }
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /* < q, \nabla\cdot u >
197c4762a1bSJed Brown    NcompI = 1, NcompJ = dim */
198c4762a1bSJed Brown static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
199c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
200c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
201c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   PetscInt d;
204c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
205c4762a1bSJed Brown }
206c4762a1bSJed Brown 
207c4762a1bSJed Brown /* -< \nabla\cdot v, p >
208c4762a1bSJed Brown     NcompI = dim, NcompJ = 1 */
209c4762a1bSJed Brown static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
210c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
211c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
212c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
213c4762a1bSJed Brown {
214c4762a1bSJed Brown   PetscInt d;
215c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T >
219c4762a1bSJed Brown    This just gives \nabla u, give the perdiagonal for the transpose */
220c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
224c4762a1bSJed Brown {
225c4762a1bSJed Brown   const PetscReal Re    = REYN;
226c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
227c4762a1bSJed Brown   PetscInt        compI, d;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   for (compI = 0; compI < Ncomp; ++compI) {
230c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
231c4762a1bSJed Brown       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0/Re;
232c4762a1bSJed Brown     }
233c4762a1bSJed Brown   }
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
237c4762a1bSJed Brown {
238c4762a1bSJed Brown   PetscErrorCode ierr;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   PetscFunctionBeginUser;
241c4762a1bSJed Brown   options->mms = 1;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");CHKERRQ(ierr);
244*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL));
245c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
246c4762a1bSJed Brown   PetscFunctionReturn(0);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
250c4762a1bSJed Brown {
251c4762a1bSJed Brown   PetscFunctionBeginUser;
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
255*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
256c4762a1bSJed Brown   PetscFunctionReturn(0);
257c4762a1bSJed Brown }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
260c4762a1bSJed Brown {
26145480ffeSMatthew G. Knepley   PetscDS        ds;
26245480ffeSMatthew G. Knepley   DMLabel        label;
263c4762a1bSJed Brown   const PetscInt id = 1;
26430602db0SMatthew G. Knepley   PetscInt       dim;
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   PetscFunctionBeginUser;
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
27030602db0SMatthew G. Knepley   switch (dim) {
27130602db0SMatthew G. Knepley   case 2:
272c4762a1bSJed Brown     switch (ctx->mms) {
273c4762a1bSJed Brown     case 1:
274*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u));
275*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 1, f0_p, f1_p));
276*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL,  g3_uu));
277*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
278*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
279*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx));
280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx));
281*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms1_u_2d, NULL, ctx, NULL));
282c4762a1bSJed Brown       break;
283c4762a1bSJed Brown     case 2:
284*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u));
285*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 1, f0_p, f1_p));
286*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL,  g3_uu));
287*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
288*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL));
289*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx));
290*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx));
291*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms2_u_2d, NULL, ctx, NULL));
292c4762a1bSJed Brown       break;
29398921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %D", ctx->mms);
294c4762a1bSJed Brown     }
295c4762a1bSJed Brown     break;
29698921bdaSJacob Faibussowitsch   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %D", dim);
297c4762a1bSJed Brown   }
298c4762a1bSJed Brown   PetscFunctionReturn(0);
299c4762a1bSJed Brown }
300c4762a1bSJed Brown 
301c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
302c4762a1bSJed Brown {
303c4762a1bSJed Brown   MPI_Comm        comm;
30430602db0SMatthew G. Knepley   DM              cdm = dm;
30530602db0SMatthew G. Knepley   PetscFE         fe[2];
30630602db0SMatthew G. Knepley   PetscInt        dim;
30730602db0SMatthew G. Knepley   PetscBool       simplex;
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   PetscFunctionBeginUser;
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexIsSimplex(dm, &simplex));
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe[0], "velocity"));
315*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
316*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECopyQuadrature(fe[0], fe[1]));
317*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe[1], "pressure"));
318c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
319*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe[0]));
320*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 1, NULL, (PetscObject) fe[1]));
321*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
322*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupProblem(dm, ctx));
323c4762a1bSJed Brown   while (cdm) {
324c4762a1bSJed Brown     PetscObject  pressure;
325c4762a1bSJed Brown     MatNullSpace nsp;
326c4762a1bSJed Brown 
327*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetField(cdm, 1, NULL, &pressure));
328*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp));
329*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectCompose(pressure, "nullspace", (PetscObject) nsp));
330*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNullSpaceDestroy(&nsp));
331c4762a1bSJed Brown 
332*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, cdm));
333*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
334c4762a1bSJed Brown   }
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe[0]));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe[1]));
337c4762a1bSJed Brown   PetscFunctionReturn(0);
338c4762a1bSJed Brown }
339c4762a1bSJed Brown 
340c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
341c4762a1bSJed Brown {
34230602db0SMatthew G. Knepley   PetscSimplePointFunc funcs[2];
34330602db0SMatthew G. Knepley   void                *ctxs[2];
344c4762a1bSJed Brown   DM                   dm;
34530602db0SMatthew G. Knepley   PetscDS              ds;
346c4762a1bSJed Brown   PetscReal            ferrors[2];
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   PetscFunctionBeginUser;
349*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
350*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
351*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
352*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
353*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors));
354*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: [%2.3g, %2.3g]\n", (int) step, (double) crtime, (double) ferrors[0], (double) ferrors[1]));
355c4762a1bSJed Brown   PetscFunctionReturn(0);
356c4762a1bSJed Brown }
357c4762a1bSJed Brown 
358c4762a1bSJed Brown int main(int argc, char **argv)
359c4762a1bSJed Brown {
360c4762a1bSJed Brown   AppCtx         ctx;
361c4762a1bSJed Brown   DM             dm;
362c4762a1bSJed Brown   TS             ts;
363c4762a1bSJed Brown   Vec            u, r;
364c4762a1bSJed Brown   PetscErrorCode ierr;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
367*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx));
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
369*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(dm, &ctx));
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, &ctx));
371*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateClosureIndex(dm, NULL));
372c4762a1bSJed Brown 
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u, &r));
375c4762a1bSJed Brown 
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts));
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorSet(ts, MonitorError, &ctx, NULL));
378*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, dm));
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
380*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
381*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSCheckFromOptions(ts, u));
385c4762a1bSJed Brown 
38630602db0SMatthew G. Knepley   {
38730602db0SMatthew G. Knepley     PetscSimplePointFunc funcs[2];
38830602db0SMatthew G. Knepley     void                *ctxs[2];
38930602db0SMatthew G. Knepley     PetscDS              ds;
39030602db0SMatthew G. Knepley 
391*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDS(dm, &ds));
392*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
393*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
394*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u));
39530602db0SMatthew G. Knepley   }
396*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts, u));
397*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(u, NULL, "-sol_vec_view"));
398c4762a1bSJed Brown 
399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
400*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
401*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
403c4762a1bSJed Brown   ierr = PetscFinalize();
404c4762a1bSJed Brown   return ierr;
405c4762a1bSJed Brown }
406c4762a1bSJed Brown 
407c4762a1bSJed Brown /*TEST
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   # Full solves
410c4762a1bSJed Brown   test:
411c4762a1bSJed Brown     suffix: 2d_p2p1_r1
412c4762a1bSJed Brown     requires: !single triangle
413c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
41430602db0SMatthew G. Knepley     args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
41530602db0SMatthew G. Knepley           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
41630602db0SMatthew G. Knepley           -snes_monitor_short -snes_converged_reason \
41730602db0SMatthew G. Knepley           -ksp_monitor_short -ksp_converged_reason \
41830602db0SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
41930602db0SMatthew G. Knepley             -fieldsplit_velocity_pc_type lu \
42030602db0SMatthew G. Knepley             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
42130602db0SMatthew G. Knepley 
422c4762a1bSJed Brown   test:
423c4762a1bSJed Brown     suffix: 2d_q2q1_r1
424c4762a1bSJed Brown     requires: !single
425c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
42630602db0SMatthew G. Knepley     args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
42730602db0SMatthew G. Knepley           -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ts_monitor -dmts_check \
42830602db0SMatthew G. Knepley           -snes_monitor_short -snes_converged_reason \
42930602db0SMatthew G. Knepley           -ksp_monitor_short -ksp_converged_reason \
43030602db0SMatthew G. Knepley           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
43130602db0SMatthew G. Knepley             -fieldsplit_velocity_pc_type lu \
43230602db0SMatthew G. Knepley             -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
433c4762a1bSJed Brown 
434c4762a1bSJed Brown TEST*/
435