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