xref: /petsc/src/ts/tutorials/ex46.c (revision 408cafa089b450f13fea5926aa80aaf73271e784)
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          dim;
22c4762a1bSJed Brown   PetscBool         simplex;
23c4762a1bSJed Brown   PetscInt          mms;
24c4762a1bSJed Brown   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
25c4762a1bSJed Brown } AppCtx;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown #define REYN 400.0
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /* MMS1
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   u = t + x^2 + y^2;
32c4762a1bSJed Brown   v = t + 2*x^2 - 2*x*y;
33c4762a1bSJed Brown   p = x + y - 1;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   f_x = -2*t*(x + y) + 2*x*y^2 - 4*x^2*y - 2*x^3 + 4.0/Re - 1.0
36c4762a1bSJed Brown   f_y = -2*t*x       + 2*y^3 - 4*x*y^2 - 2*x^2*y + 4.0/Re - 1.0
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   so that
39c4762a1bSJed Brown 
40c4762a1bSJed 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>
41c4762a1bSJed Brown                                                     + <-t (2x + 2y) + 2xy^2 - 4x^2y - 2x^3 + 4/Re - 1, -2xt + 2y^3 - 4xy^2 - 2x^2y + 4/Re - 1> = 0
42c4762a1bSJed Brown     \nabla \cdot u                                  = 2x - 2x = 0
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   where
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     <u, v> . <<u_x, v_x>, <u_y, v_y>> = <u u_x + v u_y, u v_x + v v_y>
47c4762a1bSJed Brown */
48c4762a1bSJed Brown PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
49c4762a1bSJed Brown {
50c4762a1bSJed Brown   u[0] = time + x[0]*x[0] + x[1]*x[1];
51c4762a1bSJed Brown   u[1] = time + 2.0*x[0]*x[0] - 2.0*x[0]*x[1];
52c4762a1bSJed Brown   return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
56c4762a1bSJed Brown {
57c4762a1bSJed Brown   *p = x[0] + x[1] - 1.0;
58c4762a1bSJed Brown   return 0;
59c4762a1bSJed Brown }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /* MMS 2*/
62c4762a1bSJed Brown 
63c4762a1bSJed Brown static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
64c4762a1bSJed Brown {
65c4762a1bSJed Brown   u[0] = PetscSinReal(time + x[0])*PetscSinReal(time + x[1]);
66c4762a1bSJed Brown   u[1] = PetscCosReal(time + x[0])*PetscCosReal(time + x[1]);
67c4762a1bSJed Brown   return 0;
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, void *ctx)
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   *p = PetscSinReal(time + x[0] - x[1]);
73c4762a1bSJed Brown   return 0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
77c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
78c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
79c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
80c4762a1bSJed Brown {
81c4762a1bSJed Brown   const PetscReal Re    = REYN;
82c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
83c4762a1bSJed Brown   PetscInt        c, d;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   for (c = 0; c < Ncomp; ++c) {
86c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
87c4762a1bSJed Brown       f0[c] += u[d] * u_x[c*dim+d];
88c4762a1bSJed Brown     }
89c4762a1bSJed Brown   }
90c4762a1bSJed Brown   f0[0] += u_t[0];
91c4762a1bSJed Brown   f0[1] += u_t[1];
92c4762a1bSJed Brown 
93c4762a1bSJed 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;
94c4762a1bSJed 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;
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
98c4762a1bSJed Brown                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
99c4762a1bSJed Brown                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
100c4762a1bSJed Brown                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
101c4762a1bSJed Brown {
102c4762a1bSJed Brown   const PetscReal Re    = REYN;
103c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
104c4762a1bSJed Brown   PetscInt        c, d;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   for (c = 0; c < Ncomp; ++c) {
107c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
108c4762a1bSJed Brown       f0[c] += u[d] * u_x[c*dim+d];
109c4762a1bSJed Brown     }
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown   f0[0] += u_t[0];
112c4762a1bSJed Brown   f0[1] += u_t[1];
113c4762a1bSJed Brown 
114c4762a1bSJed 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;
115c4762a1bSJed 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;
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
119c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
120c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
121c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   const PetscReal Re    = REYN;
124c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
125c4762a1bSJed Brown   PetscInt        comp, d;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   for (comp = 0; comp < Ncomp; ++comp) {
128c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
129c4762a1bSJed Brown       f1[comp*dim+d] = 1.0/Re * u_x[comp*dim+d];
130c4762a1bSJed Brown     }
131c4762a1bSJed Brown     f1[comp*dim+comp] -= u[Ncomp];
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
136c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
137c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
138c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   PetscInt d;
141c4762a1bSJed Brown   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
145c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
146c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
147c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
148c4762a1bSJed Brown {
149c4762a1bSJed Brown   PetscInt d;
150c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = 0.0;
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown /*
154c4762a1bSJed Brown   (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i)
155c4762a1bSJed Brown */
156c4762a1bSJed Brown static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
157c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
158c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
159c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
160c4762a1bSJed Brown {
161c4762a1bSJed Brown   PetscInt NcI = dim, NcJ = dim;
162c4762a1bSJed Brown   PetscInt fc, gc;
163c4762a1bSJed Brown   PetscInt d;
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
166c4762a1bSJed Brown     g0[d*dim+d] = u_tShift;
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   for (fc = 0; fc < NcI; ++fc) {
170c4762a1bSJed Brown     for (gc = 0; gc < NcJ; ++gc) {
171c4762a1bSJed Brown       g0[fc*NcJ+gc] += u_x[fc*NcJ+gc];
172c4762a1bSJed Brown     }
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown /*
177c4762a1bSJed Brown   (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
178c4762a1bSJed Brown */
179c4762a1bSJed Brown static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
180c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
181c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
182c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   PetscInt NcI = dim;
185c4762a1bSJed Brown   PetscInt NcJ = dim;
186c4762a1bSJed Brown   PetscInt fc, gc, dg;
187c4762a1bSJed Brown   for (fc = 0; fc < NcI; ++fc) {
188c4762a1bSJed Brown     for (gc = 0; gc < NcJ; ++gc) {
189c4762a1bSJed Brown       for (dg = 0; dg < dim; ++dg) {
190c4762a1bSJed Brown         /* kronecker delta */
191c4762a1bSJed Brown         if (fc == gc) {
192c4762a1bSJed Brown           g1[(fc*NcJ+gc)*dim+dg] += u[dg];
193c4762a1bSJed Brown         }
194c4762a1bSJed Brown       }
195c4762a1bSJed Brown     }
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /* < q, \nabla\cdot u >
200c4762a1bSJed Brown    NcompI = 1, NcompJ = dim */
201c4762a1bSJed Brown static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
202c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
203c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
204c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   PetscInt d;
207c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown /* -< \nabla\cdot v, p >
211c4762a1bSJed Brown     NcompI = dim, NcompJ = 1 */
212c4762a1bSJed Brown static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
213c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
214c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
215c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
216c4762a1bSJed Brown {
217c4762a1bSJed Brown   PetscInt d;
218c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T >
222c4762a1bSJed Brown    This just gives \nabla u, give the perdiagonal for the transpose */
223c4762a1bSJed Brown static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
224c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
225c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
226c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
227c4762a1bSJed Brown {
228c4762a1bSJed Brown   const PetscReal Re    = REYN;
229c4762a1bSJed Brown   const PetscInt  Ncomp = dim;
230c4762a1bSJed Brown   PetscInt        compI, d;
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   for (compI = 0; compI < Ncomp; ++compI) {
233c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
234c4762a1bSJed Brown       g3[((compI*Ncomp+compI)*dim+d)*dim+d] = 1.0/Re;
235c4762a1bSJed Brown     }
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
240c4762a1bSJed Brown {
241c4762a1bSJed Brown   PetscErrorCode ierr;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   PetscFunctionBeginUser;
244c4762a1bSJed Brown   options->dim     = 2;
245c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
246c4762a1bSJed Brown   options->mms     = 1;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex46.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex46.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
251c4762a1bSJed Brown   ierr = PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[])
257c4762a1bSJed Brown {
258*408cafa0SMatthew G. Knepley   DM             plex;
259c4762a1bSJed Brown   DMLabel        label;
260c4762a1bSJed Brown   PetscErrorCode ierr;
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   PetscFunctionBeginUser;
263c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
265*408cafa0SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
266*408cafa0SMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr);
267*408cafa0SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
268c4762a1bSJed Brown   PetscFunctionReturn(0);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
272c4762a1bSJed Brown {
273c4762a1bSJed Brown   DM             pdm = NULL;
274c4762a1bSJed Brown   const PetscInt dim = ctx->dim;
275c4762a1bSJed Brown   PetscBool      hasLabel;
276c4762a1bSJed Brown   PetscErrorCode ierr;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   PetscFunctionBeginUser;
279c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
281c4762a1bSJed Brown   /* If no boundary marker exists, mark the whole boundary */
282c4762a1bSJed Brown   ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr);
283c4762a1bSJed Brown   if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);}
284c4762a1bSJed Brown   /* Distribute mesh over processes */
285c4762a1bSJed Brown   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
286c4762a1bSJed Brown   if (pdm) {
287c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
288c4762a1bSJed Brown     *dm  = pdm;
289c4762a1bSJed Brown   }
290c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
291c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
292c4762a1bSJed Brown   PetscFunctionReturn(0);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
296c4762a1bSJed Brown {
297c4762a1bSJed Brown   PetscDS        prob;
298c4762a1bSJed Brown   const PetscInt id = 1;
299c4762a1bSJed Brown   PetscErrorCode ierr;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   PetscFunctionBeginUser;
302c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
303c4762a1bSJed Brown   switch (ctx->mms) {
304c4762a1bSJed Brown   case 1:
305c4762a1bSJed Brown     ierr = PetscDSSetResidual(prob, 0, f0_mms1_u, f1_u);CHKERRQ(ierr);break;
306c4762a1bSJed Brown   case 2:
307c4762a1bSJed Brown     ierr = PetscDSSetResidual(prob, 0, f0_mms2_u, f1_u);CHKERRQ(ierr);break;
308c4762a1bSJed Brown   }
309c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 1, f0_p, f1_p);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, g0_uu, g1_uu, NULL,  g3_uu);CHKERRQ(ierr);
311c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_up, NULL);CHKERRQ(ierr);
312c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL,  NULL);CHKERRQ(ierr);
313c4762a1bSJed Brown   switch (ctx->dim) {
314c4762a1bSJed Brown   case 2:
315c4762a1bSJed Brown     switch (ctx->mms) {
316c4762a1bSJed Brown     case 1:
317c4762a1bSJed Brown       ctx->exactFuncs[0] = mms1_u_2d;
318c4762a1bSJed Brown       ctx->exactFuncs[1] = mms1_p_2d;
319c4762a1bSJed Brown       break;
320c4762a1bSJed Brown     case 2:
321c4762a1bSJed Brown       ctx->exactFuncs[0] = mms2_u_2d;
322c4762a1bSJed Brown       ctx->exactFuncs[1] = mms2_p_2d;
323c4762a1bSJed Brown       break;
324c4762a1bSJed Brown     default:
325c4762a1bSJed Brown       SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %D", ctx->mms);
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown     break;
328c4762a1bSJed Brown   default:
329c4762a1bSJed Brown     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %D", ctx->dim);
330c4762a1bSJed Brown   }
331*408cafa0SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], 1, &id, ctx);CHKERRQ(ierr);
332c4762a1bSJed Brown   PetscFunctionReturn(0);
333c4762a1bSJed Brown }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
336c4762a1bSJed Brown {
337c4762a1bSJed Brown   DM              cdm = dm;
338c4762a1bSJed Brown   const PetscInt  dim = ctx->dim;
339c4762a1bSJed Brown   PetscFE         fe[2];
340c4762a1bSJed Brown   MPI_Comm        comm;
341c4762a1bSJed Brown   PetscErrorCode  ierr;
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   PetscFunctionBeginUser;
344c4762a1bSJed Brown   /* Create finite element */
345c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
346c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, dim, ctx->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
347c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, ctx->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
349c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
351c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
352c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
353c4762a1bSJed Brown   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
354c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
355c4762a1bSJed Brown   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
356c4762a1bSJed Brown   while (cdm) {
357c4762a1bSJed Brown     PetscObject  pressure;
358c4762a1bSJed Brown     MatNullSpace nsp;
359c4762a1bSJed Brown     PetscBool    hasLabel;
360c4762a1bSJed Brown 
361c4762a1bSJed Brown     ierr = DMGetField(cdm, 1, NULL, &pressure);CHKERRQ(ierr);
362c4762a1bSJed Brown     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp);CHKERRQ(ierr);
363c4762a1bSJed Brown     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nsp);CHKERRQ(ierr);
364c4762a1bSJed Brown     ierr = MatNullSpaceDestroy(&nsp);CHKERRQ(ierr);
365c4762a1bSJed Brown 
366c4762a1bSJed Brown     ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
367c4762a1bSJed Brown     if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
368*408cafa0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
369c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
370c4762a1bSJed Brown   }
371c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
372c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
373c4762a1bSJed Brown   PetscFunctionReturn(0);
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
376c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
377c4762a1bSJed Brown {
378c4762a1bSJed Brown   AppCtx        *user = (AppCtx *) ctx;
379c4762a1bSJed Brown   DM             dm;
380c4762a1bSJed Brown   PetscReal      ferrors[2];
381c4762a1bSJed Brown   PetscErrorCode ierr;
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   PetscFunctionBeginUser;
384c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
385c4762a1bSJed Brown   ierr = DMComputeL2FieldDiff(dm, crtime, user->exactFuncs, NULL, u, ferrors);CHKERRQ(ierr);
386c4762a1bSJed Brown   ierr = 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]);CHKERRQ(ierr);
387c4762a1bSJed Brown   PetscFunctionReturn(0);
388c4762a1bSJed Brown }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown int main(int argc, char **argv)
391c4762a1bSJed Brown {
392c4762a1bSJed Brown   AppCtx         ctx;
393c4762a1bSJed Brown   DM             dm;
394c4762a1bSJed Brown   TS             ts;
395c4762a1bSJed Brown   Vec            u, r;
396c4762a1bSJed Brown   PetscErrorCode ierr;
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
399c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
400c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr);
401c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
402c4762a1bSJed Brown   ierr = PetscMalloc1(2, &ctx.exactFuncs);CHKERRQ(ierr);
403c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
404c4762a1bSJed Brown   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
408c4762a1bSJed Brown 
409c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
410c4762a1bSJed Brown   ierr = TSMonitorSet(ts, MonitorError, &ctx, NULL);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
412c4762a1bSJed Brown   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
413c4762a1bSJed Brown   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
414c4762a1bSJed Brown   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
415c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
416c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   ierr = DMProjectFunction(dm, 0.0, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
419c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
420c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
421c4762a1bSJed Brown 
422c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
423c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
424c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
425c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
426c4762a1bSJed Brown   ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr);
427c4762a1bSJed Brown   ierr = PetscFinalize();
428c4762a1bSJed Brown   return ierr;
429c4762a1bSJed Brown }
430c4762a1bSJed Brown 
431c4762a1bSJed Brown /*TEST
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   # Full solves
434c4762a1bSJed Brown   test:
435c4762a1bSJed Brown     suffix: 2d_p2p1_r1
436c4762a1bSJed Brown     requires: !single triangle
437c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
438c4762a1bSJed Brown     args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
439c4762a1bSJed Brown   test:
440c4762a1bSJed Brown     suffix: 2d_p2p1_r2
441c4762a1bSJed Brown     requires: !single triangle
442c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
443c4762a1bSJed Brown     args: -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
444c4762a1bSJed Brown   test:
445c4762a1bSJed Brown     suffix: 2d_q2q1_r1
446c4762a1bSJed Brown     requires: !single
447c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
448c4762a1bSJed Brown     args: -simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
449c4762a1bSJed Brown   test:
450c4762a1bSJed Brown     suffix: 2d_q2q1_r2
451c4762a1bSJed Brown     requires: !single
452c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
453c4762a1bSJed Brown     args: -simplex 0 -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
454c4762a1bSJed Brown 
455c4762a1bSJed Brown TEST*/
456