xref: /petsc/src/snes/tutorials/ex24.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 static char help[] = "Poisson Problem in mixed form with 2d and 3d with finite elements.\n\
2 We solve the Poisson problem in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 This example supports automatic convergence estimation\n\
5 and Hdiv elements.\n\n\n";
6 
7 #include <petscdmplex.h>
8 #include <petscsnes.h>
9 #include <petscds.h>
10 #include <petscconvest.h>
11 
12 typedef enum {SOL_LINEAR, SOL_QUADRATIC, SOL_QUARTIC, SOL_UNKNOWN, NUM_SOLTYPE} SolType;
13 const char *SolTypeNames[NUM_SOLTYPE+3] = {"linear", "quadratic", "quartic", "unknown", "SolType", "SOL_", NULL};
14 
15 typedef struct {
16   /* Domain and mesh definition */
17   PetscInt  dim;     /* The topological mesh dimension */
18   PetscBool simplex; /* Simplicial mesh */
19   SolType   solType; /* The type of exact solution */
20 } AppCtx;
21 
22 /* 2D Dirichlet potential example
23 
24   u = x
25   q = <1, 0>
26   f = 0
27 
28   We will need a boundary integral of u over \Gamma.
29 */
30 static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
31 {
32   u[0] = x[0];
33   return 0;
34 }
35 static PetscErrorCode linear_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
36 {
37   PetscInt c;
38   for (c = 0; c < Nc; ++c) u[c] = c ? 0.0 : 1.0;
39   return 0;
40 }
41 
42 /* 2D Dirichlet potential example
43 
44   u = x^2 + y^2
45   q = <2x, 2y>
46   f = 4
47 
48   We will need a boundary integral of u over \Gamma.
49 */
50 static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
51 {
52   PetscInt d;
53 
54   u[0] = 0.0;
55   for (d = 0; d < dim; ++d) u[0] += x[d]*x[d];
56   return 0;
57 }
58 static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
59 {
60   PetscInt c;
61   for (c = 0; c < Nc; ++c) u[c] = 2.0*x[c];
62   return 0;
63 }
64 
65 /* 2D Dirichlet potential example
66 
67   u = x (1-x) y (1-y)
68   q = <(1-2x) y (1-y), x (1-x) (1-2y)>
69   f = -y (1-y) - x (1-x)
70 
71   u|_\Gamma = 0 so that the boundary integral vanishes
72 */
73 static PetscErrorCode quartic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
74 {
75   PetscInt d;
76 
77   u[0] = 1.0;
78   for (d = 0; d < dim; ++d) u[0] *= x[d]*(1.0 - x[d]);
79   return 0;
80 }
81 static PetscErrorCode quartic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
82 {
83   PetscInt c, d;
84 
85   for (c = 0; c < Nc; ++c) {
86     u[c] = 1.0;
87     for (d = 0; d < dim; ++d) {
88       if (c == d) u[c] *= 1 - 2.0*x[d];
89       else        u[c] *= x[d]*(1.0 - x[d]);
90     }
91   }
92   return 0;
93 }
94 
95 /* <v, -\nabla\cdot q> + <v, f> */
96 static void f0_linear_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
97                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
98                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
99                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
100 {
101   f0[0] = 0.0;
102 }
103 static void f0_bd_linear_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
105                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106                            PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
107 {
108   PetscScalar potential;
109   PetscInt    d;
110 
111   linear_u(dim, t, x, dim, &potential, NULL);
112   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
113 }
114 static void f0_quadratic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
115                            const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
116                            const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
117                            PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
118 {
119   PetscInt d;
120   f0[0] = 0.0;
121   for (d = 0; d < dim; ++d) {
122     f0[0] -= u_x[uOff_x[0]+d*dim+d];
123   }
124   f0[0] += 4.0;
125 }
126 static void f0_bd_quadratic_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
127                               const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
128                               const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
129                               PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
130 {
131   PetscScalar potential;
132   PetscInt    d;
133 
134   quadratic_u(dim, t, x, dim, &potential, NULL);
135   for (d = 0; d < dim; ++d) f0[d] = -potential*n[d];
136 }
137 static void f0_quartic_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
138                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
139                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
140                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
141 {
142   PetscInt d;
143   f0[0] = 0.0;
144   for (d = 0; d < dim; ++d) {
145     f0[0] -= u_x[uOff_x[0]+d*dim+d];
146     f0[0] += -2.0*x[d]*(1.0 - x[d]);
147   }
148 }
149 
150 /* <w, q> */
151 static void f0_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
152                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
153                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
154                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
155 {
156   PetscInt c;
157 
158   for (c = 0; c < dim; ++c) {
159     f0[c] = u[uOff[0]+c];
160   }
161 }
162 
163 /* <\nabla\cdot w, u> = <\nabla w, Iu> */
164 static void f1_q(PetscInt dim, PetscInt Nf, PetscInt NfAux,
165                  const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
166                  const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
167                  PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
168 {
169   PetscInt c, d;
170   for (c = 0; c < dim; ++c) {
171     for (d = 0; d < dim; ++d) {
172       if (c == d) f1[c*dim+d] = u[uOff[1]];
173     }
174   }
175 }
176 
177 /* <w, q> */
178 static void g0_qq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
179                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
180                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
181                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
182 {
183   PetscInt c;
184   for (c = 0; c < dim; ++c) g0[c*dim+c] = 1.0;
185 }
186 
187 /* <\nabla\cdot w, u> = <\nabla w, Iu> */
188 static void g2_qu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
189                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
190                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
191                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
192 {
193   PetscInt d;
194   for (d = 0; d < dim; ++d) g2[d*dim+d] = 1.0;
195 }
196 
197 /* <v, -\nabla\cdot q> */
198 static void g1_uq(PetscInt dim, PetscInt Nf, PetscInt NfAux,
199                   const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
200                   const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
201                   PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
202 {
203   PetscInt d;
204   for (d = 0; d < dim; ++d) g1[d*dim+d] = -1.0;
205 }
206 
207 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
208 {
209   PetscErrorCode ierr;
210 
211   PetscFunctionBeginUser;
212   options->dim     = 2;
213   options->simplex = PETSC_TRUE;
214   options->solType = SOL_LINEAR;
215 
216   ierr = PetscOptionsBegin(comm, "", "Poisson Problem Options", "DMPLEX");CHKERRQ(ierr);
217   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex24.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
218   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex24.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
219   ierr = PetscOptionsEnum("-sol_type", "Type of exact solution", "ex24.c", SolTypeNames, (PetscEnum) options->solType, (PetscEnum *) &options->solType, NULL);CHKERRQ(ierr);
220 
221   ierr = PetscOptionsEnd();
222   PetscFunctionReturn(0);
223 }
224 
225 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
226 {
227   PetscErrorCode ierr;
228 
229   PetscFunctionBeginUser;
230   if (0) {
231     DMLabel     label;
232     const char *name = "marker";
233 
234     ierr = DMPlexCreateReferenceCell(comm, user->dim, user->simplex, dm);CHKERRQ(ierr);
235     ierr = DMCreateLabel(*dm, name);CHKERRQ(ierr);
236     ierr = DMGetLabel(*dm, name, &label);CHKERRQ(ierr);
237     ierr = DMPlexMarkBoundaryFaces(*dm, 1, label);CHKERRQ(ierr);
238     ierr = DMPlexLabelComplete(*dm, label);CHKERRQ(ierr);
239   } else {
240     /* Create box mesh */
241     ierr = DMPlexCreateBoxMesh(comm, user->dim, user->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
242   }
243   /* Distribute mesh over processes */
244   {
245     DM               dmDist = NULL;
246     PetscPartitioner part;
247 
248     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
249     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
250     ierr = DMPlexDistribute(*dm, 0, NULL, &dmDist);CHKERRQ(ierr);
251     if (dmDist) {
252       ierr = DMDestroy(dm);CHKERRQ(ierr);
253       *dm  = dmDist;
254     }
255   }
256   /* TODO: This should be pulled into the library */
257   {
258     char      convType[256];
259     PetscBool flg;
260 
261     ierr = PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");CHKERRQ(ierr);
262     ierr = PetscOptionsFList("-dm_plex_convert_type","Convert DMPlex to another format","ex12",DMList,DMPLEX,convType,256,&flg);CHKERRQ(ierr);
263     ierr = PetscOptionsEnd();
264     if (flg) {
265       DM dmConv;
266 
267       ierr = DMConvert(*dm,convType,&dmConv);CHKERRQ(ierr);
268       if (dmConv) {
269         ierr = DMDestroy(dm);CHKERRQ(ierr);
270         *dm  = dmConv;
271       }
272     }
273   }
274   /* TODO: This should be pulled into the library */
275   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
276 
277   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
278   ierr = DMSetApplicationContext(*dm, user);CHKERRQ(ierr);
279   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
280   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
281   PetscFunctionReturn(0);
282 }
283 
284 static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user)
285 {
286   PetscDS        prob;
287   const PetscInt id = 1;
288   PetscErrorCode ierr;
289 
290   PetscFunctionBeginUser;
291   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
292   ierr = PetscDSSetResidual(prob, 0, f0_q, f1_q);CHKERRQ(ierr);
293   ierr = PetscDSSetJacobian(prob, 0, 0, g0_qq, NULL, NULL, NULL);CHKERRQ(ierr);
294   ierr = PetscDSSetJacobian(prob, 0, 1, NULL, NULL, g2_qu, NULL);CHKERRQ(ierr);
295   ierr = PetscDSSetJacobian(prob, 1, 0, NULL, g1_uq, NULL, NULL);CHKERRQ(ierr);
296   switch (user->solType)
297   {
298     case SOL_LINEAR:
299       ierr = PetscDSSetResidual(prob, 1, f0_linear_u, NULL);CHKERRQ(ierr);
300       ierr = PetscDSSetBdResidual(prob, 0, f0_bd_linear_q, NULL);CHKERRQ(ierr);
301       ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr);
302       ierr = PetscDSSetExactSolution(prob, 0, linear_q, user);CHKERRQ(ierr);
303       ierr = PetscDSSetExactSolution(prob, 1, linear_u, user);CHKERRQ(ierr);
304       break;
305     case SOL_QUADRATIC:
306       ierr = PetscDSSetResidual(prob, 1, f0_quadratic_u, NULL);CHKERRQ(ierr);
307       ierr = PetscDSSetBdResidual(prob, 0, f0_bd_quadratic_q, NULL);CHKERRQ(ierr);
308       ierr = PetscDSAddBoundary(prob, DM_BC_NATURAL, "Dirichlet Bd Integral", "marker", 0, 0, NULL, (void (*)(void)) NULL, 1, &id, user);CHKERRQ(ierr);
309       ierr = PetscDSSetExactSolution(prob, 0, quadratic_q, user);CHKERRQ(ierr);
310       ierr = PetscDSSetExactSolution(prob, 1, quadratic_u, user);CHKERRQ(ierr);
311       break;
312     case SOL_QUARTIC:
313       ierr = PetscDSSetResidual(prob, 1, f0_quartic_u, NULL);CHKERRQ(ierr);
314       ierr = PetscDSSetExactSolution(prob, 0, quartic_q, user);CHKERRQ(ierr);
315       ierr = PetscDSSetExactSolution(prob, 1, quartic_u, user);CHKERRQ(ierr);
316       break;
317     default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]);
318   }
319   PetscFunctionReturn(0);
320 }
321 
322 static PetscErrorCode SetupDiscretization(DM dm, PetscErrorCode (*setup)(DM, AppCtx *), AppCtx *user)
323 {
324   DM              cdm = dm;
325   PetscFE         feq, feu;
326   const PetscInt  dim = user->dim;
327   PetscErrorCode  ierr;
328 
329   PetscFunctionBeginUser;
330   /* Create finite element */
331   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, user->simplex, "field_",     -1, &feq);CHKERRQ(ierr);
332   ierr = PetscObjectSetName((PetscObject) feq, "field");CHKERRQ(ierr);
333   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1,   user->simplex, "potential_", -1, &feu);CHKERRQ(ierr);
334   ierr = PetscObjectSetName((PetscObject) feu, "potential");CHKERRQ(ierr);
335   ierr = PetscFECopyQuadrature(feq, feu);CHKERRQ(ierr);
336   /* Set discretization and boundary conditions for each mesh */
337   ierr = DMSetField(dm, 0, NULL, (PetscObject) feq);CHKERRQ(ierr);
338   ierr = DMSetField(dm, 1, NULL, (PetscObject) feu);CHKERRQ(ierr);
339   ierr = DMCreateDS(dm);CHKERRQ(ierr);
340   ierr = (*setup)(dm, user);CHKERRQ(ierr);
341   while (cdm) {
342     ierr = DMCopyDisc(dm,cdm);CHKERRQ(ierr);
343     /* TODO: Check whether the boundary of coarse meshes is marked */
344     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
345   }
346   ierr = PetscFEDestroy(&feq);CHKERRQ(ierr);
347   ierr = PetscFEDestroy(&feu);CHKERRQ(ierr);
348   PetscFunctionReturn(0);
349 }
350 
351 int main(int argc, char **argv)
352 {
353   DM             dm;   /* Problem specification */
354   SNES           snes; /* Nonlinear solver */
355   Vec            u;    /* Solutions */
356   AppCtx         user; /* User-defined work context */
357   PetscErrorCode ierr;
358 
359   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
360   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
361   /* Primal system */
362   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
363   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
364   ierr = SNESSetDM(snes, dm);CHKERRQ(ierr);
365   ierr = SetupDiscretization(dm, SetupPrimalProblem, &user);CHKERRQ(ierr);
366   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
367   ierr = VecSet(u, 0.0);CHKERRQ(ierr);
368   ierr = PetscObjectSetName((PetscObject) u, "potential");CHKERRQ(ierr);
369   ierr = DMPlexSetSNESLocalFEM(dm, &user, &user, &user);CHKERRQ(ierr);
370   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
371   ierr = DMSNESCheckFromOptions(snes, u, NULL, NULL);CHKERRQ(ierr);
372   ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr);
373   ierr = SNESGetSolution(snes, &u);CHKERRQ(ierr);
374   ierr = VecViewFromOptions(u, NULL, "-potential_view");CHKERRQ(ierr);
375   /* Cleanup */
376   ierr = VecDestroy(&u);CHKERRQ(ierr);
377   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
378   ierr = DMDestroy(&dm);CHKERRQ(ierr);
379   ierr = PetscFinalize();
380   return ierr;
381 }
382 
383 /*TEST
384 
385   test:
386     suffix: 2d_bdm1_p0_0
387     requires: triangle
388     args: -sol_type linear \
389           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
390           -dmsnes_check .001 -snes_error_if_not_converged \
391           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
392           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
393             -fieldsplit_field_pc_type lu \
394             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
395   test:
396     suffix: 2d_bdm1_p0_1
397     requires: triangle
398     args: -sol_type quadratic \
399           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
400           -dmsnes_check .001 -snes_error_if_not_converged \
401           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
402           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
403             -fieldsplit_field_pc_type lu \
404             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
405   test:
406     suffix: 2d_bdm1_p0_2
407     requires: triangle
408     args: -sol_type quartic \
409           -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
410           -snes_error_if_not_converged \
411           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
412           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
413             -fieldsplit_field_pc_type lu \
414             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
415 
416   test:
417     suffix: 2d_p2_p0_vtk
418     requires: triangle
419     args: -sol_type linear \
420           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
421           -dmsnes_check .001 -snes_error_if_not_converged \
422           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
423           -potential_view vtk: -exact_vec_view vtk: \
424           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
425             -fieldsplit_field_pc_type lu \
426             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
427 
428   test:
429     suffix: 2d_p2_p0_vtu
430     requires: triangle
431     args: -sol_type linear \
432           -field_petscspace_degree 2 -dm_refine 0 -convest_num_refine 1 -snes_convergence_estimate \
433           -dmsnes_check .001 -snes_error_if_not_converged \
434           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
435           -potential_view vtk:multifield.vtu -exact_vec_view vtk:exact.vtu \
436           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
437             -fieldsplit_field_pc_type lu \
438             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
439 TEST*/
440 
441 /* These tests will be active once tensor P^- is working
442 
443   test:
444     suffix: 2d_bdmq1_p0_0
445     requires: triangle
446     args: -simplex 0 -sol_type linear \
447           -field_petscspace_poly_type pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
448           -dmsnes_check .001 -snes_error_if_not_converged \
449           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
450           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
451             -fieldsplit_field_pc_type lu \
452             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
453   test:
454     suffix: 2d_bdmq1_p0_2
455     requires: triangle
456     args: -simplex 0 -sol_type quartic \
457           -field_petscspace_poly_type_no pminus_hdiv -field_petscspace_degree 1 -field_petscdualspace_type bdm -dm_refine 0 -convest_num_refine 3 -snes_convergence_estimate \
458           -dmsnes_check .001 -snes_error_if_not_converged \
459           -ksp_rtol 1e-10 -ksp_error_if_not_converged \
460           -pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_factorization_type full -pc_fieldsplit_schur_precondition full \
461             -fieldsplit_field_pc_type lu \
462             -fieldsplit_potential_ksp_rtol 1e-10 -fieldsplit_potential_pc_type lu
463 
464 */
465