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