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