xref: /petsc/src/snes/tutorials/ex71.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Poiseuille Flow in 2d and 3d channels with finite elements.\n\
2*c4762a1bSJed Brown We solve the Poiseuille flow problem in a rectangular\n\
3*c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown /*F
6*c4762a1bSJed Brown A Poiseuille flow is a steady-state isoviscous Stokes flow in a pipe of constant cross-section. We discretize using the
7*c4762a1bSJed Brown finite element method on an unstructured mesh. The weak form equations are
8*c4762a1bSJed Brown \begin{align*}
9*c4762a1bSJed Brown   < \nabla v, \nu (\nabla u + {\nabla u}^T) > - < \nabla\cdot v, p > + < v, \Delta \hat n >_{\Gamma_o} = 0
10*c4762a1bSJed Brown   < q, \nabla\cdot u >                                                                                 = 0
11*c4762a1bSJed Brown \end{align*}
12*c4762a1bSJed Brown where $\nu$ is the kinematic viscosity, $\Delta$ is the pressure drop per unit length, assuming that pressure is 0 on
13*c4762a1bSJed Brown the left edge, and $\Gamma_o$ is the outlet boundary at the right edge of the pipe. The normal velocity will be zero at
14*c4762a1bSJed Brown the wall, but we will allow a fixed tangential velocity $u_0$.
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown In order to test our global to local basis transformation, we will allow the pipe to be at an angle $\alpha$ to the
17*c4762a1bSJed Brown coordinate axes.
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown For visualization, use
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   -dm_view hdf5:$PWD/sol.h5 -sol_vec_view hdf5:$PWD/sol.h5::append -exact_vec_view hdf5:$PWD/sol.h5::append
22*c4762a1bSJed Brown F*/
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown #include <petscdmplex.h>
25*c4762a1bSJed Brown #include <petscsnes.h>
26*c4762a1bSJed Brown #include <petscds.h>
27*c4762a1bSJed Brown #include <petscbag.h>
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown typedef struct {
30*c4762a1bSJed Brown   PetscReal Delta; /* Pressure drop per unit length */
31*c4762a1bSJed Brown   PetscReal nu;    /* Kinematic viscosity */
32*c4762a1bSJed Brown   PetscReal u_0;   /* Tangential velocity at the wall */
33*c4762a1bSJed Brown   PetscReal alpha; /* Angle of pipe wall to x-axis */
34*c4762a1bSJed Brown } Parameter;
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown typedef struct {
37*c4762a1bSJed Brown   /* Domain and mesh definition */
38*c4762a1bSJed Brown   PetscInt  dim;               /* The topological mesh dimension */
39*c4762a1bSJed Brown   PetscBool simplex;           /* Use simplices or tensor product cells */
40*c4762a1bSJed Brown   PetscInt  cells[3];          /* The initial domain division */
41*c4762a1bSJed Brown   /* Problem definition */
42*c4762a1bSJed Brown   PetscBag  bag;               /* Holds problem parameters */
43*c4762a1bSJed Brown   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
44*c4762a1bSJed Brown } AppCtx;
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown /*
47*c4762a1bSJed Brown   In 2D, plane Poiseuille flow has exact solution:
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown     u = \Delta/(2 \nu) y (1 - y) + u_0
50*c4762a1bSJed Brown     v = 0
51*c4762a1bSJed Brown     p = -\Delta x
52*c4762a1bSJed Brown     f = 0
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown   so that
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown     -\nu \Delta u + \nabla p + f = <\Delta, 0> + <-\Delta, 0> + <0, 0> = 0
57*c4762a1bSJed Brown     \nabla \cdot u               = 0 + 0                               = 0
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   In 3D we use exact solution:
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown     u = \Delta/(4 \nu) (y (1 - y) + z (1 - z)) + u_0
62*c4762a1bSJed Brown     v = 0
63*c4762a1bSJed Brown     w = 0
64*c4762a1bSJed Brown     p = -\Delta x
65*c4762a1bSJed Brown     f = 0
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   so that
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown     -\nu \Delta u + \nabla p + f = <Delta, 0, 0> + <-Delta, 0, 0> + <0, 0, 0> = 0
70*c4762a1bSJed Brown     \nabla \cdot u               = 0 + 0 + 0                                  = 0
71*c4762a1bSJed Brown 
72*c4762a1bSJed Brown   Note that these functions use coordinates X in the global (rotated) frame
73*c4762a1bSJed Brown */
74*c4762a1bSJed Brown PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
75*c4762a1bSJed Brown {
76*c4762a1bSJed Brown   Parameter *param = (Parameter *) ctx;
77*c4762a1bSJed Brown   PetscReal  Delta = param->Delta;
78*c4762a1bSJed Brown   PetscReal  nu    = param->nu;
79*c4762a1bSJed Brown   PetscReal  u_0   = param->u_0;
80*c4762a1bSJed Brown   PetscReal  fac   = (PetscReal) (dim - 1);
81*c4762a1bSJed Brown   PetscInt   d;
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown   u[0] = u_0;
84*c4762a1bSJed Brown   for (d = 1; d < dim; ++d) u[0] += Delta/(fac * 2.0*nu) * X[d] * (1.0 - X[d]);
85*c4762a1bSJed Brown   for (d = 1; d < dim; ++d) u[d]  = 0.0;
86*c4762a1bSJed Brown   return 0;
87*c4762a1bSJed Brown }
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *p, void *ctx)
90*c4762a1bSJed Brown {
91*c4762a1bSJed Brown   Parameter *param = (Parameter *) ctx;
92*c4762a1bSJed Brown   PetscReal  Delta = param->Delta;
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   p[0] = -Delta * X[0];
95*c4762a1bSJed Brown   return 0;
96*c4762a1bSJed Brown }
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown PetscErrorCode wall_velocity(PetscInt dim, PetscReal time, const PetscReal X[], PetscInt Nf, PetscScalar *u, void *ctx)
99*c4762a1bSJed Brown {
100*c4762a1bSJed Brown   Parameter *param = (Parameter *) ctx;
101*c4762a1bSJed Brown   PetscReal  u_0   = param->u_0;
102*c4762a1bSJed Brown   PetscInt   d;
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown   u[0] = u_0;
105*c4762a1bSJed Brown   for (d = 1; d < dim; ++d) u[d] = 0.0;
106*c4762a1bSJed Brown   return 0;
107*c4762a1bSJed Brown }
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z}
110*c4762a1bSJed Brown    u[Ncomp]          = {p} */
111*c4762a1bSJed Brown void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
112*c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
113*c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
114*c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
115*c4762a1bSJed Brown {
116*c4762a1bSJed Brown   const PetscReal nu = PetscRealPart(constants[1]);
117*c4762a1bSJed Brown   const PetscInt  Nc = dim;
118*c4762a1bSJed Brown   PetscInt        c, d;
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
121*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
122*c4762a1bSJed Brown       /* f1[c*dim+d] = 0.5*nu*(u_x[c*dim+d] + u_x[d*dim+c]); */
123*c4762a1bSJed Brown       f1[c*dim+d] = nu*u_x[c*dim+d];
124*c4762a1bSJed Brown     }
125*c4762a1bSJed Brown     f1[c*dim+c] -= u[uOff[1]];
126*c4762a1bSJed Brown   }
127*c4762a1bSJed Brown }
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y, v_x, v_y} or {u_x, u_y, u_z, v_x, v_y, v_z, w_x, w_y, w_z} */
130*c4762a1bSJed Brown void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
131*c4762a1bSJed Brown           const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
132*c4762a1bSJed Brown           const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
133*c4762a1bSJed Brown           PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
134*c4762a1bSJed Brown {
135*c4762a1bSJed Brown   PetscInt d;
136*c4762a1bSJed Brown   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] += u_x[d*dim+d];
137*c4762a1bSJed Brown }
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown /* Residual functions are in reference coordinates */
140*c4762a1bSJed Brown static void f0_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
141*c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
142*c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
143*c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
144*c4762a1bSJed Brown {
145*c4762a1bSJed Brown   const PetscReal Delta = PetscRealPart(constants[0]);
146*c4762a1bSJed Brown   PetscReal       alpha = PetscRealPart(constants[3]);
147*c4762a1bSJed Brown   PetscReal       X     = PetscCosReal(alpha)*x[0] + PetscSinReal(alpha)*x[1];
148*c4762a1bSJed Brown   PetscInt        d;
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
151*c4762a1bSJed Brown     f0[d] = -Delta * X * n[d];
152*c4762a1bSJed Brown   }
153*c4762a1bSJed Brown }
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown /* < q, \nabla\cdot u >
156*c4762a1bSJed Brown    NcompI = 1, NcompJ = dim */
157*c4762a1bSJed Brown void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
158*c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
159*c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
160*c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
161*c4762a1bSJed Brown {
162*c4762a1bSJed Brown   PetscInt d;
163*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d*dim+d] = 1.0; /* \frac{\partial\phi^{u_d}}{\partial x_d} */
164*c4762a1bSJed Brown }
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown /* -< \nabla\cdot v, p >
167*c4762a1bSJed Brown     NcompI = dim, NcompJ = 1 */
168*c4762a1bSJed Brown void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
169*c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
170*c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
171*c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
172*c4762a1bSJed Brown {
173*c4762a1bSJed Brown   PetscInt d;
174*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* \frac{\partial\psi^{u_d}}{\partial x_d} */
175*c4762a1bSJed Brown }
176*c4762a1bSJed Brown 
177*c4762a1bSJed Brown /* < \nabla v, \nabla u + {\nabla u}^T >
178*c4762a1bSJed Brown    This just gives \nabla u, give the perdiagonal for the transpose */
179*c4762a1bSJed Brown void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
180*c4762a1bSJed Brown            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
181*c4762a1bSJed Brown            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
182*c4762a1bSJed Brown            PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
183*c4762a1bSJed Brown {
184*c4762a1bSJed Brown   const PetscReal nu = PetscRealPart(constants[1]);
185*c4762a1bSJed Brown   const PetscInt  Nc = dim;
186*c4762a1bSJed Brown   PetscInt        c, d;
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   for (c = 0; c < Nc; ++c) {
189*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
190*c4762a1bSJed Brown       g3[((c*Nc+c)*dim+d)*dim+d] = nu;
191*c4762a1bSJed Brown     }
192*c4762a1bSJed Brown   }
193*c4762a1bSJed Brown }
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
196*c4762a1bSJed Brown {
197*c4762a1bSJed Brown   PetscInt       n = 3;
198*c4762a1bSJed Brown   PetscErrorCode ierr;
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown   PetscFunctionBeginUser;
201*c4762a1bSJed Brown   options->dim      = 2;
202*c4762a1bSJed Brown   options->simplex  = PETSC_TRUE;
203*c4762a1bSJed Brown   options->cells[0] = 3;
204*c4762a1bSJed Brown   options->cells[1] = 3;
205*c4762a1bSJed Brown   options->cells[2] = 3;
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Poiseuille Flow Options", "DMPLEX");CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex62.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Use simplices or tensor product cells", "ex62.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-cells", "The initial mesh division", "ex62.c", options->cells, &n, NULL);CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
212*c4762a1bSJed Brown   PetscFunctionReturn(0);
213*c4762a1bSJed Brown }
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown static PetscErrorCode SetupParameters(AppCtx *user)
216*c4762a1bSJed Brown {
217*c4762a1bSJed Brown   PetscBag       bag;
218*c4762a1bSJed Brown   Parameter     *p;
219*c4762a1bSJed Brown   PetscErrorCode ierr;
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown   PetscFunctionBeginUser;
222*c4762a1bSJed Brown   /* setup PETSc parameter bag */
223*c4762a1bSJed Brown   ierr = PetscBagGetData(user->bag, (void **) &p);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = PetscBagSetName(user->bag, "par", "Poiseuille flow parameters");CHKERRQ(ierr);
225*c4762a1bSJed Brown   bag  = user->bag;
226*c4762a1bSJed Brown   ierr = PetscBagRegisterReal(bag, &p->Delta, 1.0, "Delta", "Pressure drop per unit length");CHKERRQ(ierr);
227*c4762a1bSJed Brown   ierr = PetscBagRegisterReal(bag, &p->nu,    1.0, "nu",    "Kinematic viscosity");CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = PetscBagRegisterReal(bag, &p->u_0,   0.0, "u_0",   "Tangential velocity at the wall");CHKERRQ(ierr);
229*c4762a1bSJed Brown   ierr = PetscBagRegisterReal(bag, &p->alpha, 0.0, "alpha", "Angle of pipe wall to x-axis");CHKERRQ(ierr);
230*c4762a1bSJed Brown   PetscFunctionReturn(0);
231*c4762a1bSJed Brown }
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
234*c4762a1bSJed Brown {
235*c4762a1bSJed Brown   PetscInt       dim = user->dim;
236*c4762a1bSJed Brown   PetscErrorCode ierr;
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown   PetscFunctionBeginUser;
239*c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, user->simplex, user->cells, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
240*c4762a1bSJed Brown   {
241*c4762a1bSJed Brown     Parameter   *param;
242*c4762a1bSJed Brown     Vec          coordinates;
243*c4762a1bSJed Brown     PetscScalar *coords;
244*c4762a1bSJed Brown     PetscReal    alpha;
245*c4762a1bSJed Brown     PetscInt     cdim, N, bs, i;
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown     ierr = DMGetCoordinateDim(*dm, &cdim);CHKERRQ(ierr);
248*c4762a1bSJed Brown     ierr = DMGetCoordinates(*dm, &coordinates);CHKERRQ(ierr);
249*c4762a1bSJed Brown     ierr = VecGetLocalSize(coordinates, &N);CHKERRQ(ierr);
250*c4762a1bSJed Brown     ierr = VecGetBlockSize(coordinates, &bs);CHKERRQ(ierr);
251*c4762a1bSJed Brown     if (bs != cdim) SETERRQ2(comm, PETSC_ERR_ARG_WRONG, "Invalid coordinate blocksize %D != embedding dimension %D", bs, cdim);
252*c4762a1bSJed Brown     ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
253*c4762a1bSJed Brown     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
254*c4762a1bSJed Brown     alpha = param->alpha;
255*c4762a1bSJed Brown     for (i = 0; i < N; i += cdim) {
256*c4762a1bSJed Brown       PetscScalar x = coords[i+0];
257*c4762a1bSJed Brown       PetscScalar y = coords[i+1];
258*c4762a1bSJed Brown 
259*c4762a1bSJed Brown       coords[i+0] = PetscCosReal(alpha)*x - PetscSinReal(alpha)*y;
260*c4762a1bSJed Brown       coords[i+1] = PetscSinReal(alpha)*x + PetscCosReal(alpha)*y;
261*c4762a1bSJed Brown     }
262*c4762a1bSJed Brown     ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
263*c4762a1bSJed Brown     ierr = DMSetCoordinates(*dm, coordinates);CHKERRQ(ierr);
264*c4762a1bSJed Brown   }
265*c4762a1bSJed Brown   {
266*c4762a1bSJed Brown     DM               pdm = NULL;
267*c4762a1bSJed Brown     PetscPartitioner part;
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
270*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
271*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
272*c4762a1bSJed Brown     if (pdm) {
273*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
274*c4762a1bSJed Brown       *dm  = pdm;
275*c4762a1bSJed Brown     }
276*c4762a1bSJed Brown   }
277*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
278*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
279*c4762a1bSJed Brown   PetscFunctionReturn(0);
280*c4762a1bSJed Brown }
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown PetscErrorCode SetupProblem(DM dm, AppCtx *user)
283*c4762a1bSJed Brown {
284*c4762a1bSJed Brown   PetscDS        prob;
285*c4762a1bSJed Brown   Parameter     *ctx;
286*c4762a1bSJed Brown   PetscInt       id;
287*c4762a1bSJed Brown   PetscErrorCode ierr;
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown   PetscFunctionBeginUser;
290*c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
291*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, NULL, f1_u);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 1, f0_p, NULL);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = PetscDSSetBdResidual(prob, 0, f0_bd_u, NULL);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, NULL, NULL,  NULL,  g3_uu);CHKERRQ(ierr);
295*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL,  g2_up, NULL);CHKERRQ(ierr);
296*c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_pu, NULL,  NULL);CHKERRQ(ierr);
297*c4762a1bSJed Brown   /* Setup constants */
298*c4762a1bSJed Brown   {
299*c4762a1bSJed Brown     Parameter  *param;
300*c4762a1bSJed Brown     PetscScalar constants[4];
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown     constants[0] = param->Delta;
305*c4762a1bSJed Brown     constants[1] = param->nu;
306*c4762a1bSJed Brown     constants[2] = param->u_0;
307*c4762a1bSJed Brown     constants[3] = param->alpha;
308*c4762a1bSJed Brown     ierr = PetscDSSetConstants(prob, 4, constants);CHKERRQ(ierr);
309*c4762a1bSJed Brown   }
310*c4762a1bSJed Brown   /* Setup Boundary Conditions */
311*c4762a1bSJed Brown   ierr = PetscBagGetData(user->bag, (void **) &ctx);CHKERRQ(ierr);
312*c4762a1bSJed Brown   id   = 3;
313*c4762a1bSJed Brown   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "top wall",    "marker", 0, 0, NULL, (void (*)(void)) wall_velocity, 1, &id, ctx);CHKERRQ(ierr);
314*c4762a1bSJed Brown   id   = 1;
315*c4762a1bSJed Brown   ierr = PetscDSAddBoundary(prob, DM_BC_ESSENTIAL, "bottom wall", "marker", 0, 0, NULL, (void (*)(void)) wall_velocity, 1, &id, ctx);CHKERRQ(ierr);
316*c4762a1bSJed Brown   id   = 2;
317*c4762a1bSJed Brown   ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL,   "right wall",  "marker", 0, 0, NULL, (void (*)(void)) NULL,          1, &id, ctx);CHKERRQ(ierr);
318*c4762a1bSJed Brown   /* Setup exact solution */
319*c4762a1bSJed Brown   user->exactFuncs[0] = quadratic_u;
320*c4762a1bSJed Brown   user->exactFuncs[1] = linear_p;
321*c4762a1bSJed Brown   ierr = PetscDSSetExactSolution(prob, 0, user->exactFuncs[0], ctx);CHKERRQ(ierr);
322*c4762a1bSJed Brown   ierr = PetscDSSetExactSolution(prob, 1, user->exactFuncs[1], ctx);CHKERRQ(ierr);
323*c4762a1bSJed Brown   PetscFunctionReturn(0);
324*c4762a1bSJed Brown }
325*c4762a1bSJed Brown 
326*c4762a1bSJed Brown PetscErrorCode SetupDiscretization(DM dm, AppCtx *user)
327*c4762a1bSJed Brown {
328*c4762a1bSJed Brown   DM              cdm   = dm;
329*c4762a1bSJed Brown   const PetscInt  dim   = user->dim;
330*c4762a1bSJed Brown   PetscFE         fe[2];
331*c4762a1bSJed Brown   Parameter      *param;
332*c4762a1bSJed Brown   MPI_Comm        comm;
333*c4762a1bSJed Brown   PetscErrorCode  ierr;
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown   PetscFunctionBeginUser;
336*c4762a1bSJed Brown   /* Create finite element */
337*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
338*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, dim, user->simplex, "vel_", PETSC_DEFAULT, &fe[0]);CHKERRQ(ierr);
339*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[0], "velocity");CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, user->simplex, "pres_", PETSC_DEFAULT, &fe[1]);CHKERRQ(ierr);
341*c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe[0], fe[1]);CHKERRQ(ierr);
342*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe[1], "pressure");CHKERRQ(ierr);
343*c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
344*c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
345*c4762a1bSJed Brown   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
346*c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
347*c4762a1bSJed Brown   ierr = SetupProblem(dm, user);CHKERRQ(ierr);
348*c4762a1bSJed Brown   ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
349*c4762a1bSJed Brown   while (cdm) {
350*c4762a1bSJed Brown     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
351*c4762a1bSJed Brown     ierr = DMPlexCreateBasisRotation(cdm, param->alpha, 0.0, 0.0);CHKERRQ(ierr);
352*c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
353*c4762a1bSJed Brown   }
354*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[0]);CHKERRQ(ierr);
355*c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe[1]);CHKERRQ(ierr);
356*c4762a1bSJed Brown   PetscFunctionReturn(0);
357*c4762a1bSJed Brown }
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown int main(int argc, char **argv)
360*c4762a1bSJed Brown {
361*c4762a1bSJed Brown   SNES           snes; /* nonlinear solver */
362*c4762a1bSJed Brown   DM             dm;   /* problem definition */
363*c4762a1bSJed Brown   Vec            u, r; /* solution and residual */
364*c4762a1bSJed Brown   AppCtx         user; /* user-defined work context */
365*c4762a1bSJed Brown   PetscErrorCode ierr;
366*c4762a1bSJed Brown 
367*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
368*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.bag);CHKERRQ(ierr);
370*c4762a1bSJed Brown   ierr = SetupParameters(&user);CHKERRQ(ierr);
371*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
372*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
373*c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
375*c4762a1bSJed Brown   /* Setup problem */
376*c4762a1bSJed Brown   ierr = PetscMalloc(2 * sizeof(void (*)(const PetscReal[], PetscScalar *, void *)), &user.exactFuncs);CHKERRQ(ierr);
377*c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &user);CHKERRQ(ierr);
378*c4762a1bSJed Brown   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
379*c4762a1bSJed Brown 
380*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
381*c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
382*c4762a1bSJed Brown 
383*c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm,&user,&user,&user);CHKERRQ(ierr);
384*c4762a1bSJed Brown 
385*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
386*c4762a1bSJed Brown 
387*c4762a1bSJed Brown   {
388*c4762a1bSJed Brown     Parameter *param;
389*c4762a1bSJed Brown     void      *ctxs[2];
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown     ierr = PetscBagGetData(user.bag, (void **) &param);CHKERRQ(ierr);
392*c4762a1bSJed Brown     ctxs[0] = ctxs[1] = param;
393*c4762a1bSJed Brown     ierr = DMProjectFunction(dm, 0.0, user.exactFuncs, ctxs, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
394*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) u, "Exact Solution");CHKERRQ(ierr);
395*c4762a1bSJed Brown     ierr = VecViewFromOptions(u, NULL, "-exact_vec_view");CHKERRQ(ierr);
396*c4762a1bSJed Brown   }
397*c4762a1bSJed Brown   ierr = DMSNESCheckFromOptions(snes, u, NULL, NULL);CHKERRQ(ierr);
398*c4762a1bSJed Brown   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
399*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
400*c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
401*c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
402*c4762a1bSJed Brown 
403*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
404*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
405*c4762a1bSJed Brown   ierr = PetscFree(user.exactFuncs);CHKERRQ(ierr);
406*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
407*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
408*c4762a1bSJed Brown   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
409*c4762a1bSJed Brown   ierr = PetscFinalize();
410*c4762a1bSJed Brown   return ierr;
411*c4762a1bSJed Brown }
412*c4762a1bSJed Brown 
413*c4762a1bSJed Brown /*TEST
414*c4762a1bSJed Brown 
415*c4762a1bSJed Brown   # Convergence
416*c4762a1bSJed Brown   test:
417*c4762a1bSJed Brown     suffix: 2d_quad_q1_p0_conv
418*c4762a1bSJed Brown     requires: !single
419*c4762a1bSJed Brown     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 \
420*c4762a1bSJed Brown       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
421*c4762a1bSJed Brown       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
422*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
423*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
424*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type lu \
425*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
426*c4762a1bSJed Brown   test:
427*c4762a1bSJed Brown     suffix: 2d_quad_q1_p0_conv_u0
428*c4762a1bSJed Brown     requires: !single
429*c4762a1bSJed Brown     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 \
430*c4762a1bSJed Brown       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
431*c4762a1bSJed Brown       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
432*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
433*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
434*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type lu \
435*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
436*c4762a1bSJed Brown   test:
437*c4762a1bSJed Brown     suffix: 2d_quad_q1_p0_conv_u0_alpha
438*c4762a1bSJed Brown     requires: !single
439*c4762a1bSJed Brown     args: -simplex 0 -dm_plex_separate_marker -dm_refine 1 -u_0 0.125 -alpha 0.3927 \
440*c4762a1bSJed Brown       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
441*c4762a1bSJed Brown       -snes_convergence_estimate -convest_num_refine 2 -snes_error_if_not_converged \
442*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
443*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
444*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type lu \
445*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
446*c4762a1bSJed Brown   test:
447*c4762a1bSJed Brown     suffix: 2d_quad_q1_p0_conv_gmg_vanka
448*c4762a1bSJed Brown     requires: !single long_runtime
449*c4762a1bSJed Brown     args: -simplex 0 -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \
450*c4762a1bSJed Brown       -vel_petscspace_degree 1 -pres_petscspace_degree 0 \
451*c4762a1bSJed Brown       -snes_convergence_estimate -convest_num_refine 1 -snes_error_if_not_converged \
452*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
453*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
454*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type mg \
455*c4762a1bSJed Brown           -fieldsplit_velocity_mg_levels_pc_type patch -fieldsplit_velocity_mg_levels_pc_patch_exclude_subspaces 1 \
456*c4762a1bSJed Brown           -fieldsplit_velocity_mg_levels_pc_patch_construct_codim 0 -fieldsplit_velocity_mg_levels_pc_patch_construct_type vanka \
457*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-5 -fieldsplit_pressure_pc_type jacobi
458*c4762a1bSJed Brown   test:
459*c4762a1bSJed Brown     suffix: 2d_tri_p2_p1_conv
460*c4762a1bSJed Brown     requires: triangle !single
461*c4762a1bSJed Brown     args: -dm_plex_separate_marker -dm_refine 1 \
462*c4762a1bSJed Brown       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
463*c4762a1bSJed Brown       -dmsnes_check .001 -snes_error_if_not_converged \
464*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
465*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
466*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type lu \
467*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
468*c4762a1bSJed Brown   test:
469*c4762a1bSJed Brown     suffix: 2d_tri_p2_p1_conv_u0_alpha
470*c4762a1bSJed Brown     requires: triangle !single
471*c4762a1bSJed Brown     args: -dm_plex_separate_marker -dm_refine 0 -u_0 0.125 -alpha 0.3927 \
472*c4762a1bSJed Brown       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
473*c4762a1bSJed Brown       -dmsnes_check .001 -snes_error_if_not_converged \
474*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
475*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
476*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type lu \
477*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
478*c4762a1bSJed Brown   test:
479*c4762a1bSJed Brown     suffix: 2d_tri_p2_p1_conv_gmg_vcycle
480*c4762a1bSJed Brown     requires: triangle !single
481*c4762a1bSJed Brown     args: -dm_plex_separate_marker -cells 2,2 -dm_refine_hierarchy 1 \
482*c4762a1bSJed Brown       -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
483*c4762a1bSJed Brown       -dmsnes_check .001 -snes_error_if_not_converged \
484*c4762a1bSJed Brown       -ksp_type fgmres -ksp_gmres_restart 10 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
485*c4762a1bSJed Brown       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
486*c4762a1bSJed Brown         -fieldsplit_velocity_pc_type mg \
487*c4762a1bSJed Brown         -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
488*c4762a1bSJed Brown TEST*/
489