xref: /petsc/src/ts/tutorials/ex45.c (revision 408cafa089b450f13fea5926aa80aaf73271e784)
1c4762a1bSJed Brown static char help[] = "Heat Equation in 2d and 3d with finite elements.\n\
2c4762a1bSJed Brown We solve the heat equation in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\
4c4762a1bSJed Brown Contributed by: Julian Andrej <juan@tf.uni-kiel.de>\n\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscdmplex.h>
7c4762a1bSJed Brown #include <petscds.h>
8c4762a1bSJed Brown #include <petscts.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown /*
11c4762a1bSJed Brown   Heat equation:
12c4762a1bSJed Brown 
13c4762a1bSJed Brown     du/dt - \Delta u = -1 * dim
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   Exact 2D solution:
16c4762a1bSJed Brown 
17c4762a1bSJed Brown     u = 2t + x^2 + y^2
18c4762a1bSJed Brown 
19c4762a1bSJed Brown     2 - (2 + 2) + 2 = 0
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   Exact 3D solution:
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     u = 3t + x^2 + y^2 + z^2
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     3 - (2 + 2 + 2) + 3 = 0
26c4762a1bSJed Brown */
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   PetscInt          dim;
30c4762a1bSJed Brown   PetscBool         simplex;
31c4762a1bSJed Brown   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown static PetscErrorCode analytic_temp(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   PetscInt d;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   *u = dim*time;
39c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
40c4762a1bSJed Brown   return 0;
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown static void f0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
44c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
45c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
46c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   f0[0] = u_t[0] + (PetscScalar) dim;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
52c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
53c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
54c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
55c4762a1bSJed Brown {
56c4762a1bSJed Brown   PetscInt d;
57c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
58c4762a1bSJed Brown     f1[d] = u_x[d];
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
63c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
64c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
65c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
66c4762a1bSJed Brown {
67c4762a1bSJed Brown   PetscInt d;
68c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
69c4762a1bSJed Brown     g3[d*dim+d] = 1.0;
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
74c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
75c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
76c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   g0[0] = u_tShift*1.0;
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   PetscErrorCode ierr;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   PetscFunctionBeginUser;
86c4762a1bSJed Brown   options->dim     = 2;
87c4762a1bSJed Brown   options->simplex = PETSC_TRUE;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex45.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex45.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[])
97c4762a1bSJed Brown {
98*408cafa0SMatthew G. Knepley   DM             plex;
99c4762a1bSJed Brown   DMLabel        label;
100c4762a1bSJed Brown   PetscErrorCode ierr;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
103c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
105*408cafa0SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
106*408cafa0SMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr);
107*408cafa0SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
112c4762a1bSJed Brown {
113c4762a1bSJed Brown   DM             pdm = NULL;
114c4762a1bSJed Brown   const PetscInt dim = ctx->dim;
115c4762a1bSJed Brown   PetscBool      hasLabel;
116c4762a1bSJed Brown   PetscErrorCode ierr;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBeginUser;
119c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
121c4762a1bSJed Brown   /* If no boundary marker exists, mark the whole boundary */
122c4762a1bSJed Brown   ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr);
123c4762a1bSJed Brown   if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);}
124c4762a1bSJed Brown   /* Distribute mesh over processes */
125c4762a1bSJed Brown   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
126c4762a1bSJed Brown   if (pdm) {
127c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
128c4762a1bSJed Brown     *dm  = pdm;
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
132c4762a1bSJed Brown   PetscFunctionReturn(0);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   PetscDS        prob;
138c4762a1bSJed Brown   const PetscInt id = 1;
139c4762a1bSJed Brown   PetscErrorCode ierr;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBeginUser;
142c4762a1bSJed Brown   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = PetscDSSetResidual(prob, 0, f0_temp, f1_temp);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = PetscDSSetJacobian(prob, 0, 0, g0_temp, NULL, NULL, g3_temp);CHKERRQ(ierr);
145c4762a1bSJed Brown   ctx->exactFuncs[0] = analytic_temp;
146*408cafa0SMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "marker", 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], 1, &id, ctx);CHKERRQ(ierr);
147c4762a1bSJed Brown   PetscFunctionReturn(0);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   DM             cdm = dm;
153c4762a1bSJed Brown   const PetscInt dim = ctx->dim;
154c4762a1bSJed Brown   PetscFE        fe;
155c4762a1bSJed Brown   PetscErrorCode ierr;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBeginUser;
158c4762a1bSJed Brown   /* Create finite element */
159c4762a1bSJed Brown   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1, ctx->simplex, "temp_", -1, &fe);CHKERRQ(ierr);
160c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, "temperature");CHKERRQ(ierr);
161c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
162c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
165c4762a1bSJed Brown   while (cdm) {
166c4762a1bSJed Brown     PetscBool hasLabel;
167c4762a1bSJed Brown 
168c4762a1bSJed Brown     ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
169c4762a1bSJed Brown     if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
170*408cafa0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
171c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
174c4762a1bSJed Brown   PetscFunctionReturn(0);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown int main(int argc, char **argv)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   AppCtx         ctx;
180c4762a1bSJed Brown   DM             dm;
181c4762a1bSJed Brown   TS             ts;
182c4762a1bSJed Brown   Vec            u, r;
183c4762a1bSJed Brown   PetscReal      t       = 0.0;
184c4762a1bSJed Brown   PetscReal      L2error = 0.0;
185c4762a1bSJed Brown   PetscErrorCode ierr;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
188c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = PetscMalloc1(1, &ctx.exactFuncs);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
195c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "temperature");CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
200c4762a1bSJed Brown   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
203c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
204c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   ierr = DMProjectFunction(dm, t, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
210c4762a1bSJed Brown   ierr = DMComputeL2Diff(dm, t, ctx.exactFuncs, NULL, u, &L2error);CHKERRQ(ierr);
211c4762a1bSJed Brown   if (L2error < 1.0e-11) {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: < 1.0e-11\n");CHKERRQ(ierr);}
212c4762a1bSJed Brown   else                   {ierr = PetscPrintf(PETSC_COMM_WORLD, "L_2 Error: %g\n", (double)L2error);CHKERRQ(ierr);}
213c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
216c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
219c4762a1bSJed Brown   ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr);
220c4762a1bSJed Brown   ierr = PetscFinalize();
221c4762a1bSJed Brown   return ierr;
222c4762a1bSJed Brown }
223c4762a1bSJed Brown 
224c4762a1bSJed Brown /*TEST
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   # Full solves
227c4762a1bSJed Brown   test:
228c4762a1bSJed Brown     suffix: 2d_p1_r1
229c4762a1bSJed Brown     requires: triangle
230c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
231c4762a1bSJed Brown     args: -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
232c4762a1bSJed Brown   test:
233c4762a1bSJed Brown     suffix: 2d_p1_r3
234c4762a1bSJed Brown     requires: triangle
235c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
236c4762a1bSJed Brown     args: -dm_refine 3 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
237c4762a1bSJed Brown   test:
238c4762a1bSJed Brown     suffix: 2d_p2_r1
239c4762a1bSJed Brown     requires: triangle
240c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
241c4762a1bSJed Brown     args: -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
242c4762a1bSJed Brown   test:
243c4762a1bSJed Brown     suffix: 2d_p2_r3
244c4762a1bSJed Brown     requires: triangle
245c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
246c4762a1bSJed Brown     args: -dm_refine 3 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
247c4762a1bSJed Brown   test:
248c4762a1bSJed Brown     suffix: 2d_q1_r1
249c4762a1bSJed Brown     requires: !single
250c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
251c4762a1bSJed Brown     args: -simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
252c4762a1bSJed Brown   test:
253c4762a1bSJed Brown     suffix: 2d_q1_r3
254c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
255c4762a1bSJed Brown     args: -simplex 0 -dm_refine 3 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
256c4762a1bSJed Brown   test:
257c4762a1bSJed Brown     suffix: 2d_q2_r1
258c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
259c4762a1bSJed Brown     args: -simplex 0 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
260c4762a1bSJed Brown   test:
261c4762a1bSJed Brown     suffix: 2d_q2_r3
262c4762a1bSJed Brown     requires: !single
263c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
264c4762a1bSJed Brown     args: -simplex 0 -dm_refine 3 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
265c4762a1bSJed Brown   test:
266c4762a1bSJed Brown     suffix: 3d_p1_r1
267c4762a1bSJed Brown     requires: ctetgen
268c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
269c4762a1bSJed Brown     args: -dim 3 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
270c4762a1bSJed Brown   test:
271c4762a1bSJed Brown     suffix: 3d_p1_r2
272c4762a1bSJed Brown     requires: ctetgen
273c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
274c4762a1bSJed Brown     args: -dim 3 -dm_refine 2 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
275c4762a1bSJed Brown   test:
276c4762a1bSJed Brown     suffix: 3d_p2_r1
277c4762a1bSJed Brown     requires: ctetgen
278c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
279c4762a1bSJed Brown     args: -dim 3 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
280c4762a1bSJed Brown   test:
281c4762a1bSJed Brown     suffix: 3d_p2_r2
282c4762a1bSJed Brown     requires: ctetgen
283c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
284c4762a1bSJed Brown     args: -dim 3 -dm_refine 2 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
285c4762a1bSJed Brown   test:
286c4762a1bSJed Brown     suffix: 3d_q1_r1
287c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
288c4762a1bSJed Brown     args: -dim 3 -simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
289c4762a1bSJed Brown   test:
290c4762a1bSJed Brown     suffix: 3d_q1_r2
291c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
292c4762a1bSJed Brown     args: -dim 3 -simplex 0 -dm_refine 2 -temp_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
293c4762a1bSJed Brown   test:
294c4762a1bSJed Brown     suffix: 3d_q2_r1
295c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
296c4762a1bSJed Brown     args: -dim 3 -simplex 0 -dm_refine 1 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
297c4762a1bSJed Brown   test:
298c4762a1bSJed Brown     suffix: 3d_q2_r2
299c4762a1bSJed Brown     filter: sed -e "s~ATOL~RTOL~g" -e "s~ABS~RELATIVE~g"
300c4762a1bSJed Brown     args: -dim 3 -simplex 0 -dm_refine 2 -temp_petscspace_degree 2 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -ksp_monitor_short -ksp_converged_reason -snes_monitor_short -snes_converged_reason -ts_monitor
301c4762a1bSJed Brown 
302c4762a1bSJed Brown TEST*/
303