xref: /petsc/src/ts/tutorials/ex45.c (revision 6a5217c03994f2d95bb2e6dbd8bed42381aeb015)
1 static char help[] = "Heat Equation in 2d and 3d with finite elements.\n\
2 We solve the heat equation in a rectangular\n\
3 domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4 Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n";
5 
6 #include <petscdmplex.h>
7 #include <petscds.h>
8 #include <petscts.h>
9 
10 /*
11   Heat equation:
12 
13     du/dt - \Delta u + f = 0
14 */
15 
16 typedef enum {SOL_QUADRATIC_LINEAR, SOL_QUADRATIC_TRIG, SOL_TRIG_LINEAR, SOL_TRIG_TRIG, NUM_SOLUTION_TYPES} SolutionType;
17 const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "trig_trig", "unknown"};
18 
19 typedef struct {
20   SolutionType solType; /* Type of exact solution */
21   /* Solver setup */
22   PetscBool expTS;      /* Flag for explicit timestepping */
23   PetscBool lumped;     /* Lump the mass matrix */
24 } AppCtx;
25 
26 /*
27 Exact 2D solution:
28   u    = 2t + x^2 + y^2
29   u_t  = 2
30   \Delta u = 2 + 2 = 4
31   f    = 2
32   F(u) = 2 - (2 + 2) + 2 = 0
33 
34 Exact 3D solution:
35   u = 3t + x^2 + y^2 + z^2
36   F(u) = 3 - (2 + 2 + 2) + 3 = 0
37 */
38 static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
39 {
40   PetscInt d;
41 
42   *u = dim*time;
43   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
44   return 0;
45 }
46 
47 static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
48 {
49   *u = dim;
50   return 0;
51 }
52 
53 static void f0_quad_lin_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
54                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
55                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
56                             PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
57 {
58   f0[0] = -(PetscScalar) dim;
59 }
60 static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
61                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
62                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
63                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
64 {
65   PetscScalar exp[1] = {0.};
66   f0_quad_lin_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp);
67   f0[0] = u_t[0] - exp[0];
68 }
69 
70 /*
71 Exact 2D solution:
72   u = 2*cos(t) + x^2 + y^2
73   F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0
74 
75 Exact 3D solution:
76   u = 3*cos(t) + x^2 + y^2 + z^2
77   F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0
78 */
79 static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
80 {
81   PetscInt d;
82 
83   *u = dim*PetscCosReal(time);
84   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
85   return 0;
86 }
87 
88 static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
89 {
90   *u = -dim*PetscSinReal(time);
91   return 0;
92 }
93 
94 static void f0_quad_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98 {
99   f0[0] = -dim*(PetscSinReal(t) + 2.0);
100 }
101 static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
102                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
103                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
104                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
105 {
106   PetscScalar exp[1] = {0.};
107   f0_quad_trig_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp);
108   f0[0] = u_t[0] - exp[0];
109 }
110 
111 /*
112 Exact 2D solution:
113   u = 2\pi^2 t + cos(\pi x) + cos(\pi y)
114   F(u) = 2\pi^2 - \pi^2 (cos(\pi x) + cos(\pi y)) + \pi^2 (cos(\pi x) + cos(\pi y)) - 2\pi^2 = 0
115 
116 Exact 3D solution:
117   u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z)
118   F(u) = 3\pi^2 - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) + \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) - 3\pi^2 = 0
119 */
120 static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
121 {
122   PetscInt d;
123 
124   *u = dim*PetscSqr(PETSC_PI)*time;
125   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]);
126   return 0;
127 }
128 
129 static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
130 {
131   *u = dim*PetscSqr(PETSC_PI);
132   return 0;
133 }
134 
135 static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
136                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
137                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
138                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
139 {
140   PetscInt d;
141   f0[0] = u_t[0];
142   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*(PetscCosReal(PETSC_PI*x[d]) - 1.0);
143 }
144 
145 /*
146 Exact 2D solution:
147   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y)
148   u_t  = -pi^2 sin(t)
149   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y))
150   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y))
151   F(u) = -\pi^2 sin(t) + \pi^2 (cos(\pi x) + cos(\pi y)) - \pi^2 (cos(\pi x) + cos(\pi y)) + \pi^2 sin(t) = 0
152 
153 Exact 3D solution:
154   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z)
155   u_t  = -pi^2 sin(t)
156   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
157   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
158   F(u) = -\pi^2 sin(t) + \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z)) + \pi^2 sin(t) = 0
159 */
160 static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
161 {
162   PetscInt d;
163 
164   *u = PetscSqr(PETSC_PI)*PetscCosReal(time);
165   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]);
166   return 0;
167 }
168 
169 static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
170 {
171   *u = -PetscSqr(PETSC_PI)*PetscSinReal(time);
172   return 0;
173 }
174 
175 static void f0_trig_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
179 {
180   PetscInt d;
181   f0[0] -= PetscSqr(PETSC_PI)*PetscSinReal(t);
182   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*PetscCosReal(PETSC_PI*x[d]);
183 }
184 static void f0_trig_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
185                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
186                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
187                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
188 {
189   PetscScalar exp[1] = {0.};
190   f0_trig_trig_exp(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, a_t, a_x, t, x, numConstants, constants, exp);
191   f0[0] = u_t[0] - exp[0];
192 }
193 
194 static void f1_temp_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
195                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
196                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
197                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
198 {
199   PetscInt d;
200   for (d = 0; d < dim; ++d) f1[d] = -u_x[d];
201 }
202 static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
203                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
204                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
205                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
206 {
207   PetscInt d;
208   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
209 }
210 
211 static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
212                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
213                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
214                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
215 {
216   PetscInt d;
217   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
218 }
219 
220 static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
224 {
225   g0[0] = u_tShift*1.0;
226 }
227 
228 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
229 {
230   PetscInt       sol;
231   PetscErrorCode ierr;
232 
233   PetscFunctionBeginUser;
234   options->solType  = SOL_QUADRATIC_LINEAR;
235   options->expTS    = PETSC_FALSE;
236   options->lumped   = PETSC_TRUE;
237 
238   ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");PetscCall(ierr);
239   PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
240   options->solType = (SolutionType) sol;
241   PetscCall(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL));
242   PetscCall(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL));
243   ierr = PetscOptionsEnd();PetscCall(ierr);
244   PetscFunctionReturn(0);
245 }
246 
247 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
248 {
249   PetscFunctionBeginUser;
250   PetscCall(DMCreate(comm, dm));
251   PetscCall(DMSetType(*dm, DMPLEX));
252   PetscCall(DMSetFromOptions(*dm));
253   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
254   PetscFunctionReturn(0);
255 }
256 
257 static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
258 {
259   PetscDS        ds;
260   DMLabel        label;
261   const PetscInt id = 1;
262 
263   PetscFunctionBeginUser;
264   PetscCall(DMGetLabel(dm, "marker", &label));
265   PetscCall(DMGetDS(dm, &ds));
266   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp));
267   switch (ctx->solType) {
268     case SOL_QUADRATIC_LINEAR:
269       PetscCall(PetscDSSetResidual(ds, 0, f0_quad_lin,  f1_temp));
270       PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp));
271       PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx));
272       PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx));
273       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_quad_lin, (void (*)(void)) mms_quad_lin_t, ctx, NULL));
274       break;
275     case SOL_QUADRATIC_TRIG:
276       PetscCall(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp));
277       PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp));
278       PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx));
279       PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx));
280       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_quad_trig, (void (*)(void)) mms_quad_trig_t, ctx, NULL));
281       break;
282     case SOL_TRIG_LINEAR:
283       PetscCall(PetscDSSetResidual(ds, 0, f0_trig_lin,  f1_temp));
284       PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx));
285       PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx));
286       PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) mms_trig_lin, (void (*)(void)) mms_trig_lin_t, ctx, NULL));
287       break;
288     case SOL_TRIG_TRIG:
289       PetscCall(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp));
290       PetscCall(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp));
291       PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx));
292       PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx));
293       break;
294     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
295   }
296   PetscFunctionReturn(0);
297 }
298 
299 static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
300 {
301   DM             cdm = dm;
302   PetscFE        fe;
303   DMPolytopeType ct;
304   PetscBool      simplex;
305   PetscInt       dim, cStart;
306 
307   PetscFunctionBeginUser;
308   PetscCall(DMGetDimension(dm, &dim));
309   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
310   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
311   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
312   /* Create finite element */
313   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe));
314   PetscCall(PetscObjectSetName((PetscObject) fe, "temperature"));
315   /* Set discretization and boundary conditions for each mesh */
316   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
317   PetscCall(DMCreateDS(dm));
318   if (ctx->expTS) {
319     PetscDS ds;
320 
321     PetscCall(DMGetDS(dm, &ds));
322     PetscCall(PetscDSSetImplicit(ds, 0, PETSC_FALSE));
323   }
324   PetscCall(SetupProblem(dm, ctx));
325   while (cdm) {
326     PetscCall(DMCopyDisc(dm, cdm));
327     PetscCall(DMGetCoarseDM(cdm, &cdm));
328   }
329   PetscCall(PetscFEDestroy(&fe));
330   PetscFunctionReturn(0);
331 }
332 
333 static PetscErrorCode SetInitialConditions(TS ts, Vec u)
334 {
335   DM             dm;
336   PetscReal      t;
337 
338   PetscFunctionBegin;
339   PetscCall(TSGetDM(ts, &dm));
340   PetscCall(TSGetTime(ts, &t));
341   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
342   PetscFunctionReturn(0);
343 }
344 
345 int main(int argc, char **argv)
346 {
347   DM             dm;
348   TS             ts;
349   Vec            u;
350   AppCtx         ctx;
351 
352   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
353   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
354   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
355   PetscCall(DMSetApplicationContext(dm, &ctx));
356   PetscCall(SetupDiscretization(dm, &ctx));
357 
358   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
359   PetscCall(TSSetDM(ts, dm));
360   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
361   if (ctx.expTS) {
362     PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx));
363     if (ctx.lumped) PetscCall(DMTSCreateRHSMassMatrixLumped(dm));
364     else            PetscCall(DMTSCreateRHSMassMatrix(dm));
365   } else {
366     PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
367     PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
368   }
369   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
370   PetscCall(TSSetFromOptions(ts));
371   PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions));
372 
373   PetscCall(DMCreateGlobalVector(dm, &u));
374   PetscCall(DMTSCheckFromOptions(ts, u));
375   PetscCall(SetInitialConditions(ts, u));
376   PetscCall(PetscObjectSetName((PetscObject) u, "temperature"));
377   PetscCall(TSSolve(ts, u));
378   PetscCall(DMTSCheckFromOptions(ts, u));
379   if (ctx.expTS) PetscCall(DMTSDestroyRHSMassMatrix(dm));
380 
381   PetscCall(VecDestroy(&u));
382   PetscCall(TSDestroy(&ts));
383   PetscCall(DMDestroy(&dm));
384   PetscCall(PetscFinalize());
385   return 0;
386 }
387 
388 /*TEST
389 
390   test:
391     suffix: 2d_p1
392     requires: triangle
393     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
394           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
395   test:
396     suffix: 2d_p1_exp
397     requires: triangle
398     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \
399           -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error
400   test:
401     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
402     suffix: 2d_p1_sconv
403     requires: triangle
404     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
405           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
406   test:
407     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1]
408     suffix: 2d_p1_sconv_2
409     requires: triangle
410     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
411           -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu
412   test:
413     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
414     suffix: 2d_p1_tconv
415     requires: triangle
416     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
417           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
418   test:
419     # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0]
420     suffix: 2d_p1_tconv_2
421     requires: triangle
422     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
423           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
424   test:
425     # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid
426     suffix: 2d_p1_exp_tconv_2
427     requires: triangle
428     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
429           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu
430   test:
431     # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid
432     suffix: 2d_p1_exp_tconv_2_lump
433     requires: triangle
434     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
435           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4
436   test:
437     suffix: 2d_p2
438     requires: triangle
439     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
440           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
441   test:
442     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
443     suffix: 2d_p2_sconv
444     requires: triangle
445     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
446           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
447   test:
448     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1]
449     suffix: 2d_p2_sconv_2
450     requires: triangle
451     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
452           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
453   test:
454     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
455     suffix: 2d_p2_tconv
456     requires: triangle
457     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
458           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
459   test:
460     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
461     suffix: 2d_p2_tconv_2
462     requires: triangle
463     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
464           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
465   test:
466     suffix: 2d_q1
467     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
468           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
469   test:
470     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
471     suffix: 2d_q1_sconv
472     args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
473           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
474   test:
475     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
476     suffix: 2d_q1_tconv
477     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
478           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
479   test:
480     suffix: 2d_q2
481     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
482           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
483   test:
484     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
485     suffix: 2d_q2_sconv
486     args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
487           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
488   test:
489     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
490     suffix: 2d_q2_tconv
491     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
492           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
493 
494   test:
495     suffix: 3d_p1
496     requires: ctetgen
497     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
498           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
499   test:
500     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
501     suffix: 3d_p1_sconv
502     requires: ctetgen
503     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
504           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
505   test:
506     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
507     suffix: 3d_p1_tconv
508     requires: ctetgen
509     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
510           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
511   test:
512     suffix: 3d_p2
513     requires: ctetgen
514     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
515           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
516   test:
517     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
518     suffix: 3d_p2_sconv
519     requires: ctetgen
520     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
521           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
522   test:
523     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
524     suffix: 3d_p2_tconv
525     requires: ctetgen
526     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
527           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
528   test:
529     suffix: 3d_q1
530     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
531           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
532   test:
533     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
534     suffix: 3d_q1_sconv
535     args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
536           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
537   test:
538     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
539     suffix: 3d_q1_tconv
540     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
541           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
542   test:
543     suffix: 3d_q2
544     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
545           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
546   test:
547     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
548     suffix: 3d_q2_sconv
549     args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
550           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
551   test:
552     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
553     suffix: 3d_q2_tconv
554     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
555           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
556 
557   test:
558     # For a nice picture, -bd_dm_refine 2 -dm_refine 1 -dm_view hdf5:${PETSC_DIR}/sol.h5 -ts_monitor_solution hdf5:${PETSC_DIR}/sol.h5::append
559     suffix: egads_sphere
560     requires: egads ctetgen
561     args: -sol_type quadratic_linear \
562           -dm_plex_boundary_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_plex_boundary_label marker -bd_dm_plex_scale 40 \
563           -temp_petscspace_degree 2 -dmts_check .0001 \
564           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
565 
566 TEST*/
567