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
169371c9d4SSatish Balay typedef enum {
179371c9d4SSatish Balay SOL_QUADRATIC_LINEAR,
189371c9d4SSatish Balay SOL_QUADRATIC_TRIG,
199371c9d4SSatish Balay SOL_TRIG_LINEAR,
209371c9d4SSatish Balay SOL_TRIG_TRIG,
219371c9d4SSatish Balay NUM_SOLUTION_TYPES
229371c9d4SSatish 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 */
30d9255fafSStefano Zampini PetscInt remesh_every; /* Remesh every number of steps */
31d9255fafSStefano Zampini DM remesh_dm; /* New DM after remeshing */
32c4762a1bSJed Brown } AppCtx;
33c4762a1bSJed Brown
34a3d0cf85SMatthew G. Knepley /*
35a3d0cf85SMatthew G. Knepley Exact 2D solution:
36a3d0cf85SMatthew G. Knepley u = 2t + x^2 + y^2
37742ee2edSMatthew G. Knepley u_t = 2
38742ee2edSMatthew G. Knepley \Delta u = 2 + 2 = 4
39742ee2edSMatthew G. Knepley f = 2
40a3d0cf85SMatthew G. Knepley F(u) = 2 - (2 + 2) + 2 = 0
41a3d0cf85SMatthew G. Knepley
42a3d0cf85SMatthew G. Knepley Exact 3D solution:
43a3d0cf85SMatthew G. Knepley u = 3t + x^2 + y^2 + z^2
44a3d0cf85SMatthew G. Knepley F(u) = 3 - (2 + 2 + 2) + 3 = 0
45a3d0cf85SMatthew G. Knepley */
mms_quad_lin(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)46*2a8381b2SBarry Smith static PetscErrorCode mms_quad_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
47d71ae5a4SJacob Faibussowitsch {
48c4762a1bSJed Brown PetscInt d;
49c4762a1bSJed Brown
50c4762a1bSJed Brown *u = dim * time;
51c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += x[d] * x[d];
523ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown
mms_quad_lin_t(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)55*2a8381b2SBarry Smith static PetscErrorCode mms_quad_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
56d71ae5a4SJacob Faibussowitsch {
57a3d0cf85SMatthew G. Knepley *u = dim;
583ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
59a3d0cf85SMatthew G. Knepley }
60a3d0cf85SMatthew G. Knepley
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[])61d71ae5a4SJacob Faibussowitsch 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[])
62d71ae5a4SJacob Faibussowitsch {
63742ee2edSMatthew G. Knepley f0[0] = -(PetscScalar)dim;
64742ee2edSMatthew G. Knepley }
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[])65d71ae5a4SJacob Faibussowitsch 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[])
66d71ae5a4SJacob Faibussowitsch {
67742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.};
68742ee2edSMatthew 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);
69742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0];
70c4762a1bSJed Brown }
71c4762a1bSJed Brown
72a3d0cf85SMatthew G. Knepley /*
73a3d0cf85SMatthew G. Knepley Exact 2D solution:
74a3d0cf85SMatthew G. Knepley u = 2*cos(t) + x^2 + y^2
75a3d0cf85SMatthew G. Knepley F(u) = -2*sint(t) - (2 + 2) + 2*sin(t) + 4 = 0
76a3d0cf85SMatthew G. Knepley
77a3d0cf85SMatthew G. Knepley Exact 3D solution:
78a3d0cf85SMatthew G. Knepley u = 3*cos(t) + x^2 + y^2 + z^2
79a3d0cf85SMatthew G. Knepley F(u) = -3*sin(t) - (2 + 2 + 2) + 3*sin(t) + 6 = 0
80a3d0cf85SMatthew G. Knepley */
mms_quad_trig(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)81*2a8381b2SBarry Smith static PetscErrorCode mms_quad_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
82d71ae5a4SJacob Faibussowitsch {
83a3d0cf85SMatthew G. Knepley PetscInt d;
84a3d0cf85SMatthew G. Knepley
85a3d0cf85SMatthew G. Knepley *u = dim * PetscCosReal(time);
86a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) *u += x[d] * x[d];
873ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
88a3d0cf85SMatthew G. Knepley }
89a3d0cf85SMatthew G. Knepley
mms_quad_trig_t(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)90*2a8381b2SBarry Smith static PetscErrorCode mms_quad_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
91d71ae5a4SJacob Faibussowitsch {
92a3d0cf85SMatthew G. Knepley *u = -dim * PetscSinReal(time);
933ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
94a3d0cf85SMatthew G. Knepley }
95a3d0cf85SMatthew G. Knepley
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[])96d71ae5a4SJacob Faibussowitsch 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[])
97d71ae5a4SJacob Faibussowitsch {
98742ee2edSMatthew G. Knepley f0[0] = -dim * (PetscSinReal(t) + 2.0);
99742ee2edSMatthew G. Knepley }
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[])100d71ae5a4SJacob Faibussowitsch 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[])
101d71ae5a4SJacob Faibussowitsch {
102742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.};
103742ee2edSMatthew 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);
104742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0];
105a3d0cf85SMatthew G. Knepley }
106a3d0cf85SMatthew G. Knepley
107a3d0cf85SMatthew G. Knepley /*
108a3d0cf85SMatthew G. Knepley Exact 2D solution:
109a3d0cf85SMatthew G. Knepley u = 2\pi^2 t + cos(\pi x) + cos(\pi y)
110a3d0cf85SMatthew 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
111a3d0cf85SMatthew G. Knepley
112a3d0cf85SMatthew G. Knepley Exact 3D solution:
113a3d0cf85SMatthew G. Knepley u = 3\pi^2 t + cos(\pi x) + cos(\pi y) + cos(\pi z)
114a3d0cf85SMatthew 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
115a3d0cf85SMatthew G. Knepley */
mms_trig_lin(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)116*2a8381b2SBarry Smith static PetscErrorCode mms_trig_lin(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
117d71ae5a4SJacob Faibussowitsch {
118a3d0cf85SMatthew G. Knepley PetscInt d;
119a3d0cf85SMatthew G. Knepley
120a3d0cf85SMatthew G. Knepley *u = dim * PetscSqr(PETSC_PI) * time;
121a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]);
1223ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
123a3d0cf85SMatthew G. Knepley }
124a3d0cf85SMatthew G. Knepley
mms_trig_lin_t(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)125*2a8381b2SBarry Smith static PetscErrorCode mms_trig_lin_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
126d71ae5a4SJacob Faibussowitsch {
127a3d0cf85SMatthew G. Knepley *u = dim * PetscSqr(PETSC_PI);
1283ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
129a3d0cf85SMatthew G. Knepley }
130a3d0cf85SMatthew G. Knepley
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[])131d71ae5a4SJacob Faibussowitsch 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[])
132d71ae5a4SJacob Faibussowitsch {
133a3d0cf85SMatthew G. Knepley PetscInt d;
134a3d0cf85SMatthew G. Knepley f0[0] = u_t[0];
135a3d0cf85SMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * (PetscCosReal(PETSC_PI * x[d]) - 1.0);
136a3d0cf85SMatthew G. Knepley }
137a3d0cf85SMatthew G. Knepley
138742ee2edSMatthew G. Knepley /*
139742ee2edSMatthew G. Knepley Exact 2D solution:
140742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y)
141742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t)
142742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y))
143742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y))
144742ee2edSMatthew 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
145742ee2edSMatthew G. Knepley
146742ee2edSMatthew G. Knepley Exact 3D solution:
147742ee2edSMatthew G. Knepley u = pi^2 cos(t) + cos(\pi x) + cos(\pi y) + cos(\pi z)
148742ee2edSMatthew G. Knepley u_t = -pi^2 sin(t)
149742ee2edSMatthew G. Knepley \Delta u = -\pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
150742ee2edSMatthew G. Knepley f = pi^2 sin(t) - \pi^2 (cos(\pi x) + cos(\pi y) + cos(\pi z))
151742ee2edSMatthew 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
152742ee2edSMatthew G. Knepley */
mms_trig_trig(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)153*2a8381b2SBarry Smith static PetscErrorCode mms_trig_trig(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
154d71ae5a4SJacob Faibussowitsch {
155742ee2edSMatthew G. Knepley PetscInt d;
156742ee2edSMatthew G. Knepley
157742ee2edSMatthew G. Knepley *u = PetscSqr(PETSC_PI) * PetscCosReal(time);
158742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) *u += PetscCosReal(PETSC_PI * x[d]);
1593ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
160742ee2edSMatthew G. Knepley }
161742ee2edSMatthew G. Knepley
mms_trig_trig_t(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)162*2a8381b2SBarry Smith static PetscErrorCode mms_trig_trig_t(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
163d71ae5a4SJacob Faibussowitsch {
164742ee2edSMatthew G. Knepley *u = -PetscSqr(PETSC_PI) * PetscSinReal(time);
1653ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
166742ee2edSMatthew G. Knepley }
167742ee2edSMatthew G. Knepley
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[])168d71ae5a4SJacob Faibussowitsch 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[])
169d71ae5a4SJacob Faibussowitsch {
170742ee2edSMatthew G. Knepley PetscInt d;
171742ee2edSMatthew G. Knepley f0[0] -= PetscSqr(PETSC_PI) * PetscSinReal(t);
172742ee2edSMatthew G. Knepley for (d = 0; d < dim; ++d) f0[0] += PetscSqr(PETSC_PI) * PetscCosReal(PETSC_PI * x[d]);
173742ee2edSMatthew G. Knepley }
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[])174d71ae5a4SJacob Faibussowitsch 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[])
175d71ae5a4SJacob Faibussowitsch {
176742ee2edSMatthew G. Knepley PetscScalar exp[1] = {0.};
177742ee2edSMatthew 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);
178742ee2edSMatthew G. Knepley f0[0] = u_t[0] - exp[0];
179742ee2edSMatthew G. Knepley }
180742ee2edSMatthew G. Knepley
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[])181d71ae5a4SJacob Faibussowitsch 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[])
182d71ae5a4SJacob Faibussowitsch {
1837279ac95SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f1[d] = -u_x[d];
184742ee2edSMatthew G. Knepley }
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[])185d71ae5a4SJacob Faibussowitsch 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[])
186d71ae5a4SJacob Faibussowitsch {
1877279ac95SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f1[d] = u_x[d];
188c4762a1bSJed Brown }
189c4762a1bSJed Brown
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[])190d71ae5a4SJacob Faibussowitsch 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[])
191d71ae5a4SJacob Faibussowitsch {
1927279ac95SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = 1.0;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
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[])195d71ae5a4SJacob Faibussowitsch 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[])
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown g0[0] = u_tShift * 1.0;
198c4762a1bSJed Brown }
199c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)200d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
201d71ae5a4SJacob Faibussowitsch {
202a3d0cf85SMatthew G. Knepley PetscInt sol;
203c4762a1bSJed Brown
204c4762a1bSJed Brown PetscFunctionBeginUser;
205a3d0cf85SMatthew G. Knepley options->solType = SOL_QUADRATIC_LINEAR;
206742ee2edSMatthew G. Knepley options->expTS = PETSC_FALSE;
207742ee2edSMatthew G. Knepley options->lumped = PETSC_TRUE;
208d9255fafSStefano Zampini options->remesh_every = 0;
209c4762a1bSJed Brown
210d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Heat Equation Options", "DMPLEX");
211d9255fafSStefano Zampini sol = options->solType;
2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-sol_type", "Type of exact solution", "ex45.c", solutionTypes, NUM_SOLUTION_TYPES, solutionTypes[options->solType], &sol, NULL));
213a3d0cf85SMatthew G. Knepley options->solType = (SolutionType)sol;
2149566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-explicit", "Use explicit timestepping", "ex45.c", options->expTS, &options->expTS, NULL));
2159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-lumped", "Lump the mass matrix", "ex45.c", options->lumped, &options->lumped, NULL));
216d9255fafSStefano Zampini PetscCall(PetscOptionsInt("-remesh_every", "Remesh every number of steps", "ex45.c", options->remesh_every, &options->remesh_every, NULL));
217d0609cedSBarry Smith PetscOptionsEnd();
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
219c4762a1bSJed Brown }
220c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * ctx)221d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
222d71ae5a4SJacob Faibussowitsch {
223c4762a1bSJed Brown PetscFunctionBeginUser;
2249566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
2259566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
2269566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
227d9255fafSStefano Zampini {
228d9255fafSStefano Zampini char convType[256];
229d9255fafSStefano Zampini PetscBool flg;
230d9255fafSStefano Zampini PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX");
231d9255fafSStefano Zampini PetscCall(PetscOptionsFList("-dm_plex_convert_type", "Convert DMPlex to another format", __FILE__, DMList, DMPLEX, convType, 256, &flg));
232d9255fafSStefano Zampini PetscOptionsEnd();
233d9255fafSStefano Zampini if (flg) {
234d9255fafSStefano Zampini DM dmConv;
235d9255fafSStefano Zampini PetscCall(DMConvert(*dm, convType, &dmConv));
236d9255fafSStefano Zampini if (dmConv) {
237d9255fafSStefano Zampini PetscCall(DMDestroy(dm));
238d9255fafSStefano Zampini *dm = dmConv;
239d9255fafSStefano Zampini PetscCall(DMSetFromOptions(*dm));
240d9255fafSStefano Zampini PetscCall(DMSetUp(*dm));
241d9255fafSStefano Zampini }
242d9255fafSStefano Zampini }
243d9255fafSStefano Zampini }
2449566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
246c4762a1bSJed Brown }
247c4762a1bSJed Brown
SetupProblem(DM dm,AppCtx * ctx)248d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
249d71ae5a4SJacob Faibussowitsch {
250a3d0cf85SMatthew G. Knepley PetscDS ds;
25145480ffeSMatthew G. Knepley DMLabel label;
252c4762a1bSJed Brown const PetscInt id = 1;
253c4762a1bSJed Brown
254c4762a1bSJed Brown PetscFunctionBeginUser;
2559566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
2569566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2579566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_temp, NULL, NULL, g3_temp));
258a3d0cf85SMatthew G. Knepley switch (ctx->solType) {
259a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_LINEAR:
2609566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_lin, f1_temp));
2619566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_lin_exp, f1_temp_exp));
2629566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_lin, ctx));
2639566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_lin_t, ctx));
26457d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms_quad_lin, (PetscVoidFn *)mms_quad_lin_t, ctx, NULL));
265a3d0cf85SMatthew G. Knepley break;
266a3d0cf85SMatthew G. Knepley case SOL_QUADRATIC_TRIG:
2679566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_quad_trig, f1_temp));
2689566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_quad_trig_exp, f1_temp_exp));
2699566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_quad_trig, ctx));
2709566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_quad_trig_t, ctx));
27157d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms_quad_trig, (PetscVoidFn *)mms_quad_trig_t, ctx, NULL));
272a3d0cf85SMatthew G. Knepley break;
273a3d0cf85SMatthew G. Knepley case SOL_TRIG_LINEAR:
2749566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_lin, f1_temp));
2759566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_lin, ctx));
2769566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_lin_t, ctx));
27757d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)mms_trig_lin, (PetscVoidFn *)mms_trig_lin_t, ctx, NULL));
278a3d0cf85SMatthew G. Knepley break;
279742ee2edSMatthew G. Knepley case SOL_TRIG_TRIG:
2809566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_trig_trig, f1_temp));
2819566063dSJacob Faibussowitsch PetscCall(PetscDSSetRHSResidual(ds, 0, f0_trig_trig_exp, f1_temp_exp));
2829566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, mms_trig_trig, ctx));
2839566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolutionTimeDerivative(ds, 0, mms_trig_trig_t, ctx));
284742ee2edSMatthew G. Knepley break;
285d71ae5a4SJacob Faibussowitsch default:
286d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid solution type: %s (%d)", solutionTypes[PetscMin(ctx->solType, NUM_SOLUTION_TYPES)], ctx->solType);
287a3d0cf85SMatthew G. Knepley }
2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown
SetupDiscretization(DM dm,AppCtx * ctx)291d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
292d71ae5a4SJacob Faibussowitsch {
293d9255fafSStefano Zampini DM plex, cdm = dm;
294c4762a1bSJed Brown PetscFE fe;
295d9255fafSStefano Zampini PetscBool simplex;
296d9255fafSStefano Zampini PetscInt dim;
297c4762a1bSJed Brown
298c4762a1bSJed Brown PetscFunctionBeginUser;
2999566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
300d9255fafSStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex));
301d9255fafSStefano Zampini PetscCall(DMPlexIsSimplex(plex, &simplex));
302d9255fafSStefano Zampini PetscCall(DMDestroy(&plex));
303d9255fafSStefano Zampini PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, "temp_", -1, &fe));
3049566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "temperature"));
305c4762a1bSJed Brown /* Set discretization and boundary conditions for each mesh */
3069566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
3079566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
308742ee2edSMatthew G. Knepley if (ctx->expTS) {
309742ee2edSMatthew G. Knepley PetscDS ds;
310742ee2edSMatthew G. Knepley
3119566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3129566063dSJacob Faibussowitsch PetscCall(PetscDSSetImplicit(ds, 0, PETSC_FALSE));
313742ee2edSMatthew G. Knepley }
3149566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx));
315c4762a1bSJed Brown while (cdm) {
3169566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
3179566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
318c4762a1bSJed Brown }
3199566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
321c4762a1bSJed Brown }
322c4762a1bSJed Brown
323d9255fafSStefano Zampini #include <petsc/private/dmpleximpl.h>
Remesh(DM dm,Vec U,DM * newdm)324d9255fafSStefano Zampini static PetscErrorCode Remesh(DM dm, Vec U, DM *newdm)
325d9255fafSStefano Zampini {
326d9255fafSStefano Zampini PetscFunctionBeginUser;
327d9255fafSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-remesh_dmin_view"));
328d9255fafSStefano Zampini *newdm = NULL;
329d9255fafSStefano Zampini
330d9255fafSStefano Zampini PetscInt dim;
331d9255fafSStefano Zampini DM plex;
332d9255fafSStefano Zampini PetscBool simplex;
333d9255fafSStefano Zampini PetscCall(DMGetDimension(dm, &dim));
334d9255fafSStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex));
335d9255fafSStefano Zampini PetscCall(DMPlexIsSimplex(plex, &simplex));
336d9255fafSStefano Zampini PetscCall(DMDestroy(&plex));
337d9255fafSStefano Zampini
338d9255fafSStefano Zampini DM dmGrad;
339d9255fafSStefano Zampini PetscFE fe;
340d9255fafSStefano Zampini PetscCall(DMClone(dm, &dmGrad));
341d9255fafSStefano Zampini PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, dim, simplex, "grad_", -1, &fe));
342d9255fafSStefano Zampini PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
343d9255fafSStefano Zampini PetscCall(PetscFEDestroy(&fe));
344d9255fafSStefano Zampini PetscCall(DMCreateDS(dmGrad));
345d9255fafSStefano Zampini
346d9255fafSStefano Zampini Vec locU, locG;
347d9255fafSStefano Zampini PetscCall(DMGetLocalVector(dm, &locU));
348d9255fafSStefano Zampini PetscCall(DMGetLocalVector(dmGrad, &locG));
349d9255fafSStefano Zampini PetscCall(DMGlobalToLocal(dm, U, INSERT_VALUES, locU));
350d9255fafSStefano Zampini PetscCall(VecViewFromOptions(locU, NULL, "-remesh_input_view"));
351d9255fafSStefano Zampini PetscCall(DMPlexComputeGradientClementInterpolant(dm, locU, locG));
352d9255fafSStefano Zampini PetscCall(VecViewFromOptions(locG, NULL, "-remesh_grad_view"));
353d9255fafSStefano Zampini
354d9255fafSStefano Zampini const PetscScalar *g;
355d9255fafSStefano Zampini PetscScalar *u;
356d9255fafSStefano Zampini PetscInt n;
357d9255fafSStefano Zampini PetscCall(VecGetLocalSize(locU, &n));
358d9255fafSStefano Zampini PetscCall(VecGetArrayWrite(locU, &u));
359d9255fafSStefano Zampini PetscCall(VecGetArrayRead(locG, &g));
360d9255fafSStefano Zampini for (PetscInt i = 0; i < n; i++) {
361d9255fafSStefano Zampini PetscReal norm = 0.0;
362d9255fafSStefano Zampini for (PetscInt d = 0; d < dim; d++) norm += PetscSqr(PetscRealPart(g[dim * i + d]));
363d9255fafSStefano Zampini u[i] = PetscSqrtReal(norm);
364d9255fafSStefano Zampini }
365d9255fafSStefano Zampini PetscCall(VecRestoreArrayRead(locG, &g));
366d9255fafSStefano Zampini PetscCall(VecRestoreArrayWrite(locU, &u));
367d9255fafSStefano Zampini
368d9255fafSStefano Zampini DM dmM;
369d9255fafSStefano Zampini Vec metric;
370d9255fafSStefano Zampini PetscCall(DMClone(dm, &dmM));
371d9255fafSStefano Zampini PetscCall(DMPlexMetricCreateIsotropic(dmM, 0, locU, &metric));
372d9255fafSStefano Zampini PetscCall(DMDestroy(&dmM));
373d9255fafSStefano Zampini PetscCall(DMRestoreLocalVector(dm, &locU));
374d9255fafSStefano Zampini PetscCall(DMRestoreLocalVector(dmGrad, &locG));
375d9255fafSStefano Zampini PetscCall(DMDestroy(&dmGrad));
376d9255fafSStefano Zampini
377d9255fafSStefano Zampini // TODO remove?
378d9255fafSStefano Zampini PetscScalar scale = 10.0;
379d9255fafSStefano Zampini PetscCall(PetscOptionsGetScalar(NULL, NULL, "-metric_scale", &scale, NULL));
380d9255fafSStefano Zampini PetscCall(VecScale(metric, scale));
381d9255fafSStefano Zampini PetscCall(VecViewFromOptions(metric, NULL, "-metric_view"));
382d9255fafSStefano Zampini
383d9255fafSStefano Zampini DMLabel label = NULL;
384d9255fafSStefano Zampini PetscCall(DMGetLabel(dm, "marker", &label));
385d9255fafSStefano Zampini PetscCall(DMAdaptMetric(dm, metric, label, NULL, newdm));
386d9255fafSStefano Zampini PetscCall(VecDestroy(&metric));
387d9255fafSStefano Zampini
388d9255fafSStefano Zampini PetscCall(DMViewFromOptions(*newdm, NULL, "-remesh_dmout_view"));
389d9255fafSStefano Zampini
390d9255fafSStefano Zampini AppCtx *ctx;
391d9255fafSStefano Zampini PetscCall(DMGetApplicationContext(dm, &ctx));
392d9255fafSStefano Zampini PetscCall(DMSetApplicationContext(*newdm, ctx));
393d9255fafSStefano Zampini PetscCall(SetupDiscretization(*newdm, ctx));
394d9255fafSStefano Zampini
395d9255fafSStefano Zampini // TODO
396d9255fafSStefano Zampini ((DM_Plex *)(*newdm)->data)->useHashLocation = ((DM_Plex *)dm->data)->useHashLocation;
397d9255fafSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
398d9255fafSStefano Zampini }
399d9255fafSStefano Zampini
SetInitialConditions(TS ts,Vec u)400d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetInitialConditions(TS ts, Vec u)
401d71ae5a4SJacob Faibussowitsch {
402a3d0cf85SMatthew G. Knepley DM dm;
403a3d0cf85SMatthew G. Knepley PetscReal t;
404a3d0cf85SMatthew G. Knepley
4057510d9b0SBarry Smith PetscFunctionBeginUser;
4069566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
4079566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t));
4089566063dSJacob Faibussowitsch PetscCall(DMComputeExactSolution(dm, t, u, NULL));
409d9255fafSStefano Zampini PetscCall(VecSetOptionsPrefix(u, NULL));
410d9255fafSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
411d9255fafSStefano Zampini }
412d9255fafSStefano Zampini
TransferSetUp(TS ts,PetscInt step,PetscReal time,Vec U,PetscBool * remesh,void * tctx)413d9255fafSStefano Zampini static PetscErrorCode TransferSetUp(TS ts, PetscInt step, PetscReal time, Vec U, PetscBool *remesh, void *tctx)
414d9255fafSStefano Zampini {
415d9255fafSStefano Zampini AppCtx *ctx = (AppCtx *)tctx;
416d9255fafSStefano Zampini
417d9255fafSStefano Zampini PetscFunctionBeginUser;
418d9255fafSStefano Zampini *remesh = PETSC_FALSE;
419d9255fafSStefano Zampini if (ctx->remesh_every > 0 && step && step % ctx->remesh_every == 0) {
420d9255fafSStefano Zampini DM dm;
421d9255fafSStefano Zampini
422d9255fafSStefano Zampini *remesh = PETSC_TRUE;
423d9255fafSStefano Zampini PetscCall(TSGetDM(ts, &dm));
424d9255fafSStefano Zampini PetscCall(Remesh(dm, U, &ctx->remesh_dm));
425d9255fafSStefano Zampini }
426d9255fafSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
427d9255fafSStefano Zampini }
428d9255fafSStefano Zampini
TransferVecs(TS ts,PetscInt nv,Vec vin[],Vec vout[],void * tctx)429d9255fafSStefano Zampini static PetscErrorCode TransferVecs(TS ts, PetscInt nv, Vec vin[], Vec vout[], void *tctx)
430d9255fafSStefano Zampini {
431d9255fafSStefano Zampini AppCtx *ctx = (AppCtx *)tctx;
432d9255fafSStefano Zampini DM dm;
433d9255fafSStefano Zampini Mat Interp;
434d9255fafSStefano Zampini
435d9255fafSStefano Zampini PetscFunctionBeginUser;
436d9255fafSStefano Zampini PetscCall(TSGetDM(ts, &dm));
437d9255fafSStefano Zampini PetscCall(DMCreateInterpolation(dm, ctx->remesh_dm, &Interp, NULL));
438d9255fafSStefano Zampini for (PetscInt i = 0; i < nv; i++) {
439d9255fafSStefano Zampini PetscCall(DMCreateGlobalVector(ctx->remesh_dm, &vout[i]));
440d9255fafSStefano Zampini PetscCall(MatMult(Interp, vin[i], vout[i]));
441d9255fafSStefano Zampini }
442d9255fafSStefano Zampini PetscCall(MatDestroy(&Interp));
443d9255fafSStefano Zampini PetscCall(TSSetDM(ts, ctx->remesh_dm));
444d9255fafSStefano Zampini PetscCall(DMDestroy(&ctx->remesh_dm));
4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
446a3d0cf85SMatthew G. Knepley }
447a3d0cf85SMatthew G. Knepley
main(int argc,char ** argv)448d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
449d71ae5a4SJacob Faibussowitsch {
450c4762a1bSJed Brown DM dm;
451c4762a1bSJed Brown TS ts;
452a3d0cf85SMatthew G. Knepley Vec u;
453a3d0cf85SMatthew G. Knepley AppCtx ctx;
454c4762a1bSJed Brown
455327415f7SBarry Smith PetscFunctionBeginUser;
4569566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
4579566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
4589566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
4599566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx));
4609566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx));
461c4762a1bSJed Brown
4629566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
4639566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm));
4649566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
465742ee2edSMatthew G. Knepley if (ctx.expTS) {
4669566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunctionLocal(dm, DMPlexTSComputeRHSFunctionFEM, &ctx));
4679566063dSJacob Faibussowitsch if (ctx.lumped) PetscCall(DMTSCreateRHSMassMatrixLumped(dm));
4689566063dSJacob Faibussowitsch else PetscCall(DMTSCreateRHSMassMatrix(dm));
469742ee2edSMatthew G. Knepley } else {
4709566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
4719566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
472742ee2edSMatthew G. Knepley }
473d9255fafSStefano Zampini PetscCall(TSSetMaxTime(ts, 1.0));
4749566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
4759566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
4769566063dSJacob Faibussowitsch PetscCall(TSSetComputeInitialCondition(ts, SetInitialConditions));
477ecc87898SStefano Zampini PetscCall(TSSetResize(ts, PETSC_FALSE, TransferSetUp, TransferVecs, &ctx));
478c4762a1bSJed Brown
4799566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
4809566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u));
4819566063dSJacob Faibussowitsch PetscCall(SetInitialConditions(ts, u));
4829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "temperature"));
483d9255fafSStefano Zampini
484d9255fafSStefano Zampini PetscCall(TSSetSolution(ts, u));
485d9255fafSStefano Zampini PetscCall(VecViewFromOptions(u, NULL, "-u0_view"));
486d9255fafSStefano Zampini PetscCall(VecDestroy(&u));
487d9255fafSStefano Zampini PetscCall(TSSolve(ts, NULL));
488d9255fafSStefano Zampini
489d9255fafSStefano Zampini PetscCall(TSGetSolution(ts, &u));
490d9255fafSStefano Zampini PetscCall(VecViewFromOptions(u, NULL, "-uf_view"));
4919566063dSJacob Faibussowitsch PetscCall(DMTSCheckFromOptions(ts, u));
4929566063dSJacob Faibussowitsch if (ctx.expTS) PetscCall(DMTSDestroyRHSMassMatrix(dm));
493c4762a1bSJed Brown
4949566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
4959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
4969566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
497b122ec5aSJacob Faibussowitsch return 0;
498c4762a1bSJed Brown }
499c4762a1bSJed Brown
500c4762a1bSJed Brown /*TEST
501c4762a1bSJed Brown
502c4762a1bSJed Brown test:
503a3d0cf85SMatthew G. Knepley suffix: 2d_p1
504c4762a1bSJed Brown requires: triangle
505a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
506188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
507c4762a1bSJed Brown test:
508742ee2edSMatthew G. Knepley suffix: 2d_p1_exp
509742ee2edSMatthew G. Knepley requires: triangle
510742ee2edSMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -explicit \
511188af4bfSBarry Smith -ts_type euler -ts_max_steps 4 -ts_time_step 1e-3 -ts_monitor_error
512742ee2edSMatthew G. Knepley test:
513a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
514a3d0cf85SMatthew G. Knepley suffix: 2d_p1_sconv
515c4762a1bSJed Brown requires: triangle
516a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
517188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00001 -snes_error_if_not_converged -pc_type lu
518c4762a1bSJed Brown test:
519742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.1]
520742ee2edSMatthew G. Knepley suffix: 2d_p1_sconv_2
521742ee2edSMatthew G. Knepley requires: triangle
522742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
523188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 1e-6 -snes_error_if_not_converged -pc_type lu
524742ee2edSMatthew G. Knepley test:
525a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
526a3d0cf85SMatthew G. Knepley suffix: 2d_p1_tconv
527c4762a1bSJed Brown requires: triangle
528a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
529188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
530c4762a1bSJed Brown test:
531742ee2edSMatthew G. Knepley # -dm_refine 6 -convest_num_refine 3 get L_2 convergence rate: [1.0]
532742ee2edSMatthew G. Knepley suffix: 2d_p1_tconv_2
533742ee2edSMatthew G. Knepley requires: triangle
534742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
535188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
536742ee2edSMatthew G. Knepley test:
537742ee2edSMatthew 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
538742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2
539742ee2edSMatthew G. Knepley requires: triangle
540742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
541188af4bfSBarry Smith -ts_type euler -ts_max_steps 4 -ts_time_step 1e-4 -lumped 0 -mass_pc_type lu
542742ee2edSMatthew G. Knepley test:
543742ee2edSMatthew 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
544742ee2edSMatthew G. Knepley suffix: 2d_p1_exp_tconv_2_lump
545742ee2edSMatthew G. Knepley requires: triangle
546742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 1 -explicit -ts_convergence_estimate -convest_num_refine 1 \
547188af4bfSBarry Smith -ts_type euler -ts_max_steps 4 -ts_time_step 1e-4
548742ee2edSMatthew G. Knepley test:
549a3d0cf85SMatthew G. Knepley suffix: 2d_p2
550c4762a1bSJed Brown requires: triangle
551a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
552188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
553c4762a1bSJed Brown test:
554a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
555a3d0cf85SMatthew G. Knepley suffix: 2d_p2_sconv
556a3d0cf85SMatthew G. Knepley requires: triangle
557a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
558188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00000001 -snes_error_if_not_converged -pc_type lu
559c4762a1bSJed Brown test:
560742ee2edSMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [3.1]
561742ee2edSMatthew G. Knepley suffix: 2d_p2_sconv_2
562742ee2edSMatthew G. Knepley requires: triangle
563742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
564188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00000001 -snes_error_if_not_converged -pc_type lu
565742ee2edSMatthew G. Knepley test:
566a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
567a3d0cf85SMatthew G. Knepley suffix: 2d_p2_tconv
568a3d0cf85SMatthew G. Knepley requires: triangle
569a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
570188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
571c4762a1bSJed Brown test:
572742ee2edSMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
573742ee2edSMatthew G. Knepley suffix: 2d_p2_tconv_2
574742ee2edSMatthew G. Knepley requires: triangle
575742ee2edSMatthew G. Knepley args: -sol_type trig_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
576188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
577742ee2edSMatthew G. Knepley test:
578a3d0cf85SMatthew G. Knepley suffix: 2d_q1
57930602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
580188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
581c4762a1bSJed Brown test:
582a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
583a3d0cf85SMatthew G. Knepley suffix: 2d_q1_sconv
58430602db0SMatthew 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 \
585188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00001 -snes_error_if_not_converged -pc_type lu
586c4762a1bSJed Brown test:
587a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
588a3d0cf85SMatthew G. Knepley suffix: 2d_q1_tconv
58930602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
590188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
591a3d0cf85SMatthew G. Knepley test:
592a3d0cf85SMatthew G. Knepley suffix: 2d_q2
59330602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
594188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
595a3d0cf85SMatthew G. Knepley test:
596a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
597a3d0cf85SMatthew G. Knepley suffix: 2d_q2_sconv
59830602db0SMatthew 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 \
599188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00000001 -snes_error_if_not_converged -pc_type lu
600a3d0cf85SMatthew G. Knepley test:
601a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
602a3d0cf85SMatthew G. Knepley suffix: 2d_q2_tconv
60330602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
604188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
605a3d0cf85SMatthew G. Knepley
606a3d0cf85SMatthew G. Knepley test:
607a3d0cf85SMatthew G. Knepley suffix: 3d_p1
608c4762a1bSJed Brown requires: ctetgen
609a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
610188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
611c4762a1bSJed Brown test:
612a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
613a3d0cf85SMatthew G. Knepley suffix: 3d_p1_sconv
614c4762a1bSJed Brown requires: ctetgen
615a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -temp_petscspace_degree 1 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
616188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00001 -snes_error_if_not_converged -pc_type lu
617c4762a1bSJed Brown test:
618a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
619a3d0cf85SMatthew G. Knepley suffix: 3d_p1_tconv
620c4762a1bSJed Brown requires: ctetgen
621a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
622188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
623c4762a1bSJed Brown test:
624a3d0cf85SMatthew G. Knepley suffix: 3d_p2
625c4762a1bSJed Brown requires: ctetgen
626a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_linear -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
627188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
628c4762a1bSJed Brown test:
629a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
630a3d0cf85SMatthew G. Knepley suffix: 3d_p2_sconv
631a3d0cf85SMatthew G. Knepley requires: ctetgen
632a3d0cf85SMatthew G. Knepley args: -sol_type trig_linear -temp_petscspace_degree 2 -ts_convergence_estimate -ts_convergence_temporal 0 -convest_num_refine 1 \
633188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00000001 -snes_error_if_not_converged -pc_type lu
634c4762a1bSJed Brown test:
635a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
636a3d0cf85SMatthew G. Knepley suffix: 3d_p2_tconv
637a3d0cf85SMatthew G. Knepley requires: ctetgen
638a3d0cf85SMatthew G. Knepley args: -sol_type quadratic_trig -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
639188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
640c4762a1bSJed Brown test:
641a3d0cf85SMatthew G. Knepley suffix: 3d_q1
64230602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 1 -temp_petscspace_degree 1 -dmts_check .0001 \
643188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
644c4762a1bSJed Brown test:
645a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [1.9]
646a3d0cf85SMatthew G. Knepley suffix: 3d_q1_sconv
64730602db0SMatthew 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 \
648188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00001 -snes_error_if_not_converged -pc_type lu
649a3d0cf85SMatthew G. Knepley test:
650a3d0cf85SMatthew G. Knepley # -dm_refine 4 -convest_num_refine 3 get L_2 convergence rate: [1.2]
651a3d0cf85SMatthew G. Knepley suffix: 3d_q1_tconv
65230602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 1 -ts_convergence_estimate -convest_num_refine 1 \
653188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
654a3d0cf85SMatthew G. Knepley test:
655a3d0cf85SMatthew G. Knepley suffix: 3d_q2
65630602db0SMatthew G. Knepley args: -sol_type quadratic_linear -dm_plex_simplex 0 -dm_refine 0 -temp_petscspace_degree 2 -dmts_check .0001 \
657188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
658a3d0cf85SMatthew G. Knepley test:
659a3d0cf85SMatthew G. Knepley # -dm_refine 2 -convest_num_refine 3 get L_2 convergence rate: [2.9]
660a3d0cf85SMatthew G. Knepley suffix: 3d_q2_sconv
66130602db0SMatthew 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 \
662188af4bfSBarry Smith -ts_type beuler -ts_max_steps 1 -ts_time_step 0.00000001 -snes_error_if_not_converged -pc_type lu
663a3d0cf85SMatthew G. Knepley test:
664a3d0cf85SMatthew G. Knepley # -dm_refine 3 -convest_num_refine 3 get L_2 convergence rate: [1.0]
665a3d0cf85SMatthew G. Knepley suffix: 3d_q2_tconv
66630602db0SMatthew G. Knepley args: -sol_type quadratic_trig -dm_plex_simplex 0 -temp_petscspace_degree 2 -ts_convergence_estimate -convest_num_refine 1 \
667188af4bfSBarry Smith -ts_type beuler -ts_max_steps 4 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
668a3d0cf85SMatthew G. Knepley
669a3d0cf85SMatthew G. Knepley test:
670a3d0cf85SMatthew 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
671a3d0cf85SMatthew G. Knepley suffix: egads_sphere
6725552b385SBrandon requires: egads ctetgen datafilespath
67330602db0SMatthew G. Knepley args: -sol_type quadratic_linear \
6745552b385SBrandon -dm_plex_boundary_filename ${DATAFILESPATH}/meshes/cad/sphere_example.egadslite -dm_plex_boundary_label marker \
675a3d0cf85SMatthew G. Knepley -temp_petscspace_degree 2 -dmts_check .0001 \
676188af4bfSBarry Smith -ts_type beuler -ts_max_steps 5 -ts_time_step 0.1 -snes_error_if_not_converged -pc_type lu
677c4762a1bSJed Brown
678d9255fafSStefano Zampini test:
679d9255fafSStefano Zampini suffix: remesh
680d9255fafSStefano Zampini requires: triangle mmg
681188af4bfSBarry Smith args: -sol_type quadratic_trig -dm_refine 2 -temp_petscspace_degree 1 -ts_type beuler -ts_time_step 0.01 -snes_error_if_not_converged -pc_type lu -grad_petscspace_degree 1 -dm_adaptor mmg -dm_plex_hash_location -remesh_every 5
682e0008caeSPierre Jolivet output_file: output/empty.out
683d9255fafSStefano Zampini
684c4762a1bSJed Brown TEST*/
685