xref: /petsc/src/snes/tutorials/ex62.c (revision 8b0e23d0053fb4bd1ea930debb90f0b8b83d8d09)
1*8b0e23d0SMatthew G. Knepley static char help[] = "Stokes Problem discretized with finite elements,\n\
2*8b0e23d0SMatthew G. Knepley using a parallel unstructured mesh (DMPLEX) to represent the domain.\n\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5*8b0e23d0SMatthew G. Knepley For the isoviscous Stokes problem, which we discretize using the finite
6*8b0e23d0SMatthew G. Knepley element method on an unstructured mesh, the weak form equations are
7c4762a1bSJed Brown 
8*8b0e23d0SMatthew G. Knepley   < \nabla v, \nabla u + {\nabla u}^T > - < \nabla\cdot v, p > - < v, f > = 0
9*8b0e23d0SMatthew G. Knepley   < q, -\nabla\cdot u >                                                   = 0
10c4762a1bSJed Brown 
11c4762a1bSJed Brown Viewing:
12c4762a1bSJed Brown 
13c4762a1bSJed Brown To produce nice output, use
14c4762a1bSJed Brown 
15*8b0e23d0SMatthew G. Knepley   -dm_refine 3 -dm_view hdf5:sol1.h5 -error_vec_view hdf5:sol1.h5::append -snes_view_solution hdf5:sol1.h5::append -exact_vec_view hdf5:sol1.h5::append
16c4762a1bSJed Brown 
17c4762a1bSJed Brown You can get a LaTeX view of the mesh, with point numbering using
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   -dm_view :mesh.tex:ascii_latex -dm_plex_view_scale 8.0
20c4762a1bSJed Brown 
21c4762a1bSJed Brown The data layout can be viewed using
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   -dm_petscsection_view
24c4762a1bSJed Brown 
25c4762a1bSJed Brown Lots of information about the FEM assembly can be printed using
26c4762a1bSJed Brown 
27*8b0e23d0SMatthew G. Knepley   -dm_plex_print_fem 3
28c4762a1bSJed Brown */
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #include <petscdmplex.h>
31c4762a1bSJed Brown #include <petscsnes.h>
32c4762a1bSJed Brown #include <petscds.h>
33*8b0e23d0SMatthew G. Knepley #include <petscbag.h>
34c4762a1bSJed Brown 
35*8b0e23d0SMatthew G. Knepley // TODO: Plot residual by fields after each smoother iterate
36c4762a1bSJed Brown 
37*8b0e23d0SMatthew G. Knepley typedef enum {SOL_QUADRATIC, SOL_TRIG, SOL_UNKNOWN} SolType;
38*8b0e23d0SMatthew G. Knepley const char *SolTypes[] = {"quadratic", "trig", "unknown", "SolType", "SOL_", 0};
39c4762a1bSJed Brown 
40c4762a1bSJed Brown typedef struct {
41*8b0e23d0SMatthew G. Knepley   PetscScalar mu; /* dynamic shear viscosity */
42*8b0e23d0SMatthew G. Knepley } Parameter;
43*8b0e23d0SMatthew G. Knepley 
44*8b0e23d0SMatthew G. Knepley typedef struct {
45c4762a1bSJed Brown   /* Domain and mesh definition */
46*8b0e23d0SMatthew G. Knepley   char     filename[PETSC_MAX_PATH_LEN];
47c4762a1bSJed Brown   /* Problem definition */
48*8b0e23d0SMatthew G. Knepley   PetscBag bag; /* Problem parameters */
49*8b0e23d0SMatthew G. Knepley   SolType  sol; /* MMS solution */
50c4762a1bSJed Brown } AppCtx;
51c4762a1bSJed Brown 
52*8b0e23d0SMatthew G. Knepley static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
53c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
54c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
55c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
56c4762a1bSJed Brown {
57*8b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
58*8b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1]-uOff[0];
59*8b0e23d0SMatthew G. Knepley   PetscInt        c, d;
60c4762a1bSJed Brown 
61*8b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
62c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
63*8b0e23d0SMatthew G. Knepley       f1[c*dim+d] = mu * (u_x[c*dim+d] + u_x[d*dim+c]);
64c4762a1bSJed Brown     }
65*8b0e23d0SMatthew G. Knepley     f1[c*dim+c] -= u[uOff[1]];
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69*8b0e23d0SMatthew G. Knepley static void f0_p(PetscInt dim, PetscInt Nf, PetscInt NfAux,
70c4762a1bSJed Brown                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
71c4762a1bSJed Brown                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
72c4762a1bSJed Brown                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
73c4762a1bSJed Brown {
74c4762a1bSJed Brown   PetscInt d;
75*8b0e23d0SMatthew G. Knepley   for (d = 0, f0[0] = 0.0; d < dim; ++d) f0[0] -= u_x[d*dim+d];
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
78*8b0e23d0SMatthew G. Knepley static void g1_pu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
79c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
80c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
81c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   PetscInt d;
84*8b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0; /* < q, -\nabla\cdot u > */
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87*8b0e23d0SMatthew G. Knepley static void g2_up(PetscInt dim, PetscInt Nf, PetscInt NfAux,
88c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
89c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
90c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
91c4762a1bSJed Brown {
92c4762a1bSJed Brown   PetscInt d;
93*8b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) g2[d*dim+d] = -1.0; /* -< \nabla\cdot v, p > */
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96*8b0e23d0SMatthew G. Knepley static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97c4762a1bSJed Brown                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
98c4762a1bSJed Brown                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
99c4762a1bSJed Brown                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
100c4762a1bSJed Brown {
101*8b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
102*8b0e23d0SMatthew G. Knepley   const PetscInt  Nc = uOff[1]-uOff[0];
103*8b0e23d0SMatthew G. Knepley   PetscInt        c, d;
104c4762a1bSJed Brown 
105*8b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) {
106c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {
107*8b0e23d0SMatthew G. Knepley       g3[((c*Nc+c)*dim+d)*dim+d] += mu; /* < \nabla v, \nabla u > */
108*8b0e23d0SMatthew G. Knepley       g3[((c*Nc+d)*dim+d)*dim+c] += mu; /* < \nabla v, {\nabla u}^T > */
109c4762a1bSJed Brown     }
110c4762a1bSJed Brown   }
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113*8b0e23d0SMatthew G. Knepley static void g0_pp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
114*8b0e23d0SMatthew G. Knepley                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
115*8b0e23d0SMatthew G. Knepley                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
116*8b0e23d0SMatthew G. Knepley                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
117*8b0e23d0SMatthew G. Knepley {
118*8b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
119*8b0e23d0SMatthew G. Knepley 
120*8b0e23d0SMatthew G. Knepley   g0[0] = 1.0/mu;
121*8b0e23d0SMatthew G. Knepley }
122*8b0e23d0SMatthew G. Knepley 
123*8b0e23d0SMatthew G. Knepley /* Quadratic MMS Solution
124*8b0e23d0SMatthew G. Knepley    2D:
125c4762a1bSJed Brown 
126c4762a1bSJed Brown      u = x^2 + y^2
127*8b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
128*8b0e23d0SMatthew G. Knepley      p = x + y - 1
129*8b0e23d0SMatthew G. Knepley      f = <1 - 4 mu, 1 - 4 mu>
130c4762a1bSJed Brown 
131c4762a1bSJed Brown    so that
132c4762a1bSJed Brown 
133*8b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 4x  4x \
134*8b0e23d0SMatthew G. Knepley                                   \ 4x -4x /
135*8b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <4, 4> - <1, 1> + <1 - 4 mu, 1 - 4 mu> = 0
136*8b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2x - 2x = 0
137*8b0e23d0SMatthew G. Knepley 
138*8b0e23d0SMatthew G. Knepley    3D:
139*8b0e23d0SMatthew G. Knepley 
140*8b0e23d0SMatthew G. Knepley      u = 2 x^2 + y^2 + z^2
141*8b0e23d0SMatthew G. Knepley      v = 2 x^2 - 2xy
142*8b0e23d0SMatthew G. Knepley      w = 2 x^2 - 2xz
143*8b0e23d0SMatthew G. Knepley      p = x + y + z - 3/2
144*8b0e23d0SMatthew G. Knepley      f = <1 - 8 mu, 1 - 4 mu, 1 - 4 mu>
145*8b0e23d0SMatthew G. Knepley 
146*8b0e23d0SMatthew G. Knepley    so that
147*8b0e23d0SMatthew G. Knepley 
148*8b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = / 8x  4x  4x \
149*8b0e23d0SMatthew G. Knepley                                   | 4x -4x  0  |
150*8b0e23d0SMatthew G. Knepley                                   \ 4x  0  -4x /
151*8b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <8, 4, 4> - <1, 1, 1> + <1 - 8 mu, 1 - 4 mu, 1 - 4 mu> = 0
152*8b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 4x - 2x - 2x = 0
153c4762a1bSJed Brown */
154*8b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
155c4762a1bSJed Brown {
156*8b0e23d0SMatthew G. Knepley   PetscInt c;
157*8b0e23d0SMatthew G. Knepley 
158*8b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSqr(x[0]);
159*8b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
160*8b0e23d0SMatthew G. Knepley     u[0] += PetscSqr(x[c]);
161*8b0e23d0SMatthew G. Knepley     u[c]  = 2.0*PetscSqr(x[0]) - 2.0*x[0]*x[c];
162*8b0e23d0SMatthew G. Knepley   }
163c4762a1bSJed Brown   return 0;
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166*8b0e23d0SMatthew G. Knepley static PetscErrorCode quadratic_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
167c4762a1bSJed Brown {
168*8b0e23d0SMatthew G. Knepley   PetscInt d;
169*8b0e23d0SMatthew G. Knepley 
170*8b0e23d0SMatthew G. Knepley   u[0] = -0.5*dim;
171*8b0e23d0SMatthew G. Knepley   for (d = 0; d < dim; ++d) u[0] += x[d];
172c4762a1bSJed Brown   return 0;
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175*8b0e23d0SMatthew G. Knepley static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176c4762a1bSJed Brown                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177c4762a1bSJed Brown                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178*8b0e23d0SMatthew G. Knepley                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
179c4762a1bSJed Brown {
180*8b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
181*8b0e23d0SMatthew G. Knepley   PetscInt        d;
182*8b0e23d0SMatthew G. Knepley 
183*8b0e23d0SMatthew G. Knepley   f0[0] = (dim-1)*4.0*mu - 1.0;
184*8b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) f0[d] = 4.0*mu - 1.0;
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187*8b0e23d0SMatthew G. Knepley /* Trigonometric MMS Solution
188*8b0e23d0SMatthew G. Knepley    2D:
189*8b0e23d0SMatthew G. Knepley 
190*8b0e23d0SMatthew G. Knepley      u = sin(pi x) + sin(pi y)
191*8b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
192*8b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y)
193*8b0e23d0SMatthew G. Knepley      f = <2pi cos(2 pi x) + mu pi^2 sin(pi x) + mu pi^2 sin(pi y), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y>
194*8b0e23d0SMatthew G. Knepley 
195*8b0e23d0SMatthew G. Knepley    so that
196*8b0e23d0SMatthew G. Knepley 
197*8b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = /        2pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y \
198*8b0e23d0SMatthew G. Knepley                                   \ pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)          /
199*8b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <-pi^2 sin(pi x) - pi^2 sin(pi y), pi^3 cos(pi x) y> - <2pi cos(2 pi x), 2pi cos(2 pi y)> + <f_x, f_y> = 0
200*8b0e23d0SMatthew G. Knepley      \nabla \cdot u             = pi cos(pi x) - pi cos(pi x) = 0
201*8b0e23d0SMatthew G. Knepley 
202*8b0e23d0SMatthew G. Knepley    3D:
203*8b0e23d0SMatthew G. Knepley 
204*8b0e23d0SMatthew G. Knepley      u = 2 sin(pi x) + sin(pi y) + sin(pi z)
205*8b0e23d0SMatthew G. Knepley      v = -pi cos(pi x) y
206*8b0e23d0SMatthew G. Knepley      w = -pi cos(pi x) z
207*8b0e23d0SMatthew G. Knepley      p = sin(2 pi x) + sin(2 pi y) + sin(2 pi z)
208*8b0e23d0SMatthew G. Knepley      f = <2pi cos(2 pi x) + mu 2pi^2 sin(pi x) + mu pi^2 sin(pi y) + mu pi^2 sin(pi z), 2pi cos(2 pi y) - mu pi^3 cos(pi x) y, 2pi cos(2 pi z) - mu pi^3 cos(pi x) z>
209*8b0e23d0SMatthew G. Knepley 
210*8b0e23d0SMatthew G. Knepley    so that
211*8b0e23d0SMatthew G. Knepley 
212*8b0e23d0SMatthew G. Knepley      e(u) = (grad u + grad u^T) = /        4pi cos(pi x)             pi cos(pi y) + pi^2 sin(pi x) y  pi cos(pi z) + pi^2 sin(pi x) z \
213*8b0e23d0SMatthew G. Knepley                                   | pi cos(pi y) + pi^2 sin(pi x) y          -2pi cos(pi x)                        0                  |
214*8b0e23d0SMatthew G. Knepley                                   \ pi cos(pi z) + pi^2 sin(pi x) z               0                         -2pi cos(pi x)            /
215*8b0e23d0SMatthew G. Knepley      div mu e(u) - \nabla p + f = mu <-2pi^2 sin(pi x) - pi^2 sin(pi y) - pi^2 sin(pi z), pi^3 cos(pi x) y, pi^3 cos(pi x) z> - <2pi cos(2 pi x), 2pi cos(2 pi y), 2pi cos(2 pi z)> + <f_x, f_y, f_z> = 0
216*8b0e23d0SMatthew G. Knepley      \nabla \cdot u             = 2 pi cos(pi x) - pi cos(pi x) - pi cos(pi x) = 0
217*8b0e23d0SMatthew G. Knepley */
218*8b0e23d0SMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
219c4762a1bSJed Brown {
220*8b0e23d0SMatthew G. Knepley   PetscInt c;
221*8b0e23d0SMatthew G. Knepley 
222*8b0e23d0SMatthew G. Knepley   u[0] = (dim-1)*PetscSinReal(PETSC_PI*x[0]);
223*8b0e23d0SMatthew G. Knepley   for (c = 1; c < Nc; ++c) {
224*8b0e23d0SMatthew G. Knepley     u[0] += PetscSinReal(PETSC_PI*x[c]);
225*8b0e23d0SMatthew G. Knepley     u[c]  = -PETSC_PI*PetscCosReal(PETSC_PI*x[0]) * x[c];
226*8b0e23d0SMatthew G. Knepley   }
227*8b0e23d0SMatthew G. Knepley   return 0;
228*8b0e23d0SMatthew G. Knepley }
229*8b0e23d0SMatthew G. Knepley 
230*8b0e23d0SMatthew G. Knepley static PetscErrorCode trig_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
231*8b0e23d0SMatthew G. Knepley {
232*8b0e23d0SMatthew G. Knepley   PetscInt d;
233*8b0e23d0SMatthew G. Knepley 
234*8b0e23d0SMatthew G. Knepley   for (d = 0, u[0] = 0.0; d < dim; ++d) u[0] += PetscSinReal(2.0*PETSC_PI*x[d]);
235*8b0e23d0SMatthew G. Knepley   return 0;
236*8b0e23d0SMatthew G. Knepley }
237*8b0e23d0SMatthew G. Knepley 
238*8b0e23d0SMatthew G. Knepley static void f0_trig_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
239*8b0e23d0SMatthew G. Knepley                       const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
240*8b0e23d0SMatthew G. Knepley                       const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
241*8b0e23d0SMatthew G. Knepley                       PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
242*8b0e23d0SMatthew G. Knepley {
243*8b0e23d0SMatthew G. Knepley   const PetscReal mu = PetscRealPart(constants[0]);
244*8b0e23d0SMatthew G. Knepley   PetscInt        d;
245*8b0e23d0SMatthew G. Knepley 
246*8b0e23d0SMatthew G. Knepley   f0[0] = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[0]) - (dim-1)*mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[0]);
247*8b0e23d0SMatthew G. Knepley   for (d = 1; d < dim; ++d) {
248*8b0e23d0SMatthew G. Knepley     f0[0] -= mu*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x[d]);
249*8b0e23d0SMatthew G. Knepley     f0[d]  = -2.0*PETSC_PI*PetscCosReal(2.0*PETSC_PI*x[d]) + mu*PetscPowRealInt(PETSC_PI, 3)*PetscCosReal(PETSC_PI*x[0])*x[d];
250*8b0e23d0SMatthew G. Knepley   }
251*8b0e23d0SMatthew G. Knepley }
252*8b0e23d0SMatthew G. Knepley 
253*8b0e23d0SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
254*8b0e23d0SMatthew G. Knepley {
255*8b0e23d0SMatthew G. Knepley   PetscInt       sol;
256c4762a1bSJed Brown   PetscErrorCode ierr;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBeginUser;
259*8b0e23d0SMatthew G. Knepley   options->filename[0] = '\0';
260*8b0e23d0SMatthew G. Knepley   options->sol         = SOL_QUADRATIC;
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Stokes Problem Options", "DMPLEX");CHKERRQ(ierr);
263*8b0e23d0SMatthew G. Knepley   sol  = options->sol;
264*8b0e23d0SMatthew G. Knepley   ierr = PetscOptionsEList("-sol", "The MMS solution", "ex64.c", SolTypes, (sizeof(SolTypes)/sizeof(SolTypes[0]))-3, SolTypes[options->sol], &sol, NULL);CHKERRQ(ierr);
265*8b0e23d0SMatthew G. Knepley   options->sol = (SolType) sol;
266c4762a1bSJed Brown   ierr = PetscOptionsEnd();
267c4762a1bSJed Brown   PetscFunctionReturn(0);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
270*8b0e23d0SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
271c4762a1bSJed Brown {
272*8b0e23d0SMatthew G. Knepley   size_t         len;
273c4762a1bSJed Brown   PetscErrorCode ierr;
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   PetscFunctionBeginUser;
276*8b0e23d0SMatthew G. Knepley   ierr = PetscStrlen(user->filename, &len);CHKERRQ(ierr);
277*8b0e23d0SMatthew G. Knepley   if (len) {ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr);}
278*8b0e23d0SMatthew G. Knepley   else     {ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);}
279c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
281c4762a1bSJed Brown   PetscFunctionReturn(0);
282c4762a1bSJed Brown }
283c4762a1bSJed Brown 
284*8b0e23d0SMatthew G. Knepley static PetscErrorCode SetupParameters(MPI_Comm comm, AppCtx *ctx)
285c4762a1bSJed Brown {
286*8b0e23d0SMatthew G. Knepley   Parameter     *p;
287*8b0e23d0SMatthew G. Knepley   PetscErrorCode ierr;
288*8b0e23d0SMatthew G. Knepley 
289*8b0e23d0SMatthew G. Knepley   PetscFunctionBeginUser;
290*8b0e23d0SMatthew G. Knepley   /* setup PETSc parameter bag */
291*8b0e23d0SMatthew G. Knepley   ierr = PetscBagCreate(PETSC_COMM_SELF, sizeof(Parameter), &ctx->bag);CHKERRQ(ierr);
292*8b0e23d0SMatthew G. Knepley   ierr = PetscBagGetData(ctx->bag, (void **) &p);CHKERRQ(ierr);
293*8b0e23d0SMatthew G. Knepley   ierr = PetscBagSetName(ctx->bag, "par", "Stokes Parameters");CHKERRQ(ierr);
294*8b0e23d0SMatthew G. Knepley   ierr = PetscBagRegisterScalar(ctx->bag, &p->mu, 1.0, "mu", "Dynamic Shear Viscosity, Pa s");CHKERRQ(ierr);
295*8b0e23d0SMatthew G. Knepley   ierr = PetscBagSetFromOptions(ctx->bag);CHKERRQ(ierr);
296*8b0e23d0SMatthew G. Knepley   {
297*8b0e23d0SMatthew G. Knepley     PetscViewer       viewer;
298*8b0e23d0SMatthew G. Knepley     PetscViewerFormat format;
299*8b0e23d0SMatthew G. Knepley     PetscBool         flg;
300*8b0e23d0SMatthew G. Knepley 
301*8b0e23d0SMatthew G. Knepley     ierr = PetscOptionsGetViewer(comm, NULL, NULL, "-param_view", &viewer, &format, &flg);CHKERRQ(ierr);
302*8b0e23d0SMatthew G. Knepley     if (flg) {
303*8b0e23d0SMatthew G. Knepley       ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
304*8b0e23d0SMatthew G. Knepley       ierr = PetscBagView(ctx->bag, viewer);CHKERRQ(ierr);
305*8b0e23d0SMatthew G. Knepley       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
306*8b0e23d0SMatthew G. Knepley       ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
307*8b0e23d0SMatthew G. Knepley       ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
308*8b0e23d0SMatthew G. Knepley     }
309*8b0e23d0SMatthew G. Knepley   }
310*8b0e23d0SMatthew G. Knepley   PetscFunctionReturn(0);
311*8b0e23d0SMatthew G. Knepley }
312*8b0e23d0SMatthew G. Knepley 
313*8b0e23d0SMatthew G. Knepley static PetscErrorCode SetupEqn(DM dm, AppCtx *user)
314*8b0e23d0SMatthew G. Knepley {
315*8b0e23d0SMatthew G. Knepley   PetscErrorCode (*exactFuncs[2])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *);
316*8b0e23d0SMatthew G. Knepley   PetscDS          ds;
317c4762a1bSJed Brown   const PetscInt   id = 1;
318c4762a1bSJed Brown   PetscErrorCode   ierr;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   PetscFunctionBeginUser;
321*8b0e23d0SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
322*8b0e23d0SMatthew G. Knepley   switch (user->sol) {
323c4762a1bSJed Brown     case SOL_QUADRATIC:
324*8b0e23d0SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_quadratic_u, f1_u);CHKERRQ(ierr);
325*8b0e23d0SMatthew G. Knepley       exactFuncs[0] = quadratic_u;
326*8b0e23d0SMatthew G. Knepley       exactFuncs[1] = quadratic_p;
327c4762a1bSJed Brown       break;
328c4762a1bSJed Brown     case SOL_TRIG:
329*8b0e23d0SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_trig_u, f1_u);CHKERRQ(ierr);
330*8b0e23d0SMatthew G. Knepley       exactFuncs[0] = trig_u;
331*8b0e23d0SMatthew G. Knepley       exactFuncs[1] = trig_p;
332c4762a1bSJed Brown       break;
333*8b0e23d0SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unsupported solution type: %s (%D)", SolTypes[PetscMin(user->sol, SOL_UNKNOWN)], user->sol);
334c4762a1bSJed Brown   }
335*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetResidual(ds, 1, f0_p, NULL);CHKERRQ(ierr);
336*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, NULL, NULL,  NULL,  g3_uu);CHKERRQ(ierr);
337*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 1, NULL, NULL,  g2_up, NULL);CHKERRQ(ierr);
338*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 1, 0, NULL, g1_pu, NULL,  NULL);CHKERRQ(ierr);
339*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobianPreconditioner(ds, 0, 0, NULL, NULL, NULL, g3_uu);CHKERRQ(ierr);
340*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetJacobianPreconditioner(ds, 1, 1, g0_pp, NULL, NULL, NULL);CHKERRQ(ierr);
341c4762a1bSJed Brown 
342*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 0, exactFuncs[0], user);CHKERRQ(ierr);
343*8b0e23d0SMatthew G. Knepley   ierr = PetscDSSetExactSolution(ds, 1, exactFuncs[1], user);CHKERRQ(ierr);
344c4762a1bSJed Brown 
345*8b0e23d0SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) exactFuncs[0], NULL, 1, &id, user);CHKERRQ(ierr);
346*8b0e23d0SMatthew G. Knepley 
347c4762a1bSJed Brown   {
348*8b0e23d0SMatthew G. Knepley     Parameter  *param;
349*8b0e23d0SMatthew G. Knepley     PetscScalar constants[1];
350c4762a1bSJed Brown 
351*8b0e23d0SMatthew G. Knepley     ierr = PetscBagGetData(user->bag, (void **) &param);CHKERRQ(ierr);
352*8b0e23d0SMatthew G. Knepley     constants[0] = param->mu; /* dynamic shear viscosity, Pa s */
353*8b0e23d0SMatthew G. Knepley     ierr = PetscDSSetConstants(ds, 1, constants);CHKERRQ(ierr);
354c4762a1bSJed Brown   }
355c4762a1bSJed Brown   PetscFunctionReturn(0);
356c4762a1bSJed Brown }
357c4762a1bSJed Brown 
358*8b0e23d0SMatthew G. Knepley static PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
359*8b0e23d0SMatthew G. Knepley {
360*8b0e23d0SMatthew G. Knepley   PetscInt c;
361*8b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 0.0;
362*8b0e23d0SMatthew G. Knepley   return 0;
363*8b0e23d0SMatthew G. Knepley }
364*8b0e23d0SMatthew G. Knepley static PetscErrorCode one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
365*8b0e23d0SMatthew G. Knepley {
366*8b0e23d0SMatthew G. Knepley   PetscInt c;
367*8b0e23d0SMatthew G. Knepley   for (c = 0; c < Nc; ++c) u[c] = 1.0;
368*8b0e23d0SMatthew G. Knepley   return 0;
369*8b0e23d0SMatthew G. Knepley }
370*8b0e23d0SMatthew G. Knepley 
371*8b0e23d0SMatthew G. Knepley static PetscErrorCode CreatePressureNullSpace(DM dm, PetscInt origField, PetscInt field, MatNullSpace *nullspace)
372c4762a1bSJed Brown {
373c4762a1bSJed Brown   Vec              vec;
374*8b0e23d0SMatthew G. Knepley   PetscErrorCode (*funcs[2])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void* ctx) = {zero, zero};
375c4762a1bSJed Brown   PetscErrorCode   ierr;
376c4762a1bSJed Brown 
377c4762a1bSJed Brown   PetscFunctionBeginUser;
378*8b0e23d0SMatthew G. Knepley   if (origField != 1) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Field %D should be 1 for pressure", origField);
379*8b0e23d0SMatthew G. Knepley   funcs[field] = one;
380*8b0e23d0SMatthew G. Knepley   {
381*8b0e23d0SMatthew G. Knepley     PetscDS ds;
382*8b0e23d0SMatthew G. Knepley     ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
383*8b0e23d0SMatthew G. Knepley     ierr = PetscObjectViewFromOptions((PetscObject) ds, NULL, "-ds_view");CHKERRQ(ierr);
384*8b0e23d0SMatthew G. Knepley   }
385c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &vec);CHKERRQ(ierr);
386c4762a1bSJed Brown   ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_ALL_VALUES, vec);CHKERRQ(ierr);
387c4762a1bSJed Brown   ierr = VecNormalize(vec, NULL);CHKERRQ(ierr);
388c4762a1bSJed Brown   ierr = MatNullSpaceCreate(PetscObjectComm((PetscObject)dm), PETSC_FALSE, 1, &vec, nullspace);CHKERRQ(ierr);
389c4762a1bSJed Brown   ierr = VecDestroy(&vec);CHKERRQ(ierr);
390c4762a1bSJed Brown   /* New style for field null spaces */
391c4762a1bSJed Brown   {
392c4762a1bSJed Brown     PetscObject  pressure;
393c4762a1bSJed Brown     MatNullSpace nullspacePres;
394c4762a1bSJed Brown 
395*8b0e23d0SMatthew G. Knepley     ierr = DMGetField(dm, field, NULL, &pressure);CHKERRQ(ierr);
396c4762a1bSJed Brown     ierr = MatNullSpaceCreate(PetscObjectComm(pressure), PETSC_TRUE, 0, NULL, &nullspacePres);CHKERRQ(ierr);
397c4762a1bSJed Brown     ierr = PetscObjectCompose(pressure, "nullspace", (PetscObject) nullspacePres);CHKERRQ(ierr);
398c4762a1bSJed Brown     ierr = MatNullSpaceDestroy(&nullspacePres);CHKERRQ(ierr);
399c4762a1bSJed Brown   }
400c4762a1bSJed Brown   PetscFunctionReturn(0);
401c4762a1bSJed Brown }
402c4762a1bSJed Brown 
403*8b0e23d0SMatthew G. Knepley static PetscErrorCode SetupProblem(DM dm, PetscErrorCode (*setupEqn)(DM, AppCtx *), AppCtx *user)
404c4762a1bSJed Brown {
405*8b0e23d0SMatthew G. Knepley   DM              cdm = dm;
406*8b0e23d0SMatthew G. Knepley   PetscQuadrature q   = NULL;
407*8b0e23d0SMatthew G. Knepley   DMPolytopeType  ct;
408*8b0e23d0SMatthew G. Knepley   PetscBool       simplex;
409*8b0e23d0SMatthew G. Knepley   PetscInt        dim, Nf = 2, f, Nc[2], cStart;
410*8b0e23d0SMatthew G. Knepley   const char     *name[2]   = {"velocity", "pressure"};
411*8b0e23d0SMatthew G. Knepley   const char     *prefix[2] = {"vel_",     "pres_"};
412c4762a1bSJed Brown   PetscErrorCode  ierr;
413c4762a1bSJed Brown 
414*8b0e23d0SMatthew G. Knepley   PetscFunctionBegin;
415*8b0e23d0SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
416*8b0e23d0SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
417*8b0e23d0SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
418*8b0e23d0SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
419*8b0e23d0SMatthew G. Knepley   Nc[0] = dim;
420*8b0e23d0SMatthew G. Knepley   Nc[1] = 1;
421*8b0e23d0SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
422*8b0e23d0SMatthew G. Knepley     PetscFE fe;
423*8b0e23d0SMatthew G. Knepley 
424*8b0e23d0SMatthew G. Knepley     ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc[f], simplex, prefix[f], -1, &fe);CHKERRQ(ierr);
425*8b0e23d0SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) fe, name[f]);CHKERRQ(ierr);
426*8b0e23d0SMatthew G. Knepley     if (!q) {ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);}
427*8b0e23d0SMatthew G. Knepley     ierr = PetscFESetQuadrature(fe, q);CHKERRQ(ierr);
428*8b0e23d0SMatthew G. Knepley     ierr = DMSetField(dm, f, NULL, (PetscObject) fe);CHKERRQ(ierr);
429*8b0e23d0SMatthew G. Knepley     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
430c4762a1bSJed Brown   }
431*8b0e23d0SMatthew G. Knepley   ierr = DMCreateDS(dm);CHKERRQ(ierr);
432*8b0e23d0SMatthew G. Knepley   ierr = (*setupEqn)(dm, user);CHKERRQ(ierr);
433*8b0e23d0SMatthew G. Knepley   while (cdm) {
434*8b0e23d0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
435*8b0e23d0SMatthew G. Knepley     ierr = DMSetNullSpaceConstructor(cdm, 1, CreatePressureNullSpace);CHKERRQ(ierr);
436*8b0e23d0SMatthew G. Knepley     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
437c4762a1bSJed Brown   }
438c4762a1bSJed Brown   PetscFunctionReturn(0);
439c4762a1bSJed Brown }
440c4762a1bSJed Brown 
441c4762a1bSJed Brown int main(int argc, char **argv)
442c4762a1bSJed Brown {
443*8b0e23d0SMatthew G. Knepley   SNES           snes;
444*8b0e23d0SMatthew G. Knepley   DM             dm;
445*8b0e23d0SMatthew G. Knepley   Vec            u;
446*8b0e23d0SMatthew G. Knepley   AppCtx         user;
447c4762a1bSJed Brown   PetscErrorCode ierr;
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
450c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
451c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
452*8b0e23d0SMatthew G. Knepley   ierr = SNESCreate(PetscObjectComm((PetscObject) dm), &snes);CHKERRQ(ierr);
453c4762a1bSJed Brown   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
454c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &user);CHKERRQ(ierr);
455c4762a1bSJed Brown 
456*8b0e23d0SMatthew G. Knepley   ierr = SetupParameters(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
457*8b0e23d0SMatthew G. Knepley   ierr = SetupProblem(dm, SetupEqn, &user);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = DMPlexCreateClosureIndex(dm, NULL);CHKERRQ(ierr);
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
461c4762a1bSJed Brown   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
462c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
463*8b0e23d0SMatthew G. Knepley   ierr = DMSNESCheckFromOptions(snes, u);CHKERRQ(ierr);
464c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "Solution");CHKERRQ(ierr);
465*8b0e23d0SMatthew G. Knepley   {
466*8b0e23d0SMatthew G. Knepley     Mat          J;
467*8b0e23d0SMatthew G. Knepley     MatNullSpace sp;
468c4762a1bSJed Brown 
469*8b0e23d0SMatthew G. Knepley     ierr = SNESSetUp(snes);CHKERRQ(ierr);
470*8b0e23d0SMatthew G. Knepley     ierr = CreatePressureNullSpace(dm, 1, 1, &sp);CHKERRQ(ierr);
471*8b0e23d0SMatthew G. Knepley     ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr);
472*8b0e23d0SMatthew G. Knepley     ierr = MatSetNullSpace(J, sp);CHKERRQ(ierr);
473*8b0e23d0SMatthew G. Knepley     ierr = MatNullSpaceDestroy(&sp);CHKERRQ(ierr);
474*8b0e23d0SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) J, "Jacobian");CHKERRQ(ierr);
475*8b0e23d0SMatthew G. Knepley     ierr = MatViewFromOptions(J, NULL, "-J_view");CHKERRQ(ierr);
476c4762a1bSJed Brown   }
477c4762a1bSJed Brown   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
481c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
482*8b0e23d0SMatthew G. Knepley   ierr = PetscBagDestroy(&user.bag);CHKERRQ(ierr);
483c4762a1bSJed Brown   ierr = PetscFinalize();
484c4762a1bSJed Brown   return ierr;
485c4762a1bSJed Brown }
486c4762a1bSJed Brown /*TEST
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   test:
489*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check
490c4762a1bSJed Brown     requires: triangle
491*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
492*8b0e23d0SMatthew G. Knepley 
493c4762a1bSJed Brown   test:
494*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_check_parallel
495*8b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
496c4762a1bSJed Brown     requires: triangle
497*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
498*8b0e23d0SMatthew G. Knepley 
499c4762a1bSJed Brown   test:
500*8b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check
501c4762a1bSJed Brown     requires: ctetgen
502*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
503*8b0e23d0SMatthew G. Knepley 
504c4762a1bSJed Brown   test:
505*8b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_check_parallel
506*8b0e23d0SMatthew G. Knepley     nsize: {{2 3 5}}
507c4762a1bSJed Brown     requires: ctetgen
508*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
509*8b0e23d0SMatthew G. Knepley 
510c4762a1bSJed Brown   test:
511*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv
512*8b0e23d0SMatthew G. Knepley     requires: triangle
513*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [3.0, 2.1]
514*8b0e23d0SMatthew G. Knepley     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
515*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
516*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
517*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
518*8b0e23d0SMatthew G. Knepley 
519c4762a1bSJed Brown   test:
520*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_conv_gamg
521*8b0e23d0SMatthew G. Knepley     requires: triangle
522*8b0e23d0SMatthew G. Knepley     args: -sol trig -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
523*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
524*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
525*8b0e23d0SMatthew G. Knepley 
526c4762a1bSJed Brown   test:
527*8b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_conv
528*8b0e23d0SMatthew G. Knepley     requires: ctetgen !single
529*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.8]
530*8b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
531*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
532*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
533*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
534*8b0e23d0SMatthew G. Knepley 
535c4762a1bSJed Brown   test:
536*8b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_check
537*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
538*8b0e23d0SMatthew G. Knepley 
539c4762a1bSJed Brown   test:
540*8b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_check
541*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001
542*8b0e23d0SMatthew G. Knepley 
543c4762a1bSJed Brown   test:
544*8b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_conv
545*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 -convest_num_refine 1 gives L_2 convergence rate: [3.0, 2.1]
546*8b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_simplex 0 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 -ksp_error_if_not_converged \
547*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
548*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
549*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
550*8b0e23d0SMatthew G. Knepley 
551c4762a1bSJed Brown   test:
552*8b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_conv
553c4762a1bSJed Brown     requires: !single
554*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [2.8, 2.4]
555*8b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -snes_convergence_estimate -convest_num_refine 1 \
556*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
557*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
558*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
559*8b0e23d0SMatthew G. Knepley 
560c4762a1bSJed Brown   test:
561*8b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_check
562*8b0e23d0SMatthew G. Knepley     requires: triangle
563*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
564*8b0e23d0SMatthew G. Knepley 
565*8b0e23d0SMatthew G. Knepley   test:
566*8b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_check
567*8b0e23d0SMatthew G. Knepley     requires: ctetgen !single
568*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -dmsnes_check 0.0001
569*8b0e23d0SMatthew G. Knepley 
570*8b0e23d0SMatthew G. Knepley   test:
571*8b0e23d0SMatthew G. Knepley     suffix: 2d_p3_p2_conv
572*8b0e23d0SMatthew G. Knepley     requires: triangle
573*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 gives L_2 convergence rate: [3.8, 3.0]
574*8b0e23d0SMatthew G. Knepley     args: -sol trig -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 -ksp_error_if_not_converged \
575*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
576*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
577*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
578*8b0e23d0SMatthew G. Knepley 
579*8b0e23d0SMatthew G. Knepley   test:
580*8b0e23d0SMatthew G. Knepley     suffix: 3d_p3_p2_conv
581*8b0e23d0SMatthew G. Knepley     requires: ctetgen long_runtime
582*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 1 -convest_num_refine 2 gives L_2 convergence rate: [3.6, 3.9]
583*8b0e23d0SMatthew G. Knepley     args: -sol trig -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 3 -pres_petscspace_degree 2 -snes_convergence_estimate -convest_num_refine 2 \
584*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -pc_use_amat \
585*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition a11 -pc_fieldsplit_off_diag_use_amat \
586*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type lu
587*8b0e23d0SMatthew G. Knepley 
588*8b0e23d0SMatthew G. Knepley   test:
589*8b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_conv
590c4762a1bSJed Brown     requires: !single
591*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 3 gives L_2 convergence rate: [1.9, 1.0]
592*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 2 \
593*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -petscds_jac_pre 0 \
594*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
595*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd
596*8b0e23d0SMatthew G. Knepley 
597c4762a1bSJed Brown   test:
598*8b0e23d0SMatthew G. Knepley     suffix: 3d_q1_p0_conv
599*8b0e23d0SMatthew G. Knepley     requires: !single
600*8b0e23d0SMatthew G. Knepley     # Using -dm_refine 2 -convest_num_refine 2 gives L_2 convergence rate: [1.7, 1.0]
601*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_refine 1 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -snes_convergence_estimate -convest_num_refine 1 \
602*8b0e23d0SMatthew G. Knepley       -ksp_atol 1e-10 -ksp_error_if_not_converged -petscds_jac_pre 0 \
603*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_schur_precondition full \
604*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_levels_pc_type jacobi -fieldsplit_pressure_mg_coarse_pc_type svd
605*8b0e23d0SMatthew G. Knepley 
606*8b0e23d0SMatthew G. Knepley   # Stokes preconditioners
607c4762a1bSJed Brown   #   Block diagonal \begin{pmatrix} A & 0 \\ 0 & I \end{pmatrix}
608c4762a1bSJed Brown   test:
609*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_diagonal
610*8b0e23d0SMatthew G. Knepley     requires: triangle
611*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
612*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
613*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-4 -ksp_error_if_not_converged \
614*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type additive -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
615c4762a1bSJed Brown   #   Block triangular \begin{pmatrix} A & B \\ 0 & I \end{pmatrix}
616c4762a1bSJed Brown   test:
617*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_block_triangular
618*8b0e23d0SMatthew G. Knepley     requires: triangle
619*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
620*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
621*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
622*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_pc_type jacobi
623c4762a1bSJed Brown   #   Diagonal Schur complement \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix}
624c4762a1bSJed Brown   test:
625*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_diagonal
626*8b0e23d0SMatthew G. Knepley     requires: triangle
627*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
628*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
629*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
630*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type diag -pc_fieldsplit_off_diag_use_amat \
631*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
632c4762a1bSJed Brown   #   Upper triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
633c4762a1bSJed Brown   test:
634*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_upper
635*8b0e23d0SMatthew G. Knepley     requires: triangle
636*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -dmsnes_check 0.0001 \
637*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
638*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type upper -pc_fieldsplit_off_diag_use_amat \
639*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
640c4762a1bSJed Brown   #   Lower triangular Schur complement \begin{pmatrix} A & B \\ 0 & S \end{pmatrix}
641c4762a1bSJed Brown   test:
642*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_lower
643*8b0e23d0SMatthew G. Knepley     requires: triangle
644*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
645*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
646*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
647*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type lower -pc_fieldsplit_off_diag_use_amat \
648*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
649c4762a1bSJed Brown   #   Full Schur complement \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & S \end{pmatrix} \begin{pmatrix} I & A^{-1} B \\ 0 & I \end{pmatrix}
650c4762a1bSJed Brown   test:
651*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_schur_full
652*8b0e23d0SMatthew G. Knepley     requires: triangle
653*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
654*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
655*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged -pc_use_amat \
656*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_off_diag_use_amat \
657*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi
658*8b0e23d0SMatthew G. Knepley   #   Full Schur + Velocity GMG
659*8b0e23d0SMatthew G. Knepley   test:
660*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_gmg_vcycle
661*8b0e23d0SMatthew G. Knepley     requires: triangle
662*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine_hierarchy 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 \
663*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
664*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-9 -ksp_error_if_not_converged -pc_use_amat \
665*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full -pc_fieldsplit_off_diag_use_amat \
666*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type mg -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type gamg -fieldsplit_pressure_mg_coarse_pc_type svd
667c4762a1bSJed Brown   #   SIMPLE \begin{pmatrix} I & 0 \\ B^T A^{-1} & I \end{pmatrix} \begin{pmatrix} A & 0 \\ 0 & B^T diag(A)^{-1} B \end{pmatrix} \begin{pmatrix} I & diag(A)^{-1} B \\ 0 & I \end{pmatrix}
668c4762a1bSJed Brown   test:
669*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_simple
670c4762a1bSJed Brown     requires: triangle
671*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
672*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
673*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
674*8b0e23d0SMatthew G. Knepley       -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full \
675*8b0e23d0SMatthew G. Knepley         -fieldsplit_velocity_pc_type lu -fieldsplit_pressure_ksp_rtol 1e-10 -fieldsplit_pressure_pc_type jacobi \
676*8b0e23d0SMatthew G. Knepley         -fieldsplit_pressure_inner_ksp_type preonly -fieldsplit_pressure_inner_pc_type jacobi -fieldsplit_pressure_upper_ksp_type preonly -fieldsplit_pressure_upper_pc_type jacobi
677c4762a1bSJed Brown   #   FETI-DP solvers (these solvers are quite inefficient, they are here to exercise the code)
678c4762a1bSJed Brown   test:
679*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_fetidp
680c4762a1bSJed Brown     requires: triangle mumps
681c4762a1bSJed Brown     nsize: 5
682*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_refine 2 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
683*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
684*8b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
685*8b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
686*8b0e23d0SMatthew G. Knepley         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \
687*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly
688c4762a1bSJed Brown   test:
689*8b0e23d0SMatthew G. Knepley     suffix: 2d_q2_q1_fetidp
690*8b0e23d0SMatthew G. Knepley     requires: mumps
691c4762a1bSJed Brown     nsize: 5
692*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
693*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
694*8b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
695*8b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
696*8b0e23d0SMatthew G. Knepley         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 200 -fetidp_fieldsplit_p_pc_type none \
697*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type mumps -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type mumps -fetidp_fieldsplit_lag_ksp_type preonly
698c4762a1bSJed Brown   test:
699*8b0e23d0SMatthew G. Knepley     suffix: 3d_p2_p1_fetidp
700*8b0e23d0SMatthew G. Knepley     requires: ctetgen mumps suitesparse
701c4762a1bSJed Brown     nsize: 5
702*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
703*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
704*8b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-9 -ksp_error_if_not_converged \
705*8b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
706*8b0e23d0SMatthew G. Knepley         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 1000 -fetidp_fieldsplit_p_pc_type none \
707*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_benign_trick -fetidp_bddc_pc_bddc_deluxe_singlemat \
708*8b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
709*8b0e23d0SMatthew G. Knepley         -fetidp_bddelta_pc_factor_mat_solver_type umfpack -fetidp_fieldsplit_lag_ksp_type preonly -fetidp_bddc_sub_schurs_mat_solver_type mumps -fetidp_bddc_sub_schurs_mat_mumps_icntl_14 100000 \
710*8b0e23d0SMatthew G. Knepley         -fetidp_bddelta_pc_factor_mat_ordering_type external \
711*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
712*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
713c4762a1bSJed Brown   test:
714*8b0e23d0SMatthew G. Knepley     suffix: 3d_q2_q1_fetidp
715c4762a1bSJed Brown     requires: suitesparse
716c4762a1bSJed Brown     nsize: 5
717*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_plex_box_dim 3 -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
718*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
719*8b0e23d0SMatthew G. Knepley       -ksp_type fetidp -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
720*8b0e23d0SMatthew G. Knepley       -ksp_fetidp_saddlepoint -fetidp_ksp_type cg \
721*8b0e23d0SMatthew G. Knepley         -fetidp_fieldsplit_p_ksp_max_it 1 -fetidp_fieldsplit_p_ksp_type richardson -fetidp_fieldsplit_p_ksp_richardson_scale 2000 -fetidp_fieldsplit_p_pc_type none \
722*8b0e23d0SMatthew G. Knepley         -fetidp_pc_discrete_harmonic -fetidp_harmonic_pc_factor_mat_solver_type petsc -fetidp_harmonic_pc_type cholesky \
723*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_symmetric -fetidp_fieldsplit_lag_ksp_type preonly \
724*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_solver_type umfpack -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_solver_type umfpack \
725*8b0e23d0SMatthew G. Knepley         -fetidp_bddc_pc_bddc_dirichlet_pc_factor_mat_ordering_type external -fetidp_bddc_pc_bddc_neumann_pc_factor_mat_ordering_type external
726*8b0e23d0SMatthew G. Knepley   #   BDDC solvers (these solvers are quite inefficient, they are here to exercise the code)
727c4762a1bSJed Brown   test:
728*8b0e23d0SMatthew G. Knepley     suffix: 2d_p2_p1_bddc
729*8b0e23d0SMatthew G. Knepley     nsize: 2
730c4762a1bSJed Brown     requires: triangle !single
731*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_faces 2,2,2 -dm_refine 1 -dm_mat_type is -dm_distribute -petscpartitioner_type simple -vel_petscspace_degree 2 -pres_petscspace_degree 1 -petscds_jac_pre 0 \
732*8b0e23d0SMatthew G. Knepley       -snes_error_if_not_converged \
733*8b0e23d0SMatthew G. Knepley       -ksp_type gmres -ksp_gmres_restart 100 -ksp_rtol 1.0e-8 -ksp_error_if_not_converged \
734*8b0e23d0SMatthew G. Knepley         -pc_type bddc -pc_bddc_corner_selection -pc_bddc_dirichlet_pc_type svd -pc_bddc_neumann_pc_type svd -pc_bddc_coarse_redundant_pc_type svd
735*8b0e23d0SMatthew G. Knepley   #   Vanka
736c4762a1bSJed Brown   test:
737*8b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka
738c4762a1bSJed Brown     requires: double !complex
739*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
740*8b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
741*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
742c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
743c4762a1bSJed Brown         -sub_ksp_type preonly -sub_pc_type lu
744c4762a1bSJed Brown   test:
745*8b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_vanka_denseinv
746c4762a1bSJed Brown     requires: double !complex
747*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
748*8b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
749*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
750c4762a1bSJed Brown       -pc_type patch -pc_patch_partition_of_unity 0 -pc_patch_construct_codim 0 -pc_patch_construct_type vanka \
751c4762a1bSJed Brown         -pc_patch_dense_inverse -pc_patch_sub_mat_type seqdense
752c4762a1bSJed Brown   #   Vanka smoother
753c4762a1bSJed Brown   test:
754*8b0e23d0SMatthew G. Knepley     suffix: 2d_q1_p0_gmg_vanka
755*8b0e23d0SMatthew G. Knepley     requires: double !complex
756*8b0e23d0SMatthew G. Knepley     args: -sol quadratic -dm_plex_box_simplex 0 -dm_refine_hierarchy 2 -vel_petscspace_degree 1 -pres_petscspace_degree 0 -petscds_jac_pre 0 \
757*8b0e23d0SMatthew G. Knepley       -snes_rtol 1.0e-4 \
758*8b0e23d0SMatthew G. Knepley       -ksp_type fgmres -ksp_atol 1e-5 -ksp_error_if_not_converged \
759*8b0e23d0SMatthew G. Knepley       -pc_type mg \
760*8b0e23d0SMatthew G. Knepley         -mg_levels_ksp_type gmres -mg_levels_ksp_max_it 30 \
761c4762a1bSJed Brown         -mg_levels_pc_type patch -mg_levels_pc_patch_partition_of_unity 0 -mg_levels_pc_patch_construct_codim 0 -mg_levels_pc_patch_construct_type vanka \
762c4762a1bSJed Brown           -mg_levels_sub_ksp_type preonly -mg_levels_sub_pc_type lu \
763c4762a1bSJed Brown         -mg_coarse_pc_type svd
764c4762a1bSJed Brown 
765c4762a1bSJed Brown TEST*/
766