xref: /petsc/src/ts/tutorials/ex45.c (revision 45480ffeba491aca5d3f3b8c78679069fb317d32)
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 
13a3d0cf85SMatthew G. Knepley     du/dt - \Delta u + f = 0
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16a3d0cf85SMatthew G. Knepley typedef enum {SOL_QUADRATIC_LINEAR, SOL_QUADRATIC_TRIG, SOL_TRIG_LINEAR, NUM_SOLUTION_TYPES} SolutionType;
17a3d0cf85SMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "unknown"};
18a3d0cf85SMatthew G. Knepley 
19c4762a1bSJed Brown typedef struct {
20a3d0cf85SMatthew G. Knepley   char         filename[PETSC_MAX_PATH_LEN];   /* Mesh filename */
21a3d0cf85SMatthew G. Knepley   char         bdfilename[PETSC_MAX_PATH_LEN]; /* Mesh boundary filename */
22a3d0cf85SMatthew G. Knepley   PetscReal    scale;                          /* Scale factor for mesh */
23a3d0cf85SMatthew G. Knepley   SolutionType solType;                        /* Type of exact solution */
24c4762a1bSJed Brown } AppCtx;
25c4762a1bSJed Brown 
26a3d0cf85SMatthew G. Knepley /*
27a3d0cf85SMatthew G. Knepley Exact 2D solution:
28a3d0cf85SMatthew G. Knepley   u = 2t + x^2 + y^2
29a3d0cf85SMatthew G. Knepley   F(u) = 2 - (2 + 2) + 2 = 0
30a3d0cf85SMatthew G. Knepley 
31a3d0cf85SMatthew G. Knepley Exact 3D solution:
32a3d0cf85SMatthew G. Knepley   u = 3t + x^2 + y^2 + z^2
33a3d0cf85SMatthew G. Knepley   F(u) = 3 - (2 + 2 + 2) + 3 = 0
34a3d0cf85SMatthew G. Knepley */
35a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
36c4762a1bSJed Brown {
37c4762a1bSJed Brown   PetscInt d;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   *u = dim*time;
40c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
41c4762a1bSJed Brown   return 0;
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
45a3d0cf85SMatthew G. Knepley {
46a3d0cf85SMatthew G. Knepley   *u = dim;
47a3d0cf85SMatthew G. Knepley   return 0;
48a3d0cf85SMatthew G. Knepley }
49a3d0cf85SMatthew G. Knepley 
50a3d0cf85SMatthew G. Knepley static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
51c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
52c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
53c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   f0[0] = u_t[0] + (PetscScalar) dim;
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58a3d0cf85SMatthew G. Knepley /*
59a3d0cf85SMatthew G. Knepley Exact 2D solution:
60a3d0cf85SMatthew G. Knepley   u = 2*cos(t) + x^2 + y^2
61a3d0cf85SMatthew G. Knepley   F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0
62a3d0cf85SMatthew G. Knepley 
63a3d0cf85SMatthew G. Knepley Exact 3D solution:
64a3d0cf85SMatthew G. Knepley   u = 3*cos(t) + x^2 + y^2 + z^2
65a3d0cf85SMatthew G. Knepley   F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0
66a3d0cf85SMatthew G. Knepley */
67a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
68a3d0cf85SMatthew G. Knepley {
69a3d0cf85SMatthew G. Knepley   PetscInt d;
70a3d0cf85SMatthew G. Knepley 
71a3d0cf85SMatthew G. Knepley   *u = dim*PetscCosReal(time);
72a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
73a3d0cf85SMatthew G. Knepley   return 0;
74a3d0cf85SMatthew G. Knepley }
75a3d0cf85SMatthew G. Knepley 
76a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
77a3d0cf85SMatthew G. Knepley {
78a3d0cf85SMatthew G. Knepley   *u = -dim*PetscSinReal(time);
79a3d0cf85SMatthew G. Knepley   return 0;
80a3d0cf85SMatthew G. Knepley }
81a3d0cf85SMatthew G. Knepley 
82a3d0cf85SMatthew G. Knepley static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
83a3d0cf85SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
84a3d0cf85SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
85a3d0cf85SMatthew G. Knepley                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
86a3d0cf85SMatthew G. Knepley {
87a3d0cf85SMatthew G. Knepley   f0[0] = u_t[0] + dim*(PetscSinReal(t) + 2.0);
88a3d0cf85SMatthew G. Knepley }
89a3d0cf85SMatthew G. Knepley 
90a3d0cf85SMatthew G. Knepley /*
91a3d0cf85SMatthew G. Knepley Exact 2D solution:
92a3d0cf85SMatthew G. Knepley   u = 2\pi^2 t + cos(\pi x) + cos(\pi y)
93a3d0cf85SMatthew G. Knepley   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
94a3d0cf85SMatthew G. Knepley 
95a3d0cf85SMatthew G. Knepley Exact 3D solution:
96a3d0cf85SMatthew G. Knepley   u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z)
97a3d0cf85SMatthew G. Knepley   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
98a3d0cf85SMatthew G. Knepley */
99a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
100a3d0cf85SMatthew G. Knepley {
101a3d0cf85SMatthew G. Knepley   PetscInt d;
102a3d0cf85SMatthew G. Knepley 
103a3d0cf85SMatthew G. Knepley   *u = dim*PetscSqr(PETSC_PI)*time;
104a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]);
105a3d0cf85SMatthew G. Knepley   return 0;
106a3d0cf85SMatthew G. Knepley }
107a3d0cf85SMatthew G. Knepley 
108a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
109a3d0cf85SMatthew G. Knepley {
110a3d0cf85SMatthew G. Knepley   *u = dim*PetscSqr(PETSC_PI);
111a3d0cf85SMatthew G. Knepley   return 0;
112a3d0cf85SMatthew G. Knepley }
113a3d0cf85SMatthew G. Knepley 
114a3d0cf85SMatthew G. Knepley static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
115a3d0cf85SMatthew G. Knepley                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
116a3d0cf85SMatthew G. Knepley                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
117a3d0cf85SMatthew G. Knepley                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
118a3d0cf85SMatthew G. Knepley {
119a3d0cf85SMatthew G. Knepley   PetscInt d;
120a3d0cf85SMatthew G. Knepley   f0[0] = u_t[0];
121a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*(PetscCosReal(PETSC_PI*x[d]) - 1.0);
122a3d0cf85SMatthew G. Knepley }
123a3d0cf85SMatthew G. Knepley 
124c4762a1bSJed Brown static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
125c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
126c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
127c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
128c4762a1bSJed Brown {
129c4762a1bSJed Brown   PetscInt d;
130a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
134c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
135c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
136c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   PetscInt d;
139a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
143c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
144c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
145c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
146c4762a1bSJed Brown {
147c4762a1bSJed Brown   g0[0] = u_tShift*1.0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
151c4762a1bSJed Brown {
152a3d0cf85SMatthew G. Knepley   PetscInt       sol;
153c4762a1bSJed Brown   PetscErrorCode ierr;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   PetscFunctionBeginUser;
156a3d0cf85SMatthew G. Knepley   options->filename[0]   = '\0';
157a3d0cf85SMatthew G. Knepley   options->bdfilename[0] = '\0';
158a3d0cf85SMatthew G. Knepley   options->scale         = 0.0;
159a3d0cf85SMatthew G. Knepley   options->solType       = SOL_QUADRATIC_LINEAR;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");CHKERRQ(ierr);
162a3d0cf85SMatthew G. Knepley   ierr = PetscOptionsString("-filename", "The mesh file", "ex45.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
163a3d0cf85SMatthew G. Knepley   ierr = PetscOptionsString("-bd_filename", "The mesh boundary file", "ex45.c", options->bdfilename, options->bdfilename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
164a3d0cf85SMatthew G. Knepley   ierr = PetscOptionsReal("-scale", "Scale factor for the mesh", "ex45.c", options->scale, &options->scale, NULL);CHKERRQ(ierr);
165a3d0cf85SMatthew G. Knepley   sol  = options->solType;
166a3d0cf85SMatthew G. Knepley   ierr = PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL);CHKERRQ(ierr);
167a3d0cf85SMatthew G. Knepley   options->solType = (SolutionType) sol;
168c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
169c4762a1bSJed Brown   PetscFunctionReturn(0);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[])
173c4762a1bSJed Brown {
174408cafa0SMatthew G. Knepley   DM             plex;
175c4762a1bSJed Brown   DMLabel        label;
176a3d0cf85SMatthew G. Knepley   PetscBool      hasLabel;
177c4762a1bSJed Brown   PetscErrorCode ierr;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   PetscFunctionBeginUser;
180a3d0cf85SMatthew G. Knepley   ierr = DMHasLabel(dm, name, &hasLabel);CHKERRQ(ierr);
181a3d0cf85SMatthew G. Knepley   if (hasLabel) PetscFunctionReturn(0);
182c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
183c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
184408cafa0SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
185408cafa0SMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(plex, 1, label);CHKERRQ(ierr);
186408cafa0SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
187c4762a1bSJed Brown   PetscFunctionReturn(0);
188c4762a1bSJed Brown }
189c4762a1bSJed Brown 
190c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
191c4762a1bSJed Brown {
192a3d0cf85SMatthew G. Knepley   size_t         len, lenbd;
193c4762a1bSJed Brown   PetscErrorCode ierr;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   PetscFunctionBeginUser;
196a3d0cf85SMatthew G. Knepley   ierr = PetscStrlen(ctx->filename,   &len);CHKERRQ(ierr);
197a3d0cf85SMatthew G. Knepley   ierr = PetscStrlen(ctx->bdfilename, &lenbd);CHKERRQ(ierr);
198a3d0cf85SMatthew G. Knepley   if (lenbd) {
199a3d0cf85SMatthew G. Knepley     DM bdm;
200a3d0cf85SMatthew G. Knepley 
201a3d0cf85SMatthew G. Knepley     ierr = DMPlexCreateFromFile(comm, ctx->bdfilename, PETSC_TRUE, &bdm);CHKERRQ(ierr);
202a3d0cf85SMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) bdm, "bd_");CHKERRQ(ierr);
203a3d0cf85SMatthew G. Knepley     ierr = DMSetFromOptions(bdm);CHKERRQ(ierr);
204a3d0cf85SMatthew G. Knepley     if (ctx->scale != 0.0) {
205a3d0cf85SMatthew G. Knepley       Vec coordinates, coordinatesLocal;
206a3d0cf85SMatthew G. Knepley 
207a3d0cf85SMatthew G. Knepley       ierr = DMGetCoordinates(bdm, &coordinates);CHKERRQ(ierr);
208a3d0cf85SMatthew G. Knepley       ierr = DMGetCoordinatesLocal(bdm, &coordinatesLocal);CHKERRQ(ierr);
209a3d0cf85SMatthew G. Knepley       ierr = VecScale(coordinates, ctx->scale);CHKERRQ(ierr);
210a3d0cf85SMatthew G. Knepley       ierr = VecScale(coordinatesLocal, ctx->scale);CHKERRQ(ierr);
211a3d0cf85SMatthew G. Knepley     }
212a3d0cf85SMatthew G. Knepley     ierr = DMViewFromOptions(bdm, NULL, "-dm_view");CHKERRQ(ierr);
213a3d0cf85SMatthew G. Knepley     ierr = DMPlexGenerate(bdm, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
214a3d0cf85SMatthew G. Knepley     ierr = DMDestroy(&bdm);CHKERRQ(ierr);
215a3d0cf85SMatthew G. Knepley   } else if (len) {
216a3d0cf85SMatthew G. Knepley     ierr = DMPlexCreateFromFile(comm, ctx->filename, PETSC_TRUE, dm);CHKERRQ(ierr);
217a3d0cf85SMatthew G. Knepley   } else {
218a3d0cf85SMatthew G. Knepley     ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
221a3d0cf85SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
222a3d0cf85SMatthew G. Knepley   ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
224c4762a1bSJed Brown   PetscFunctionReturn(0);
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
227c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
228c4762a1bSJed Brown {
229a3d0cf85SMatthew G. Knepley   PetscDS        ds;
230*45480ffeSMatthew G. Knepley   DMLabel        label;
231c4762a1bSJed Brown   const PetscInt id = 1;
232c4762a1bSJed Brown   PetscErrorCode ierr;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   PetscFunctionBeginUser;
235*45480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
236a3d0cf85SMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
237a3d0cf85SMatthew G. Knepley   ierr = PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp);CHKERRQ(ierr);
238a3d0cf85SMatthew G. Knepley   switch (ctx->solType) {
239a3d0cf85SMatthew G. Knepley     case SOL_QUADRATIC_LINEAR:
240a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_quad_lin,  f1_temp);CHKERRQ(ierr);
241a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx);CHKERRQ(ierr);
242a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx);CHKERRQ(ierr);
243*45480ffeSMatthew G. Knepley       ierr = 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);CHKERRQ(ierr);
244a3d0cf85SMatthew G. Knepley       break;
245a3d0cf85SMatthew G. Knepley     case SOL_QUADRATIC_TRIG:
246a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp);CHKERRQ(ierr);
247a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx);CHKERRQ(ierr);
248a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx);CHKERRQ(ierr);
249*45480ffeSMatthew G. Knepley       ierr = 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);CHKERRQ(ierr);
250a3d0cf85SMatthew G. Knepley       break;
251a3d0cf85SMatthew G. Knepley     case SOL_TRIG_LINEAR:
252a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetResidual(ds, 0, f0_trig_lin,  f1_temp);CHKERRQ(ierr);
253a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx);CHKERRQ(ierr);
254a3d0cf85SMatthew G. Knepley       ierr = PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx);CHKERRQ(ierr);
255*45480ffeSMatthew G. Knepley       ierr = 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);CHKERRQ(ierr);
256a3d0cf85SMatthew G. Knepley       break;
257a3d0cf85SMatthew G. Knepley     default: SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
258a3d0cf85SMatthew G. Knepley   }
259c4762a1bSJed Brown   PetscFunctionReturn(0);
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
263c4762a1bSJed Brown {
264c4762a1bSJed Brown   DM             cdm = dm;
265c4762a1bSJed Brown   PetscFE        fe;
266a3d0cf85SMatthew G. Knepley   DMPolytopeType ct;
267a3d0cf85SMatthew G. Knepley   PetscBool      simplex;
268a3d0cf85SMatthew G. Knepley   PetscInt       dim, cStart;
269c4762a1bSJed Brown   PetscErrorCode ierr;
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   PetscFunctionBeginUser;
272a3d0cf85SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
273a3d0cf85SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, NULL);CHKERRQ(ierr);
274a3d0cf85SMatthew G. Knepley   ierr = DMPlexGetCellType(dm, cStart, &ct);CHKERRQ(ierr);
275a3d0cf85SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
276c4762a1bSJed Brown   /* Create finite element */
277a3d0cf85SMatthew G. Knepley   ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe);CHKERRQ(ierr);
278c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, "temperature");CHKERRQ(ierr);
279c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
280c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
281c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
282c4762a1bSJed Brown   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
283c4762a1bSJed Brown   while (cdm) {
284a3d0cf85SMatthew G. Knepley     ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);
285408cafa0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
286c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
287c4762a1bSJed Brown   }
288c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
289c4762a1bSJed Brown   PetscFunctionReturn(0);
290c4762a1bSJed Brown }
291c4762a1bSJed Brown 
292a3d0cf85SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u)
293a3d0cf85SMatthew G. Knepley {
294a3d0cf85SMatthew G. Knepley   DM             dm;
295a3d0cf85SMatthew G. Knepley   PetscReal      t;
296a3d0cf85SMatthew G. Knepley   PetscErrorCode ierr;
297a3d0cf85SMatthew G. Knepley 
298a3d0cf85SMatthew G. Knepley   PetscFunctionBegin;
299a3d0cf85SMatthew G. Knepley   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
300a3d0cf85SMatthew G. Knepley   ierr = TSGetTime(ts, &t);CHKERRQ(ierr);
301a3d0cf85SMatthew G. Knepley   ierr = DMComputeExactSolution(dm, t, u, NULL);CHKERRQ(ierr);
302a3d0cf85SMatthew G. Knepley   PetscFunctionReturn(0);
303a3d0cf85SMatthew G. Knepley }
304a3d0cf85SMatthew G. Knepley 
305c4762a1bSJed Brown int main(int argc, char **argv)
306c4762a1bSJed Brown {
307c4762a1bSJed Brown   DM             dm;
308c4762a1bSJed Brown   TS             ts;
309a3d0cf85SMatthew G. Knepley   Vec            u;
310a3d0cf85SMatthew G. Knepley   AppCtx         ctx;
311c4762a1bSJed Brown   PetscErrorCode ierr;
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
314c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
315c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr);
316c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
317c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
321c4762a1bSJed Brown   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
322c4762a1bSJed Brown   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
323c4762a1bSJed Brown   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
324a3d0cf85SMatthew G. Knepley   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
325c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
326a3d0cf85SMatthew G. Knepley   ierr = TSSetComputeInitialCondition(ts, SetInitialConditions);CHKERRQ(ierr);
327c4762a1bSJed Brown 
328a3d0cf85SMatthew G. Knepley   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
329a3d0cf85SMatthew G. Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
330a3d0cf85SMatthew G. Knepley   ierr = SetInitialConditions(ts, u);CHKERRQ(ierr);
331a3d0cf85SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) u, "temperature");CHKERRQ(ierr);
332c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
333a3d0cf85SMatthew G. Knepley   ierr = DMTSCheckFromOptions(ts, u);CHKERRQ(ierr);
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
336c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
337c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
338c4762a1bSJed Brown   ierr = PetscFinalize();
339c4762a1bSJed Brown   return ierr;
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown /*TEST
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   test:
345a3d0cf85SMatthew G. Knepley     suffix: 2d_p1
346c4762a1bSJed Brown     requires: triangle
347a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
348a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
349c4762a1bSJed Brown   test:
350a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
351a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_sconv
352c4762a1bSJed Brown     requires: triangle
353a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
354a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
355c4762a1bSJed Brown   test:
356a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
357a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_tconv
358c4762a1bSJed Brown     requires: triangle
359a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
360a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
361c4762a1bSJed Brown   test:
362a3d0cf85SMatthew G. Knepley     suffix: 2d_p2
363c4762a1bSJed Brown     requires: triangle
364a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
365a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
366c4762a1bSJed Brown   test:
367a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
368a3d0cf85SMatthew G. Knepley     suffix: 2d_p2_sconv
369a3d0cf85SMatthew G. Knepley     requires: triangle
370a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
371a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
372c4762a1bSJed Brown   test:
373a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
374a3d0cf85SMatthew G. Knepley     suffix: 2d_p2_tconv
375a3d0cf85SMatthew G. Knepley     requires: triangle
376a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
377a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
378c4762a1bSJed Brown   test:
379a3d0cf85SMatthew G. Knepley     suffix: 2d_q1
380a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
381a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
382c4762a1bSJed Brown   test:
383a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
384a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_sconv
385a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
386a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
387c4762a1bSJed Brown   test:
388a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
389a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_tconv
390a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_box_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
391a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
392a3d0cf85SMatthew G. Knepley   test:
393a3d0cf85SMatthew G. Knepley     suffix: 2d_q2
394a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
395a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
396a3d0cf85SMatthew G. Knepley   test:
397a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
398a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_sconv
399a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -dm_plex_box_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
400a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
401a3d0cf85SMatthew G. Knepley   test:
402a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
403a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_tconv
404a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_box_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
405a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
406a3d0cf85SMatthew G. Knepley 
407a3d0cf85SMatthew G. Knepley   test:
408a3d0cf85SMatthew G. Knepley     suffix: 3d_p1
409c4762a1bSJed Brown     requires: ctetgen
410a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
411a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
412c4762a1bSJed Brown   test:
413a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
414a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_sconv
415c4762a1bSJed Brown     requires: ctetgen
416a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
417a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
418c4762a1bSJed Brown   test:
419a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
420a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_tconv
421c4762a1bSJed Brown     requires: ctetgen
422a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
423a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
424c4762a1bSJed Brown   test:
425a3d0cf85SMatthew G. Knepley     suffix: 3d_p2
426c4762a1bSJed Brown     requires: ctetgen
427a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
428a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
429c4762a1bSJed Brown   test:
430a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
431a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_sconv
432a3d0cf85SMatthew G. Knepley     requires: ctetgen
433a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
434a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
435c4762a1bSJed Brown   test:
436a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
437a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_tconv
438a3d0cf85SMatthew G. Knepley     requires: ctetgen
439a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
440a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
441c4762a1bSJed Brown   test:
442a3d0cf85SMatthew G. Knepley     suffix: 3d_q1
443a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
444a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
445c4762a1bSJed Brown   test:
446a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
447a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_sconv
448a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
449a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
450a3d0cf85SMatthew G. Knepley   test:
451a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
452a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_tconv
453a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_box_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
454a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
455a3d0cf85SMatthew G. Knepley   test:
456a3d0cf85SMatthew G. Knepley     suffix: 3d_q2
457a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_box_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
458a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
459a3d0cf85SMatthew G. Knepley   test:
460a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
461a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_sconv
462a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -dm_plex_box_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
463a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
464a3d0cf85SMatthew G. Knepley   test:
465a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
466a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_tconv
467a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_box_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
468a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
469a3d0cf85SMatthew G. Knepley 
470a3d0cf85SMatthew G. Knepley   test:
471a3d0cf85SMatthew G. Knepley     # 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
472a3d0cf85SMatthew G. Knepley     suffix: egads_sphere
473a3d0cf85SMatthew G. Knepley     requires: egads ctetgen
474a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -bd_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -scale 40 \
475a3d0cf85SMatthew G. Knepley           -temp_petscspace_degree 2 -dmts_check .0001 \
476a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
477c4762a1bSJed Brown 
478c4762a1bSJed Brown TEST*/
479