xref: /petsc/src/ts/tutorials/ex45.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
16*9371c9d4SSatish Balay typedef enum {
17*9371c9d4SSatish Balay   SOL_QUADRATIC_LINEAR,
18*9371c9d4SSatish Balay   SOL_QUADRATIC_TRIG,
19*9371c9d4SSatish Balay   SOL_TRIG_LINEAR,
20*9371c9d4SSatish Balay   SOL_TRIG_TRIG,
21*9371c9d4SSatish Balay   NUM_SOLUTION_TYPES
22*9371c9d4SSatish Balay } SolutionType;
23742ee2edSMatthew G. Knepley const char *solutionTypes[NUM_SOLUTION_TYPES + 1] = {"quadratic_linear", "quadratic_trig", "trig_linear", "trig_trig", "unknown"};
24a3d0cf85SMatthew G. Knepley 
25c4762a1bSJed Brown typedef struct {
26a3d0cf85SMatthew G. Knepley   SolutionType solType; /* Type of exact solution */
27742ee2edSMatthew G. Knepley   /* Solver setup */
28742ee2edSMatthew G. Knepley   PetscBool    expTS;  /* Flag for explicit timestepping */
29742ee2edSMatthew G. Knepley   PetscBool    lumped; /* Lump the mass matrix */
30c4762a1bSJed Brown } AppCtx;
31c4762a1bSJed Brown 
32a3d0cf85SMatthew G. Knepley /*
33a3d0cf85SMatthew G. Knepley Exact 2D solution:
34a3d0cf85SMatthew G. Knepley   u    = 2t + x^2 + y^2
35742ee2edSMatthew G. Knepley   u_t  = 2
36742ee2edSMatthew G. Knepley   \Delta u = 2 + 2 = 4
37742ee2edSMatthew G. Knepley   f    = 2
38a3d0cf85SMatthew G. Knepley   F(u) = 2 - (2 + 2) + 2 = 0
39a3d0cf85SMatthew G. Knepley 
40a3d0cf85SMatthew G. Knepley Exact 3D solution:
41a3d0cf85SMatthew G. Knepley   u = 3t + x^2 + y^2 + z^2
42a3d0cf85SMatthew G. Knepley   F(u) = 3 - (2 + 2 + 2) + 3 = 0
43a3d0cf85SMatthew G. Knepley */
44*9371c9d4SSatish Balay static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
45c4762a1bSJed Brown   PetscInt d;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   *u = dim * time;
48c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d] * x[d];
49c4762a1bSJed Brown   return 0;
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52*9371c9d4SSatish Balay static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
53a3d0cf85SMatthew G. Knepley   *u = dim;
54a3d0cf85SMatthew G. Knepley   return 0;
55a3d0cf85SMatthew G. Knepley }
56a3d0cf85SMatthew G. Knepley 
57*9371c9d4SSatish Balay static void f0_quad_lin_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
58742ee2edSMatthew G. Knepley   f0[0] = -(PetscScalar)dim;
59742ee2edSMatthew G. Knepley }
60*9371c9d4SSatish Balay static void f0_quad_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
61742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
62742ee2edSMatthew 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);
63742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66a3d0cf85SMatthew G. Knepley /*
67a3d0cf85SMatthew G. Knepley Exact 2D solution:
68a3d0cf85SMatthew G. Knepley   u = 2*cos(t) + x^2 + y^2
69a3d0cf85SMatthew G. Knepley   F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0
70a3d0cf85SMatthew G. Knepley 
71a3d0cf85SMatthew G. Knepley Exact 3D solution:
72a3d0cf85SMatthew G. Knepley   u = 3*cos(t) + x^2 + y^2 + z^2
73a3d0cf85SMatthew G. Knepley   F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0
74a3d0cf85SMatthew G. Knepley */
75*9371c9d4SSatish Balay static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
76a3d0cf85SMatthew G. Knepley   PetscInt d;
77a3d0cf85SMatthew G. Knepley 
78a3d0cf85SMatthew G. Knepley   *u = dim * PetscCosReal(time);
79a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += x[d] * x[d];
80a3d0cf85SMatthew G. Knepley   return 0;
81a3d0cf85SMatthew G. Knepley }
82a3d0cf85SMatthew G. Knepley 
83*9371c9d4SSatish Balay static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
84a3d0cf85SMatthew G. Knepley   *u = -dim * PetscSinReal(time);
85a3d0cf85SMatthew G. Knepley   return 0;
86a3d0cf85SMatthew G. Knepley }
87a3d0cf85SMatthew G. Knepley 
88*9371c9d4SSatish Balay static void f0_quad_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
89742ee2edSMatthew G. Knepley   f0[0] = -dim * (PetscSinReal(t) + 2.0);
90742ee2edSMatthew G. Knepley }
91*9371c9d4SSatish Balay static void f0_quad_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
92742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
93742ee2edSMatthew 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);
94742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
95a3d0cf85SMatthew G. Knepley }
96a3d0cf85SMatthew G. Knepley 
97a3d0cf85SMatthew G. Knepley /*
98a3d0cf85SMatthew G. Knepley Exact 2D solution:
99a3d0cf85SMatthew G. Knepley   u = 2\pi^2 t + cos(\pi x) + cos(\pi y)
100a3d0cf85SMatthew 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
101a3d0cf85SMatthew G. Knepley 
102a3d0cf85SMatthew G. Knepley Exact 3D solution:
103a3d0cf85SMatthew G. Knepley   u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z)
104a3d0cf85SMatthew 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
105a3d0cf85SMatthew G. Knepley */
106*9371c9d4SSatish Balay static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
107a3d0cf85SMatthew G. Knepley   PetscInt d;
108a3d0cf85SMatthew G. Knepley 
109a3d0cf85SMatthew G. Knepley   *u = dim * PetscSqr(PETSC_PI) * time;
110a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]);
111a3d0cf85SMatthew G. Knepley   return 0;
112a3d0cf85SMatthew G. Knepley }
113a3d0cf85SMatthew G. Knepley 
114*9371c9d4SSatish Balay static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
115a3d0cf85SMatthew G. Knepley   *u = dim * PetscSqr(PETSC_PI);
116a3d0cf85SMatthew G. Knepley   return 0;
117a3d0cf85SMatthew G. Knepley }
118a3d0cf85SMatthew G. Knepley 
119*9371c9d4SSatish Balay static void f0_trig_lin(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
120a3d0cf85SMatthew G. Knepley   PetscInt d;
121a3d0cf85SMatthew G. Knepley   f0[0] = u_t[0];
122a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * (PetscCosReal(PETSC_PI * x[d]) - 1.0);
123a3d0cf85SMatthew G. Knepley }
124a3d0cf85SMatthew G. Knepley 
125742ee2edSMatthew G. Knepley /*
126742ee2edSMatthew G. Knepley Exact 2D solution:
127742ee2edSMatthew G. Knepley   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y)
128742ee2edSMatthew G. Knepley   u_t  = -pi^2 sin(t)
129742ee2edSMatthew G. Knepley   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y))
130742ee2edSMatthew G. Knepley   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y))
131742ee2edSMatthew 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
132742ee2edSMatthew G. Knepley 
133742ee2edSMatthew G. Knepley Exact 3D solution:
134742ee2edSMatthew G. Knepley   u    = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z)
135742ee2edSMatthew G. Knepley   u_t  = -pi^2 sin(t)
136742ee2edSMatthew G. Knepley   \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
137742ee2edSMatthew G. Knepley   f    = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
138742ee2edSMatthew 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
139742ee2edSMatthew G. Knepley */
140*9371c9d4SSatish Balay static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
141742ee2edSMatthew G. Knepley   PetscInt d;
142742ee2edSMatthew G. Knepley 
143742ee2edSMatthew G. Knepley   *u = PetscSqr(PETSC_PI) * PetscCosReal(time);
144742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]);
145742ee2edSMatthew G. Knepley   return 0;
146742ee2edSMatthew G. Knepley }
147742ee2edSMatthew G. Knepley 
148*9371c9d4SSatish Balay static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
149742ee2edSMatthew G. Knepley   *u = -PetscSqr(PETSC_PI) * PetscSinReal(time);
150742ee2edSMatthew G. Knepley   return 0;
151742ee2edSMatthew G. Knepley }
152742ee2edSMatthew G. Knepley 
153*9371c9d4SSatish Balay static void f0_trig_trig_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
154742ee2edSMatthew G. Knepley   PetscInt d;
155742ee2edSMatthew G. Knepley   f0[0] -= PetscSqr(PETSC_PI) * PetscSinReal(t);
156742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * PetscCosReal(PETSC_PI * x[d]);
157742ee2edSMatthew G. Knepley }
158*9371c9d4SSatish Balay static void f0_trig_trig(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) {
159742ee2edSMatthew G. Knepley   PetscScalar exp[1] = {0.};
160742ee2edSMatthew 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);
161742ee2edSMatthew G. Knepley   f0[0] = u_t[0] - exp[0];
162742ee2edSMatthew G. Knepley }
163742ee2edSMatthew G. Knepley 
164*9371c9d4SSatish Balay static void f1_temp_exp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
165742ee2edSMatthew G. Knepley   PetscInt d;
166742ee2edSMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = -u_x[d];
167742ee2edSMatthew G. Knepley }
168*9371c9d4SSatish Balay static void f1_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) {
169c4762a1bSJed Brown   PetscInt d;
170a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) f1[d] = u_x[d];
171c4762a1bSJed Brown }
172c4762a1bSJed Brown 
173*9371c9d4SSatish Balay static void g3_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) {
174c4762a1bSJed Brown   PetscInt d;
175a3d0cf85SMatthew G. Knepley   for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178*9371c9d4SSatish Balay static void g0_temp(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) {
179c4762a1bSJed Brown   g0[0] = u_tShift * 1.0;
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
183a3d0cf85SMatthew G. Knepley   PetscInt sol;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   PetscFunctionBeginUser;
186a3d0cf85SMatthew G. Knepley   options->solType = SOL_QUADRATIC_LINEAR;
187742ee2edSMatthew G. Knepley   options->expTS   = PETSC_FALSE;
188742ee2edSMatthew G. Knepley   options->lumped  = PETSC_TRUE;
189c4762a1bSJed Brown 
190d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");
1919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
192a3d0cf85SMatthew G. Knepley   options->solType = (SolutionType)sol;
1939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL));
1949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL));
195d0609cedSBarry Smith   PetscOptionsEnd();
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) {
200c4762a1bSJed Brown   PetscFunctionBeginUser;
2019566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2029566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
2039566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
2049566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) {
209a3d0cf85SMatthew G. Knepley   PetscDS        ds;
21045480ffeSMatthew G. Knepley   DMLabel        label;
211c4762a1bSJed Brown   const PetscInt id = 1;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   PetscFunctionBeginUser;
2149566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
2159566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2169566063dSJacob Faibussowitsch   PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp));
217a3d0cf85SMatthew G. Knepley   switch (ctx->solType) {
218a3d0cf85SMatthew G. Knepley   case SOL_QUADRATIC_LINEAR:
2199566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp));
2209566063dSJacob Faibussowitsch     PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp));
2219566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx));
2229566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx));
2239566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms_quad_lin, (void (*)(void))mms_quad_lin_t, ctx, NULL));
224a3d0cf85SMatthew G. Knepley     break;
225a3d0cf85SMatthew G. Knepley   case SOL_QUADRATIC_TRIG:
2269566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp));
2279566063dSJacob Faibussowitsch     PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp));
2289566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx));
2299566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx));
2309566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms_quad_trig, (void (*)(void))mms_quad_trig_t, ctx, NULL));
231a3d0cf85SMatthew G. Knepley     break;
232a3d0cf85SMatthew G. Knepley   case SOL_TRIG_LINEAR:
2339566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp));
2349566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx));
2359566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx));
2369566063dSJacob Faibussowitsch     PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))mms_trig_lin, (void (*)(void))mms_trig_lin_t, ctx, NULL));
237a3d0cf85SMatthew G. Knepley     break;
238742ee2edSMatthew G. Knepley   case SOL_TRIG_TRIG:
2399566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp));
2409566063dSJacob Faibussowitsch     PetscCall(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp));
2419566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx));
2429566063dSJacob Faibussowitsch     PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx));
243742ee2edSMatthew G. Knepley     break;
24463a3b9bcSJacob Faibussowitsch   default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
245a3d0cf85SMatthew G. Knepley   }
246c4762a1bSJed Brown   PetscFunctionReturn(0);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) {
250c4762a1bSJed Brown   DM             cdm = dm;
251c4762a1bSJed Brown   PetscFE        fe;
252a3d0cf85SMatthew G. Knepley   DMPolytopeType ct;
253a3d0cf85SMatthew G. Knepley   PetscBool      simplex;
254a3d0cf85SMatthew G. Knepley   PetscInt       dim, cStart;
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   PetscFunctionBeginUser;
2579566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
2599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCellType(dm, cStart, &ct));
260a3d0cf85SMatthew G. Knepley   simplex = DMPolytopeTypeGetNumVertices(ct) == DMPolytopeTypeGetDim(ct) + 1 ? PETSC_TRUE : PETSC_FALSE;
261c4762a1bSJed Brown   /* Create finite element */
2629566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe));
2639566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "temperature"));
264c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
2659566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2669566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
267742ee2edSMatthew G. Knepley   if (ctx->expTS) {
268742ee2edSMatthew G. Knepley     PetscDS ds;
269742ee2edSMatthew G. Knepley 
2709566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
2719566063dSJacob Faibussowitsch     PetscCall(PetscDSSetImplicit(ds, 0, PETSC_FALSE));
272742ee2edSMatthew G. Knepley   }
2739566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, ctx));
274c4762a1bSJed Brown   while (cdm) {
2759566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
2769566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
277c4762a1bSJed Brown   }
2789566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
279c4762a1bSJed Brown   PetscFunctionReturn(0);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282*9371c9d4SSatish Balay static PetscErrorCode SetInitialConditions(TS ts, Vec u) {
283a3d0cf85SMatthew G. Knepley   DM        dm;
284a3d0cf85SMatthew G. Knepley   PetscReal t;
285a3d0cf85SMatthew G. Knepley 
2867510d9b0SBarry Smith   PetscFunctionBeginUser;
2879566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2889566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
2899566063dSJacob Faibussowitsch   PetscCall(DMComputeExactSolution(dm, t, u, NULL));
290a3d0cf85SMatthew G. Knepley   PetscFunctionReturn(0);
291a3d0cf85SMatthew G. Knepley }
292a3d0cf85SMatthew G. Knepley 
293*9371c9d4SSatish Balay int main(int argc, char **argv) {
294c4762a1bSJed Brown   DM     dm;
295c4762a1bSJed Brown   TS     ts;
296a3d0cf85SMatthew G. Knepley   Vec    u;
297a3d0cf85SMatthew G. Knepley   AppCtx ctx;
298c4762a1bSJed Brown 
299327415f7SBarry Smith   PetscFunctionBeginUser;
3009566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3019566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
3029566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
3039566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &ctx));
3049566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &ctx));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3079566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
3089566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
309742ee2edSMatthew G. Knepley   if (ctx.expTS) {
3109566063dSJacob Faibussowitsch     PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx));
3119566063dSJacob Faibussowitsch     if (ctx.lumped) PetscCall(DMTSCreateRHSMassMatrixLumped(dm));
3129566063dSJacob Faibussowitsch     else PetscCall(DMTSCreateRHSMassMatrix(dm));
313742ee2edSMatthew G. Knepley   } else {
3149566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
3159566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
316742ee2edSMatthew G. Knepley   }
3179566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
3189566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
3199566063dSJacob Faibussowitsch   PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions));
320c4762a1bSJed Brown 
3219566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
3229566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
3239566063dSJacob Faibussowitsch   PetscCall(SetInitialConditions(ts, u));
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "temperature"));
3259566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
3269566063dSJacob Faibussowitsch   PetscCall(DMTSCheckFromOptions(ts, u));
3279566063dSJacob Faibussowitsch   if (ctx.expTS) PetscCall(DMTSDestroyRHSMassMatrix(dm));
328c4762a1bSJed Brown 
3299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3309566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3319566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
333b122ec5aSJacob Faibussowitsch   return 0;
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown /*TEST
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   test:
339a3d0cf85SMatthew G. Knepley     suffix: 2d_p1
340c4762a1bSJed Brown     requires: triangle
341a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
342a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
343c4762a1bSJed Brown   test:
344742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp
345742ee2edSMatthew G. Knepley     requires: triangle
346742ee2edSMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \
347742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-3 -ts_monitor_error
348742ee2edSMatthew G. Knepley   test:
349a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
350a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_sconv
351c4762a1bSJed Brown     requires: triangle
352a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
353a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
354c4762a1bSJed Brown   test:
355742ee2edSMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1]
356742ee2edSMatthew G. Knepley     suffix: 2d_p1_sconv_2
357742ee2edSMatthew G. Knepley     requires: triangle
358742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
359742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 1e-6 -snes_error_if_not_converged -pc_type lu
360742ee2edSMatthew G. Knepley   test:
361a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
362a3d0cf85SMatthew G. Knepley     suffix: 2d_p1_tconv
363c4762a1bSJed Brown     requires: triangle
364a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
365a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
366c4762a1bSJed Brown   test:
367742ee2edSMatthew G. Knepley     # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0]
368742ee2edSMatthew G. Knepley     suffix: 2d_p1_tconv_2
369742ee2edSMatthew G. Knepley     requires: triangle
370742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
371742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
372742ee2edSMatthew G. Knepley   test:
373742ee2edSMatthew 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
374742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp_tconv_2
375742ee2edSMatthew G. Knepley     requires: triangle
376742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
377742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4 -lumped 0 -mass_pc_type lu
378742ee2edSMatthew G. Knepley   test:
379742ee2edSMatthew 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
380742ee2edSMatthew G. Knepley     suffix: 2d_p1_exp_tconv_2_lump
381742ee2edSMatthew G. Knepley     requires: triangle
382742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
383742ee2edSMatthew G. Knepley           -ts_type euler -ts_max_steps 4 -ts_dt 1e-4
384742ee2edSMatthew G. Knepley   test:
385a3d0cf85SMatthew G. Knepley     suffix: 2d_p2
386c4762a1bSJed Brown     requires: triangle
387a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
388a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
389c4762a1bSJed Brown   test:
390a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
391a3d0cf85SMatthew G. Knepley     suffix: 2d_p2_sconv
392a3d0cf85SMatthew G. Knepley     requires: triangle
393a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
394a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
395c4762a1bSJed Brown   test:
396742ee2edSMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1]
397742ee2edSMatthew G. Knepley     suffix: 2d_p2_sconv_2
398742ee2edSMatthew G. Knepley     requires: triangle
399742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
400742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
401742ee2edSMatthew 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_p2_tconv
404a3d0cf85SMatthew G. Knepley     requires: triangle
405a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
406a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
407c4762a1bSJed Brown   test:
408742ee2edSMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
409742ee2edSMatthew G. Knepley     suffix: 2d_p2_tconv_2
410742ee2edSMatthew G. Knepley     requires: triangle
411742ee2edSMatthew G. Knepley     args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
412742ee2edSMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
413742ee2edSMatthew G. Knepley   test:
414a3d0cf85SMatthew G. Knepley     suffix: 2d_q1
41530602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
416a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
417c4762a1bSJed Brown   test:
418a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
419a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_sconv
42030602db0SMatthew 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 \
421a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
422c4762a1bSJed Brown   test:
423a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
424a3d0cf85SMatthew G. Knepley     suffix: 2d_q1_tconv
42530602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
426a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
427a3d0cf85SMatthew G. Knepley   test:
428a3d0cf85SMatthew G. Knepley     suffix: 2d_q2
42930602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
430a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
431a3d0cf85SMatthew G. Knepley   test:
432a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
433a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_sconv
43430602db0SMatthew 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 \
435a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
436a3d0cf85SMatthew G. Knepley   test:
437a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
438a3d0cf85SMatthew G. Knepley     suffix: 2d_q2_tconv
43930602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -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
441a3d0cf85SMatthew G. Knepley 
442a3d0cf85SMatthew G. Knepley   test:
443a3d0cf85SMatthew G. Knepley     suffix: 3d_p1
444c4762a1bSJed Brown     requires: ctetgen
445a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
446a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
447c4762a1bSJed Brown   test:
448a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
449a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_sconv
450c4762a1bSJed Brown     requires: ctetgen
451a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
452a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
453c4762a1bSJed Brown   test:
454a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
455a3d0cf85SMatthew G. Knepley     suffix: 3d_p1_tconv
456c4762a1bSJed Brown     requires: ctetgen
457a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
458a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
459c4762a1bSJed Brown   test:
460a3d0cf85SMatthew G. Knepley     suffix: 3d_p2
461c4762a1bSJed Brown     requires: ctetgen
462a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
463a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
464c4762a1bSJed Brown   test:
465a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
466a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_sconv
467a3d0cf85SMatthew G. Knepley     requires: ctetgen
468a3d0cf85SMatthew G. Knepley     args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
469a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
470c4762a1bSJed Brown   test:
471a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
472a3d0cf85SMatthew G. Knepley     suffix: 3d_p2_tconv
473a3d0cf85SMatthew G. Knepley     requires: ctetgen
474a3d0cf85SMatthew G. Knepley     args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
475a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
476c4762a1bSJed Brown   test:
477a3d0cf85SMatthew G. Knepley     suffix: 3d_q1
47830602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
479a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
480c4762a1bSJed Brown   test:
481a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
482a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_sconv
48330602db0SMatthew 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 \
484a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00001 -snes_error_if_not_converged -pc_type lu
485a3d0cf85SMatthew G. Knepley   test:
486a3d0cf85SMatthew G. Knepley     # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
487a3d0cf85SMatthew G. Knepley     suffix: 3d_q1_tconv
48830602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
489a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
490a3d0cf85SMatthew G. Knepley   test:
491a3d0cf85SMatthew G. Knepley     suffix: 3d_q2
49230602db0SMatthew G. Knepley     args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
493a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
494a3d0cf85SMatthew G. Knepley   test:
495a3d0cf85SMatthew G. Knepley     # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
496a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_sconv
49730602db0SMatthew 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 \
498a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 1 -ts_dt 0.00000001 -snes_error_if_not_converged -pc_type lu
499a3d0cf85SMatthew G. Knepley   test:
500a3d0cf85SMatthew G. Knepley     # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
501a3d0cf85SMatthew G. Knepley     suffix: 3d_q2_tconv
50230602db0SMatthew G. Knepley     args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
503a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 4 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
504a3d0cf85SMatthew G. Knepley 
505a3d0cf85SMatthew G. Knepley   test:
506a3d0cf85SMatthew 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
507a3d0cf85SMatthew G. Knepley     suffix: egads_sphere
508a3d0cf85SMatthew G. Knepley     requires: egads ctetgen
50930602db0SMatthew G. Knepley     args: -sol_type quadratic_linear \
51030602db0SMatthew 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 \
511a3d0cf85SMatthew G. Knepley           -temp_petscspace_degree 2 -dmts_check .0001 \
512a3d0cf85SMatthew G. Knepley           -ts_type beuler -ts_max_steps 5 -ts_dt 0.1 -snes_error_if_not_converged -pc_type lu
513c4762a1bSJed Brown 
514c4762a1bSJed Brown TEST*/
515