xref: /petsc/src/ts/tutorials/ex45.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
16742ee2edSMatthew G. Knepley typedef enum {SOL_QUADRATIC_LINEAR, SOL_QUADRATIC_TRIG, SOL_TRIG_LINEAR, SOL_TRIG_TRIG, NUM_SOLUTION_TYPES} SolutionType;
17742ee2edSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES+1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "trig_trig", "unknown"};
18a3d0cf85SMatthew G. Knepley 
19c4762a1bSJed Brown typedef struct {
20a3d0cf85SMatthew G. Knepley   SolutionType solType; /* Type of exact solution */
21742ee2edSMatthew G. Knepley   /* Solver setup */
22742ee2edSMatthew G. Knepley   PetscBool expTS;      /* Flag for explicit timestepping */
23742ee2edSMatthew G. Knepley   PetscBool lumped;     /* Lump the mass matrix */
24c4762a1bSJed Brown } AppCtx;
25c4762a1bSJed Brown 
26a3d0cf85SMatthew G. Knepley /*
27a3d0cf85SMatthew G. Knepley Exact 2D solution:
28a3d0cf85SMatthew G. Knepley   u    = 2t + x^2 + y^2
29742ee2edSMatthew G. Knepley   u_t  = 2
30742ee2edSMatthew G. Knepley   \Delta u = 2 + 2 = 4
31742ee2edSMatthew G. Knepley   f    = 2
32a3d0cf85SMatthew G. Knepley   F(u) = 2 - (2 + 2) + 2 = 0
33a3d0cf85SMatthew G. Knepley 
34a3d0cf85SMatthew G. Knepley Exact 3D solution:
35a3d0cf85SMatthew G. Knepley   u = 3t + x^2 + y^2 + z^2
36a3d0cf85SMatthew G. Knepley   F(u) = 3 - (2 + 2 + 2) + 3 = 0
37a3d0cf85SMatthew G. Knepley */
38a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
39c4762a1bSJed Brown {
40c4762a1bSJed Brown   PetscInt d;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   *u = dim*time;
43c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
44c4762a1bSJed Brown   return 0;
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
48a3d0cf85SMatthew G. Knepley {
49a3d0cf85SMatthew G. Knepley   *u = dim;
50a3d0cf85SMatthew G. Knepley   return 0;
51a3d0cf85SMatthew G. Knepley }
52a3d0cf85SMatthew G. Knepley 
53742ee2edSMatthew G. Knepley static void f0_quad_lin_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
54742ee2edSMatthew G. Knepley                             const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
55742ee2edSMatthew G. Knepley                             const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
56742ee2edSMatthew G. Knepley                             PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
57742ee2edSMatthew G. Knepley {
58742ee2edSMatthew G. Knepley   f0[0] = -(PetscScalar) dim;
59742ee2edSMatthew G. Knepley }
60a3d0cf85SMatthew G. Knepley static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
61c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
62c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
63c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
64c4762a1bSJed Brown {
65742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
66742ee2edSMatthew G. Knepley   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);
67742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70a3d0cf85SMatthew G. Knepley /*
71a3d0cf85SMatthew G. Knepley Exact 2D solution:
72a3d0cf85SMatthew G. Knepley   u = 2*cos(t) + x^2 + y^2
73a3d0cf85SMatthew G. Knepley   F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0
74a3d0cf85SMatthew G. Knepley 
75a3d0cf85SMatthew G. Knepley Exact 3D solution:
76a3d0cf85SMatthew G. Knepley   u = 3*cos(t) + x^2 + y^2 + z^2
77a3d0cf85SMatthew G. Knepley   F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0
78a3d0cf85SMatthew G. Knepley */
79a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
80a3d0cf85SMatthew G. Knepley {
81a3d0cf85SMatthew G. Knepley   PetscInt d;
82a3d0cf85SMatthew G. Knepley 
83a3d0cf85SMatthew G. Knepley   *u = dim*PetscCosReal(time);
84a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += x[d]*x[d];
85a3d0cf85SMatthew G. Knepley   return 0;
86a3d0cf85SMatthew G. Knepley }
87a3d0cf85SMatthew G. Knepley 
88a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
89a3d0cf85SMatthew G. Knepley {
90a3d0cf85SMatthew G. Knepley   *u = -dim*PetscSinReal(time);
91a3d0cf85SMatthew G. Knepley   return 0;
92a3d0cf85SMatthew G. Knepley }
93a3d0cf85SMatthew G. Knepley 
94742ee2edSMatthew G. Knepley static void f0_quad_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
95742ee2edSMatthew G. Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
96742ee2edSMatthew G. Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
97742ee2edSMatthew G. Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
98742ee2edSMatthew G. Knepley {
99742ee2edSMatthew G. Knepley   f0[0] = -dim*(PetscSinReal(t) + 2.0);
100742ee2edSMatthew G. Knepley }
101a3d0cf85SMatthew G. Knepley static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
102a3d0cf85SMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
103a3d0cf85SMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
104a3d0cf85SMatthew G. Knepley                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
105a3d0cf85SMatthew G. Knepley {
106742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
107742ee2edSMatthew G. Knepley   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);
108742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
109a3d0cf85SMatthew G. Knepley }
110a3d0cf85SMatthew G. Knepley 
111a3d0cf85SMatthew G. Knepley /*
112a3d0cf85SMatthew G. Knepley Exact 2D solution:
113a3d0cf85SMatthew G. Knepley   u = 2\pi^2 t + cos(\pi x) + cos(\pi y)
114a3d0cf85SMatthew 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
115a3d0cf85SMatthew G. Knepley 
116a3d0cf85SMatthew G. Knepley Exact 3D solution:
117a3d0cf85SMatthew G. Knepley   u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z)
118a3d0cf85SMatthew 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
119a3d0cf85SMatthew G. Knepley */
120a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
121a3d0cf85SMatthew G. Knepley {
122a3d0cf85SMatthew G. Knepley   PetscInt d;
123a3d0cf85SMatthew G. Knepley 
124a3d0cf85SMatthew G. Knepley   *u = dim*PetscSqr(PETSC_PI)*time;
125a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]);
126a3d0cf85SMatthew G. Knepley   return 0;
127a3d0cf85SMatthew G. Knepley }
128a3d0cf85SMatthew G. Knepley 
129a3d0cf85SMatthew G. Knepley static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
130a3d0cf85SMatthew G. Knepley {
131a3d0cf85SMatthew G. Knepley   *u = dim*PetscSqr(PETSC_PI);
132a3d0cf85SMatthew G. Knepley   return 0;
133a3d0cf85SMatthew G. Knepley }
134a3d0cf85SMatthew G. Knepley 
135a3d0cf85SMatthew G. Knepley static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux,
136a3d0cf85SMatthew G. Knepley                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
137a3d0cf85SMatthew G. Knepley                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
138a3d0cf85SMatthew G. Knepley                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
139a3d0cf85SMatthew G. Knepley {
140a3d0cf85SMatthew G. Knepley   PetscInt d;
141a3d0cf85SMatthew G. Knepley   f0[0] = u_t[0];
142a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*(PetscCosReal(PETSC_PI*x[d]) - 1.0);
143a3d0cf85SMatthew G. Knepley }
144a3d0cf85SMatthew G. Knepley 
145742ee2edSMatthew G. Knepley /*
146742ee2edSMatthew G. Knepley Exact 2D solution:
147742ee2edSMatthew G. Knepley   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y)
148742ee2edSMatthew G. Knepley   u_t  = -pi^2 sin(t)
149742ee2edSMatthew G. Knepley   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y))
150742ee2edSMatthew G. Knepley   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y))
151742ee2edSMatthew G. Knepley   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
152742ee2edSMatthew G. Knepley 
153742ee2edSMatthew G. Knepley Exact 3D solution:
154742ee2edSMatthew G. Knepley   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z)
155742ee2edSMatthew G. Knepley   u_t  = -pi^2 sin(t)
156742ee2edSMatthew G. Knepley   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
157742ee2edSMatthew G. Knepley   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
158742ee2edSMatthew G. Knepley   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
159742ee2edSMatthew G. Knepley */
160742ee2edSMatthew G. Knepley static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
161742ee2edSMatthew G. Knepley {
162742ee2edSMatthew G. Knepley   PetscInt d;
163742ee2edSMatthew G. Knepley 
164742ee2edSMatthew G. Knepley   *u = PetscSqr(PETSC_PI)*PetscCosReal(time);
165742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI*x[d]);
166742ee2edSMatthew G. Knepley   return 0;
167742ee2edSMatthew G. Knepley }
168742ee2edSMatthew G. Knepley 
169742ee2edSMatthew G. Knepley static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
170742ee2edSMatthew G. Knepley {
171742ee2edSMatthew G. Knepley   *u = -PetscSqr(PETSC_PI)*PetscSinReal(time);
172742ee2edSMatthew G. Knepley   return 0;
173742ee2edSMatthew G. Knepley }
174742ee2edSMatthew G. Knepley 
175742ee2edSMatthew G. Knepley static void f0_trig_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176742ee2edSMatthew G. Knepley                              const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177742ee2edSMatthew G. Knepley                              const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178742ee2edSMatthew G. Knepley                              PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
179742ee2edSMatthew G. Knepley {
180742ee2edSMatthew G. Knepley   PetscInt d;
181742ee2edSMatthew G. Knepley   f0[0] -= PetscSqr(PETSC_PI)*PetscSinReal(t);
182742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI)*PetscCosReal(PETSC_PI*x[d]);
183742ee2edSMatthew G. Knepley }
184742ee2edSMatthew G. Knepley static void f0_trig_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux,
185742ee2edSMatthew G. Knepley                          const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
186742ee2edSMatthew G. Knepley                          const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
187742ee2edSMatthew G. Knepley                          PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
188742ee2edSMatthew G. Knepley {
189742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
190742ee2edSMatthew G. Knepley   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);
191742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
192742ee2edSMatthew G. Knepley }
193742ee2edSMatthew G. Knepley 
194742ee2edSMatthew G. Knepley static void f1_temp_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
195742ee2edSMatthew G. Knepley                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
196742ee2edSMatthew G. Knepley                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
197742ee2edSMatthew G. Knepley                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
198742ee2edSMatthew G. Knepley {
199742ee2edSMatthew G. Knepley   PetscInt d;
200742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = -u_x[d];
201742ee2edSMatthew G. Knepley }
202c4762a1bSJed Brown static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
203c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
204c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
205c4762a1bSJed Brown                     PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
206c4762a1bSJed Brown {
207c4762a1bSJed Brown   PetscInt d;
208a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
211c4762a1bSJed Brown static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
212c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
213c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
214c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
215c4762a1bSJed Brown {
216c4762a1bSJed Brown   PetscInt d;
217a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux,
221c4762a1bSJed Brown                     const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
222c4762a1bSJed Brown                     const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
223c4762a1bSJed Brown                     PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
224c4762a1bSJed Brown {
225c4762a1bSJed Brown   g0[0] = u_tShift*1.0;
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
229c4762a1bSJed Brown {
230a3d0cf85SMatthew G. Knepley   PetscInt       sol;
231c4762a1bSJed Brown   PetscErrorCode ierr;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBeginUser;
234a3d0cf85SMatthew G. Knepley   options->solType  = SOL_QUADRATIC_LINEAR;
235742ee2edSMatthew G. Knepley   options->expTS    = PETSC_FALSE;
236742ee2edSMatthew G. Knepley   options->lumped   = PETSC_TRUE;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");CHKERRQ(ierr);
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
240a3d0cf85SMatthew G. Knepley   options->solType = (SolutionType) sol;
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL));
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL));
243c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
244c4762a1bSJed Brown   PetscFunctionReturn(0);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
248c4762a1bSJed Brown {
249c4762a1bSJed Brown   PetscFunctionBeginUser;
250*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
252*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
254c4762a1bSJed Brown   PetscFunctionReturn(0);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
258c4762a1bSJed Brown {
259a3d0cf85SMatthew G. Knepley   PetscDS        ds;
26045480ffeSMatthew G. Knepley   DMLabel        label;
261c4762a1bSJed Brown   const PetscInt id = 1;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBeginUser;
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLabel(dm, "marker", &label));
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &ds));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp));
267a3d0cf85SMatthew G. Knepley   switch (ctx->solType) {
268a3d0cf85SMatthew G. Knepley     case SOL_QUADRATIC_LINEAR:
269*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_quad_lin,  f1_temp));
270*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp));
271*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx));
272*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx));
273*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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));
274a3d0cf85SMatthew G. Knepley       break;
275a3d0cf85SMatthew G. Knepley     case SOL_QUADRATIC_TRIG:
276*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp));
277*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp));
278*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx));
279*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx));
280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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));
281a3d0cf85SMatthew G. Knepley       break;
282a3d0cf85SMatthew G. Knepley     case SOL_TRIG_LINEAR:
283*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_lin,  f1_temp));
284*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx));
285*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx));
286*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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));
287a3d0cf85SMatthew G. Knepley       break;
288742ee2edSMatthew G. Knepley     case SOL_TRIG_TRIG:
289*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp));
290*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp));
291*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx));
292*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx));
293742ee2edSMatthew G. Knepley       break;
29498921bdaSJacob Faibussowitsch     default: SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%D)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
295a3d0cf85SMatthew G. Knepley   }
296c4762a1bSJed Brown   PetscFunctionReturn(0);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   DM             cdm = dm;
302c4762a1bSJed Brown   PetscFE        fe;
303a3d0cf85SMatthew G. Knepley   DMPolytopeType ct;
304a3d0cf85SMatthew G. Knepley   PetscBool      simplex;
305a3d0cf85SMatthew G. Knepley   PetscInt       dim, cStart;
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   PetscFunctionBeginUser;
308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetCellType(dm, cStart, &ct));
311a3d0cf85SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct)+1 ? PETSC_TRUE : PETSC_FALSE;
312c4762a1bSJed Brown   /* Create finite element */
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe));
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) fe, "temperature"));
315c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
316*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetField(dm, 0, NULL, (PetscObject) fe));
317*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateDS(dm));
318742ee2edSMatthew G. Knepley   if (ctx->expTS) {
319742ee2edSMatthew G. Knepley     PetscDS ds;
320742ee2edSMatthew G. Knepley 
321*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetDS(dm, &ds));
322*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDSSetImplicit(ds, 0, PETSC_FALSE));
323742ee2edSMatthew G. Knepley   }
324*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupProblem(dm, ctx));
325c4762a1bSJed Brown   while (cdm) {
326*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, cdm));
327*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoarseDM(cdm, &cdm));
328c4762a1bSJed Brown   }
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFEDestroy(&fe));
330c4762a1bSJed Brown   PetscFunctionReturn(0);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
333a3d0cf85SMatthew G. Knepley static PetscErrorCode SetInitialConditions(TS ts, Vec u)
334a3d0cf85SMatthew G. Knepley {
335a3d0cf85SMatthew G. Knepley   DM             dm;
336a3d0cf85SMatthew G. Knepley   PetscReal      t;
337a3d0cf85SMatthew G. Knepley 
338a3d0cf85SMatthew G. Knepley   PetscFunctionBegin;
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts, &dm));
340*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts, &t));
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMComputeExactSolution(dm, t, u, NULL));
342a3d0cf85SMatthew G. Knepley   PetscFunctionReturn(0);
343a3d0cf85SMatthew G. Knepley }
344a3d0cf85SMatthew G. Knepley 
345c4762a1bSJed Brown int main(int argc, char **argv)
346c4762a1bSJed Brown {
347c4762a1bSJed Brown   DM             dm;
348c4762a1bSJed Brown   TS             ts;
349a3d0cf85SMatthew G. Knepley   Vec            u;
350a3d0cf85SMatthew G. Knepley   AppCtx         ctx;
351c4762a1bSJed Brown   PetscErrorCode ierr;
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
354*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx));
355*5f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
356*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(dm, &ctx));
357*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetupDiscretization(dm, &ctx));
358c4762a1bSJed Brown 
359*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD, &ts));
360*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts, dm));
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
362742ee2edSMatthew G. Knepley   if (ctx.expTS) {
363*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx));
364*5f80ce2aSJacob Faibussowitsch     if (ctx.lumped) CHKERRQ(DMTSCreateRHSMassMatrixLumped(dm));
365*5f80ce2aSJacob Faibussowitsch     else            CHKERRQ(DMTSCreateRHSMassMatrix(dm));
366742ee2edSMatthew G. Knepley   } else {
367*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
368*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
369742ee2edSMatthew G. Knepley   }
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
371*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetComputeInitialCondition(ts, SetInitialConditions));
373c4762a1bSJed Brown 
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(dm, &u));
375*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSCheckFromOptions(ts, u));
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SetInitialConditions(ts, u));
377*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) u, "temperature"));
378*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts, u));
379*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMTSCheckFromOptions(ts, u));
380*5f80ce2aSJacob Faibussowitsch   if (ctx.expTS) CHKERRQ(DMTSDestroyRHSMassMatrix(dm));
381c4762a1bSJed Brown 
382*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
385c4762a1bSJed Brown   ierr = PetscFinalize();
386c4762a1bSJed Brown   return ierr;
387c4762a1bSJed Brown }
388c4762a1bSJed Brown 
389c4762a1bSJed Brown /*TEST
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   test:
392a3d0cf85SMatthew G. Knepley     suffix: 2d_p1
393c4762a1bSJed Brown     requires: triangle
394a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -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
396c4762a1bSJed Brown   test:
397742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp
398742ee2edSMatthew G. Knepley     requires: triangle
399742ee2edSMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \
400742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error
401742ee2edSMatthew G. Knepley   test:
402a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
403a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_sconv
404c4762a1bSJed Brown     requires: triangle
405a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
406a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
407c4762a1bSJed Brown   test:
408742ee2edSMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1]
409742ee2edSMatthew G. Knepley     suffix: 2d_p1_sconv_2
410742ee2edSMatthew G. Knepley     requires: triangle
411742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
412742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu
413742ee2edSMatthew G. Knepley   test:
414a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
415a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_tconv
416c4762a1bSJed Brown     requires: triangle
417a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
418a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
419c4762a1bSJed Brown   test:
420742ee2edSMatthew G. Knepley     # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0]
421742ee2edSMatthew G. Knepley     suffix: 2d_p1_tconv_2
422742ee2edSMatthew G. Knepley     requires: triangle
423742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
424742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
425742ee2edSMatthew G. Knepley   test:
426742ee2edSMatthew G. Knepley     # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid
427742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp_tconv_2
428742ee2edSMatthew G. Knepley     requires: triangle
429742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
430742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu
431742ee2edSMatthew G. Knepley   test:
432742ee2edSMatthew G. Knepley     # The L_2 convergence rate cannot be seen since stability of the explicit integrator requires that is be more accurate than the grid
433742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp_tconv_2_lump
434742ee2edSMatthew G. Knepley     requires: triangle
435742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
436742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4
437742ee2edSMatthew G. Knepley   test:
438a3d0cf85SMatthew G. Knepley     suffix: 2d_p2
439c4762a1bSJed Brown     requires: triangle
440a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
441a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
442c4762a1bSJed Brown   test:
443a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
444a3d0cf85SMatthew G. Knepley     suffix: 2d_p2_sconv
445a3d0cf85SMatthew G. Knepley     requires: triangle
446a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
447a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
448c4762a1bSJed Brown   test:
449742ee2edSMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1]
450742ee2edSMatthew G. Knepley     suffix: 2d_p2_sconv_2
451742ee2edSMatthew G. Knepley     requires: triangle
452742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
453742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
454742ee2edSMatthew G. Knepley   test:
455a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
456a3d0cf85SMatthew G. Knepley     suffix: 2d_p2_tconv
457a3d0cf85SMatthew G. Knepley     requires: triangle
458a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
459a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
460c4762a1bSJed Brown   test:
461742ee2edSMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
462742ee2edSMatthew G. Knepley     suffix: 2d_p2_tconv_2
463742ee2edSMatthew G. Knepley     requires: triangle
464742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
465742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
466742ee2edSMatthew G. Knepley   test:
467a3d0cf85SMatthew G. Knepley     suffix: 2d_q1
46830602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
469a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
470c4762a1bSJed Brown   test:
471a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
472a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_sconv
47330602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
474a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
475c4762a1bSJed Brown   test:
476a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
477a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_tconv
47830602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
479a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
480a3d0cf85SMatthew G. Knepley   test:
481a3d0cf85SMatthew G. Knepley     suffix: 2d_q2
48230602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
483a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
484a3d0cf85SMatthew G. Knepley   test:
485a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
486a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_sconv
48730602db0SMatthew G. Knepley     args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
488a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
489a3d0cf85SMatthew G. Knepley   test:
490a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
491a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_tconv
49230602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
493a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
494a3d0cf85SMatthew G. Knepley 
495a3d0cf85SMatthew G. Knepley   test:
496a3d0cf85SMatthew G. Knepley     suffix: 3d_p1
497c4762a1bSJed Brown     requires: ctetgen
498a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
499a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
500c4762a1bSJed Brown   test:
501a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
502a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_sconv
503c4762a1bSJed Brown     requires: ctetgen
504a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
505a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
506c4762a1bSJed Brown   test:
507a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
508a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_tconv
509c4762a1bSJed Brown     requires: ctetgen
510a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
511a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
512c4762a1bSJed Brown   test:
513a3d0cf85SMatthew G. Knepley     suffix: 3d_p2
514c4762a1bSJed Brown     requires: ctetgen
515a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
516a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
517c4762a1bSJed Brown   test:
518a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
519a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_sconv
520a3d0cf85SMatthew G. Knepley     requires: ctetgen
521a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
522a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
523c4762a1bSJed Brown   test:
524a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
525a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_tconv
526a3d0cf85SMatthew G. Knepley     requires: ctetgen
527a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
528a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
529c4762a1bSJed Brown   test:
530a3d0cf85SMatthew G. Knepley     suffix: 3d_q1
53130602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
532a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
533c4762a1bSJed Brown   test:
534a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
535a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_sconv
53630602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
537a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
538a3d0cf85SMatthew G. Knepley   test:
539a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
540a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_tconv
54130602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
542a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
543a3d0cf85SMatthew G. Knepley   test:
544a3d0cf85SMatthew G. Knepley     suffix: 3d_q2
54530602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
546a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
547a3d0cf85SMatthew G. Knepley   test:
548a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
549a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_sconv
55030602db0SMatthew G. Knepley     args: -sol_type trig_linear -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
551a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
552a3d0cf85SMatthew G. Knepley   test:
553a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
554a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_tconv
55530602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
556a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
557a3d0cf85SMatthew G. Knepley 
558a3d0cf85SMatthew G. Knepley   test:
559a3d0cf85SMatthew 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
560a3d0cf85SMatthew G. Knepley     suffix: egads_sphere
561a3d0cf85SMatthew G. Knepley     requires: egads ctetgen
56230602db0SMatthew G. Knepley     args: -sol_type quadratic_linear \
56330602db0SMatthew G. Knepley           -dm_plex_boundary_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_plex_boundary_label marker -bd_dm_plex_scale 40 \
564a3d0cf85SMatthew G. Knepley           -temp_petscspace_degree 2 -dmts_check .0001 \
565a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
566c4762a1bSJed Brown 
567c4762a1bSJed Brown TEST*/
568