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 */
mms1_u_2d(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)45*2a8381b2SBarry Smith PetscErrorCode mms1_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
46d71ae5a4SJacob Faibussowitsch {
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];
493ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
50c4762a1bSJed Brown }
51c4762a1bSJed Brown
mms1_p_2d(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * p,PetscCtx ctx)52*2a8381b2SBarry Smith PetscErrorCode mms1_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, PetscCtx ctx)
53d71ae5a4SJacob Faibussowitsch {
54c4762a1bSJed Brown *p = x[0] + x[1] - 1.0;
553ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
56c4762a1bSJed Brown }
57c4762a1bSJed Brown
58c4762a1bSJed Brown /* MMS 2*/
59c4762a1bSJed Brown
mms2_u_2d(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)60*2a8381b2SBarry Smith static PetscErrorCode mms2_u_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
61d71ae5a4SJacob Faibussowitsch {
62c4762a1bSJed Brown u[0] = PetscSinReal(time + x[0]) * PetscSinReal(time + x[1]);
63c4762a1bSJed Brown u[1] = PetscCosReal(time + x[0]) * PetscCosReal(time + x[1]);
643ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown
mms2_p_2d(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * p,PetscCtx ctx)67*2a8381b2SBarry Smith static PetscErrorCode mms2_p_2d(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *p, PetscCtx ctx)
68d71ae5a4SJacob Faibussowitsch {
69c4762a1bSJed Brown *p = PetscSinReal(time + x[0] - x[1]);
703ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown
f0_mms1_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])73d71ae5a4SJacob Faibussowitsch static void f0_mms1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
74d71ae5a4SJacob Faibussowitsch {
75c4762a1bSJed Brown const PetscReal Re = REYN;
76c4762a1bSJed Brown const PetscInt Ncomp = dim;
77c4762a1bSJed Brown PetscInt c, d;
78c4762a1bSJed Brown
79c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) {
80ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
81c4762a1bSJed Brown }
82c4762a1bSJed Brown f0[0] += u_t[0];
83c4762a1bSJed Brown f0[1] += u_t[1];
84c4762a1bSJed Brown
85c4762a1bSJed 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;
86c4762a1bSJed 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;
87c4762a1bSJed Brown }
88c4762a1bSJed Brown
f0_mms2_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])89d71ae5a4SJacob Faibussowitsch static void f0_mms2_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
90d71ae5a4SJacob Faibussowitsch {
91c4762a1bSJed Brown const PetscReal Re = REYN;
92c4762a1bSJed Brown const PetscInt Ncomp = dim;
93c4762a1bSJed Brown PetscInt c, d;
94c4762a1bSJed Brown
95c4762a1bSJed Brown for (c = 0; c < Ncomp; ++c) {
96ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f0[c] += u[d] * u_x[c * dim + d];
97c4762a1bSJed Brown }
98c4762a1bSJed Brown f0[0] += u_t[0];
99c4762a1bSJed Brown f0[1] += u_t[1];
100c4762a1bSJed Brown
101c4762a1bSJed 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;
102c4762a1bSJed 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;
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
f1_u(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])105d71ae5a4SJacob Faibussowitsch static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
106d71ae5a4SJacob Faibussowitsch {
107c4762a1bSJed Brown const PetscReal Re = REYN;
108c4762a1bSJed Brown const PetscInt Ncomp = dim;
109c4762a1bSJed Brown PetscInt comp, d;
110c4762a1bSJed Brown
111c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) {
112ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = 1.0 / Re * u_x[comp * dim + d];
113c4762a1bSJed Brown f1[comp * dim + comp] -= u[Ncomp];
114c4762a1bSJed Brown }
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
f0_p(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])117d71ae5a4SJacob Faibussowitsch static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
118d71ae5a4SJacob Faibussowitsch {
119c4762a1bSJed Brown PetscInt d;
120c4762a1bSJed Brown for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d * dim + d];
121c4762a1bSJed Brown }
122c4762a1bSJed Brown
f1_p(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])123d71ae5a4SJacob Faibussowitsch static void f1_p(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown PetscInt d;
126c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = 0.0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown
129c4762a1bSJed Brown /*
130c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \phi_j grad_j u_i)
131c4762a1bSJed Brown */
g0_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g0[])132d71ae5a4SJacob Faibussowitsch static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
133d71ae5a4SJacob Faibussowitsch {
134c4762a1bSJed Brown PetscInt NcI = dim, NcJ = dim;
135c4762a1bSJed Brown PetscInt fc, gc;
136c4762a1bSJed Brown PetscInt d;
137c4762a1bSJed Brown
138ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g0[d * dim + d] = u_tShift;
139c4762a1bSJed Brown
140c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) {
141ad540459SPierre Jolivet for (gc = 0; gc < NcJ; ++gc) g0[fc * NcJ + gc] += u_x[fc * NcJ + gc];
142c4762a1bSJed Brown }
143c4762a1bSJed Brown }
144c4762a1bSJed Brown
145c4762a1bSJed Brown /*
146c4762a1bSJed Brown (psi_i, u_j grad_j u_i) ==> (\psi_i, \u_j grad_j \phi_i)
147c4762a1bSJed Brown */
g1_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g1[])148d71ae5a4SJacob Faibussowitsch static void g1_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
149d71ae5a4SJacob Faibussowitsch {
150c4762a1bSJed Brown PetscInt NcI = dim;
151c4762a1bSJed Brown PetscInt NcJ = dim;
152c4762a1bSJed Brown PetscInt fc, gc, dg;
153c4762a1bSJed Brown for (fc = 0; fc < NcI; ++fc) {
154c4762a1bSJed Brown for (gc = 0; gc < NcJ; ++gc) {
155c4762a1bSJed Brown for (dg = 0; dg < dim; ++dg) {
156c4762a1bSJed Brown /* kronecker delta */
157ad540459SPierre Jolivet if (fc == gc) g1[(fc * NcJ + gc) * dim + dg] += u[dg];
158c4762a1bSJed Brown }
159c4762a1bSJed Brown }
160c4762a1bSJed Brown }
161c4762a1bSJed Brown }
162c4762a1bSJed Brown
163c4762a1bSJed Brown /* < q, \nabla\cdot u >
164c4762a1bSJed Brown NcompI = 1, NcompJ = dim */
g1_pu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g1[])165d71ae5a4SJacob Faibussowitsch static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
166d71ae5a4SJacob Faibussowitsch {
167c4762a1bSJed Brown PetscInt d;
168c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d * dim + d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
169c4762a1bSJed Brown }
170c4762a1bSJed Brown
171c4762a1bSJed Brown /* -< \nabla\cdot v, p >
172c4762a1bSJed Brown NcompI = dim, NcompJ = 1 */
g2_up(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g2[])173d71ae5a4SJacob Faibussowitsch static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
174d71ae5a4SJacob Faibussowitsch {
175c4762a1bSJed Brown PetscInt d;
176c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d * dim + d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
177c4762a1bSJed Brown }
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T >
180c4762a1bSJed Brown This just gives \nabla u, give the perdiagonal for the transpose */
g3_uu(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar g3[])181d71ae5a4SJacob Faibussowitsch static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
182d71ae5a4SJacob Faibussowitsch {
183c4762a1bSJed Brown const PetscReal Re = REYN;
184c4762a1bSJed Brown const PetscInt Ncomp = dim;
185c4762a1bSJed Brown PetscInt compI, d;
186c4762a1bSJed Brown
187c4762a1bSJed Brown for (compI = 0; compI < Ncomp; ++compI) {
188ad540459SPierre Jolivet for (d = 0; d < dim; ++d) g3[((compI * Ncomp + compI) * dim + d) * dim + d] = 1.0 / Re;
189c4762a1bSJed Brown }
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)192d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
193d71ae5a4SJacob Faibussowitsch {
194c4762a1bSJed Brown PetscFunctionBeginUser;
195c4762a1bSJed Brown options->mms = 1;
196c4762a1bSJed Brown
197d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Navier-Stokes Equation Options", "DMPLEX");
1989566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mms", "The manufactured solution to use", "ex46.c", options->mms, &options->mms, NULL));
199d0609cedSBarry Smith PetscOptionsEnd();
2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * ctx)203d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
204d71ae5a4SJacob Faibussowitsch {
205c4762a1bSJed Brown PetscFunctionBeginUser;
2069566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
2079566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
2089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
2099566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
211c4762a1bSJed Brown }
212c4762a1bSJed Brown
SetupProblem(DM dm,AppCtx * ctx)213d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
214d71ae5a4SJacob Faibussowitsch {
21545480ffeSMatthew G. Knepley PetscDS ds;
21645480ffeSMatthew G. Knepley DMLabel label;
217c4762a1bSJed Brown const PetscInt id = 1;
21830602db0SMatthew G. Knepley PetscInt dim;
219c4762a1bSJed Brown
220c4762a1bSJed Brown PetscFunctionBeginUser;
2219566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
2229566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2239566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
22430602db0SMatthew G. Knepley switch (dim) {
22530602db0SMatthew G. Knepley case 2:
226c4762a1bSJed Brown switch (ctx->mms) {
227c4762a1bSJed Brown case 1:
2289566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms1_u, f1_u));
2299566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
2309566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
2319566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
2329566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
2339566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms1_u_2d, ctx));
2349566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms1_p_2d, ctx));
23557d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms1_u_2d, NULL, ctx, NULL));
236c4762a1bSJed Brown break;
237c4762a1bSJed Brown case 2:
2389566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_mms2_u, f1_u));
2399566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 1, f0_p, f1_p));
2409566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_uu, g1_uu, NULL, g3_uu));
2419566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_up, NULL));
2429566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL, NULL));
2439566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms2_u_2d, ctx));
2449566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 1, mms2_p_2d, ctx));
24557d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms2_u_2d, NULL, ctx, NULL));
246c4762a1bSJed Brown break;
247d71ae5a4SJacob Faibussowitsch default:
248d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid MMS %" PetscInt_FMT, ctx->mms);
249c4762a1bSJed Brown }
250c4762a1bSJed Brown break;
251d71ae5a4SJacob Faibussowitsch default:
252d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension %" PetscInt_FMT, dim);
253c4762a1bSJed Brown }
2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown
SetupDiscretization(DM dm,AppCtx * ctx)257d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
258d71ae5a4SJacob Faibussowitsch {
259c4762a1bSJed Brown MPI_Comm comm;
26030602db0SMatthew G. Knepley DM cdm = dm;
26130602db0SMatthew G. Knepley PetscFE fe[2];
26230602db0SMatthew G. Knepley PetscInt dim;
26330602db0SMatthew G. Knepley PetscBool simplex;
264c4762a1bSJed Brown
265c4762a1bSJed Brown PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2679566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
2689566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
2699566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", PETSC_DEFAULT, &fe[0]));
2709566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[0], "velocity"));
2719566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "pres_", PETSC_DEFAULT, &fe[1]));
2729566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe[0], fe[1]));
2739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe[1], "pressure"));
274c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */
2759566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0]));
2769566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1]));
2779566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
2789566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx));
279c4762a1bSJed Brown while (cdm) {
280c4762a1bSJed Brown PetscObject pressure;
281c4762a1bSJed Brown MatNullSpace nsp;
282c4762a1bSJed Brown
2839566063dSJacob Faibussowitsch PetscCall(DMGetField(cdm, 1, NULL, &pressure));
2849566063dSJacob Faibussowitsch PetscCall(MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nsp));
2859566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(pressure, "nullspace", (PetscObject)nsp));
2869566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nsp));
287c4762a1bSJed Brown
2889566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
2899566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
290c4762a1bSJed Brown }
2919566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[0]));
2929566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe[1]));
2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
294c4762a1bSJed Brown }
295c4762a1bSJed Brown
MonitorError(TS ts,PetscInt step,PetscReal crtime,Vec u,PetscCtx ctx)296*2a8381b2SBarry Smith static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, PetscCtx ctx)
297d71ae5a4SJacob Faibussowitsch {
2988434afd1SBarry Smith PetscSimplePointFn *funcs[2];
29930602db0SMatthew G. Knepley void *ctxs[2];
300c4762a1bSJed Brown DM dm;
30130602db0SMatthew G. Knepley PetscDS ds;
302c4762a1bSJed Brown PetscReal ferrors[2];
303c4762a1bSJed Brown
304c4762a1bSJed Brown PetscFunctionBeginUser;
3059566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
3069566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3079566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
3089566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
3099566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, crtime, funcs, ctxs, u, ferrors));
3109566063dSJacob Faibussowitsch PetscCall(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]));
3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
312c4762a1bSJed Brown }
313c4762a1bSJed Brown
main(int argc,char ** argv)314d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
315d71ae5a4SJacob Faibussowitsch {
316c4762a1bSJed Brown AppCtx ctx;
317c4762a1bSJed Brown DM dm;
318c4762a1bSJed Brown TS ts;
319c4762a1bSJed Brown Vec u, r;
320c4762a1bSJed Brown
321327415f7SBarry Smith PetscFunctionBeginUser;
3229566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3239566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3249566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
3259566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx));
3269566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx));
3279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateClosureIndex(dm, NULL));
328c4762a1bSJed Brown
3299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
3309566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r));
331c4762a1bSJed Brown
3329566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3339566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MonitorError, &ctx, NULL));
3349566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm));
3359566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
3369566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
3379566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
3389566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3399566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
3409566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u));
341c4762a1bSJed Brown
34230602db0SMatthew G. Knepley {
3438434afd1SBarry Smith PetscSimplePointFn *funcs[2];
34430602db0SMatthew G. Knepley void *ctxs[2];
34530602db0SMatthew G. Knepley PetscDS ds;
34630602db0SMatthew G. Knepley
3479566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3489566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &funcs[0], &ctxs[0]));
3499566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 1, &funcs[1], &ctxs[1]));
3509566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, funcs, ctxs, INSERT_ALL_VALUES, u));
35130602db0SMatthew G. Knepley }
3529566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u));
3539566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
354c4762a1bSJed Brown
3559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
3569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
3579566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
3589566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
3599566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
360b122ec5aSJacob Faibussowitsch return 0;
361c4762a1bSJed Brown }
362c4762a1bSJed Brown
363c4762a1bSJed Brown /*TEST
364c4762a1bSJed Brown
365c4762a1bSJed Brown # Full solves
366c4762a1bSJed Brown test:
367c4762a1bSJed Brown suffix: 2d_p2p1_r1
368c4762a1bSJed Brown requires: !single triangle
369c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
37030602db0SMatthew G. Knepley args: -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
371188af4bfSBarry Smith -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ts_monitor -dmts_check \
37230602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \
37330602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \
37430602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
37530602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \
37630602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
37730602db0SMatthew G. Knepley
378c4762a1bSJed Brown test:
379c4762a1bSJed Brown suffix: 2d_q2q1_r1
380c4762a1bSJed Brown requires: !single
381c4762a1bSJed Brown filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g" -e "s~ 0\]~ 0.0\]~g"
38230602db0SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
383188af4bfSBarry Smith -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ts_monitor -dmts_check \
38430602db0SMatthew G. Knepley -snes_monitor_short -snes_converged_reason \
38530602db0SMatthew G. Knepley -ksp_monitor_short -ksp_converged_reason \
38630602db0SMatthew G. Knepley -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \
38730602db0SMatthew G. Knepley -fieldsplit_velocity_pc_type lu \
38830602db0SMatthew G. Knepley -fieldsplit_pressure_ksp_rtol 1.0e-10 -fieldsplit_pressure_pc_type jacobi
389c4762a1bSJed Brown
390c4762a1bSJed Brown TEST*/
391