xref: /petsc/src/ts/tutorials/ex47.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Pure advection with finite elements.\n\
2c4762a1bSJed Brown We solve the hyperbolic problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown The continuity equation (https://en.wikipedia.org/wiki/Continuity_equation) for advection
7c4762a1bSJed Brown (https://en.wikipedia.org/wiki/Advection) of a conserved scalar quantity phi, with source q,
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   phi_t + div (phi u) = q
10c4762a1bSJed Brown 
11c4762a1bSJed Brown if used with a solenoidal velocity field u (div u = 0) is given by
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   phi_t + u . grad phi = q
14c4762a1bSJed Brown 
15c4762a1bSJed Brown For a vector quantity a, we likewise have
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   a_t + u . grad a = q
18c4762a1bSJed Brown */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /*
21c4762a1bSJed Brown   r1: 8 SOR
22c4762a1bSJed Brown   r2: 1128 SOR
23c4762a1bSJed Brown   r3: > 10000 SOR
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   SOR is completely unreliable as a smoother, use Jacobi
26c4762a1bSJed Brown   r1: 8 MG
27c4762a1bSJed Brown   r2:
28c4762a1bSJed Brown */
29c4762a1bSJed Brown 
30c4762a1bSJed Brown #include <petscdmplex.h>
31c4762a1bSJed Brown #include <petscts.h>
32c4762a1bSJed Brown #include <petscds.h>
33c4762a1bSJed Brown 
34*9371c9d4SSatish Balay typedef enum {
35*9371c9d4SSatish Balay   PRIMITIVE,
36*9371c9d4SSatish Balay   INT_BY_PARTS
37*9371c9d4SSatish Balay } WeakFormType;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown typedef struct {
40c4762a1bSJed Brown   WeakFormType formType;
41c4762a1bSJed Brown } AppCtx;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /* MMS1:
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   2D:
46c4762a1bSJed Brown   u   = <1, 1>
47c4762a1bSJed Brown   phi = x + y - 2t
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   phi_t + u . grad phi = -2 + <1, 1> . <1, 1> = 0
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   3D:
52c4762a1bSJed Brown   u   = <1, 1, 1>
53c4762a1bSJed Brown   phi = x + y + z - 3t
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   phi_t + u . grad phi = -3 + <1, 1, 1> . <1, 1, 1> = 0
56c4762a1bSJed Brown */
57c4762a1bSJed Brown 
58*9371c9d4SSatish Balay static PetscErrorCode analytic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) {
59c4762a1bSJed Brown   PetscInt d;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   *u = -dim * time;
62c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d];
63c4762a1bSJed Brown   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66*9371c9d4SSatish Balay static PetscErrorCode velocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) {
67c4762a1bSJed Brown   PetscInt d;
68c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = 1.0;
69c4762a1bSJed Brown   return 0;
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /* <psi, phi_t> + <psi, u . grad phi> */
73*9371c9d4SSatish Balay static void f0_prim_phi(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[]) {
74c4762a1bSJed Brown   PetscInt d;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   f0[0] = u_t[0];
77c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += a[d] * u_x[d];
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown /* <psi, phi_t> */
81*9371c9d4SSatish Balay static void f0_ibp_phi(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[]) {
82c4762a1bSJed Brown   f0[0] = u_t[0];
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /* <grad psi, u phi> */
86*9371c9d4SSatish Balay static void f1_ibp_phi(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[]) {
87c4762a1bSJed Brown   PetscInt d;
88c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = a[d] * u[0];
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /* <psi, phi_t> */
92*9371c9d4SSatish Balay static void g0_prim_phi(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[]) {
93c4762a1bSJed Brown   g0[0] = u_tShift * 1.0;
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /* <psi, u . grad phi> */
97*9371c9d4SSatish Balay static void g1_prim_phi(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 g1[]) {
98c4762a1bSJed Brown   PetscInt d;
99c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d] = a[d];
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /* <grad psi, u phi> */
103*9371c9d4SSatish Balay static void g2_ibp_phi(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 g2[]) {
104c4762a1bSJed Brown   PetscInt d;
105c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d] = a[d];
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
109c4762a1bSJed Brown   const char *formTypes[2] = {"primitive", "int_by_parts"};
110c4762a1bSJed Brown   PetscInt    form;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   PetscFunctionBeginUser;
113c4762a1bSJed Brown   options->formType = PRIMITIVE;
114d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX");
115c4762a1bSJed Brown   form = options->formType;
1169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL));
117c4762a1bSJed Brown   options->formType = (WeakFormType)form;
118d0609cedSBarry Smith   PetscOptionsEnd();
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx) {
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1249566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1259566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1269566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1279566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
128c4762a1bSJed Brown   PetscFunctionReturn(0);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131*9371c9d4SSatish Balay static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx) {
13245480ffeSMatthew G. Knepley   PetscDS        ds;
13345480ffeSMatthew G. Knepley   DMLabel        label;
134c4762a1bSJed Brown   const PetscInt id = 1;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBeginUser;
1379566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
138c4762a1bSJed Brown   switch (ctx->formType) {
139c4762a1bSJed Brown   case PRIMITIVE:
1409566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_prim_phi, NULL));
1419566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL));
142c4762a1bSJed Brown     break;
143c4762a1bSJed Brown   case INT_BY_PARTS:
1449566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi));
1459566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL));
146c4762a1bSJed Brown     break;
147c4762a1bSJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, analytic_phi, ctx));
1499566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1509566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))analytic_phi, NULL, ctx, NULL));
151c4762a1bSJed Brown   PetscFunctionReturn(0);
152c4762a1bSJed Brown }
153c4762a1bSJed Brown 
154*9371c9d4SSatish Balay static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user) {
15530602db0SMatthew G. Knepley   PetscSimplePointFunc funcs[1] = {velocity};
156c4762a1bSJed Brown   Vec                  v;
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   PetscFunctionBeginUser;
1599566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmAux, &v));
1609566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v));
1619566063dSJacob Faibussowitsch   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, v));
1629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166*9371c9d4SSatish Balay static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user) {
167c4762a1bSJed Brown   DM dmAux, coordDM;
168c4762a1bSJed Brown 
1697510d9b0SBarry Smith   PetscFunctionBeginUser;
170c4762a1bSJed Brown   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
1719566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &coordDM));
172c4762a1bSJed Brown   if (!feAux) PetscFunctionReturn(0);
1739566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmAux));
1749566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
1759566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject)feAux));
1769566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmAux));
1779566063dSJacob Faibussowitsch   PetscCall(SetupVelocity(dm, dmAux, user));
1789566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAux));
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182*9371c9d4SSatish Balay static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx) {
183c4762a1bSJed Brown   DM        cdm = dm;
184c4762a1bSJed Brown   PetscFE   fe, feAux;
185c4762a1bSJed Brown   MPI_Comm  comm;
18630602db0SMatthew G. Knepley   PetscInt  dim;
18730602db0SMatthew G. Knepley   PetscBool simplex;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   PetscFunctionBeginUser;
1909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1919566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
1929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1939566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "phi_", -1, &fe));
1949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)fe, "phi"));
1959566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", -1, &feAux));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe, feAux));
1979566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1989566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
1999566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, ctx));
200c4762a1bSJed Brown   while (cdm) {
2019566063dSJacob Faibussowitsch     PetscCall(SetupAuxDM(cdm, feAux, ctx));
2029566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
2039566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
204c4762a1bSJed Brown   }
2059566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
2069566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&feAux));
207c4762a1bSJed Brown   PetscFunctionReturn(0);
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210*9371c9d4SSatish Balay static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx) {
211c4762a1bSJed Brown   DM                   dm;
21230602db0SMatthew G. Knepley   PetscDS              ds;
21330602db0SMatthew G. Knepley   PetscSimplePointFunc func[1];
21430602db0SMatthew G. Knepley   void                *ctxs[1];
215c4762a1bSJed Brown   Vec                  u, r, error;
216c4762a1bSJed Brown   PetscReal            time = 0.5, res;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(ksp, &dm));
2209566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, it, time));
221c4762a1bSJed Brown   /* Calculate residual */
2229566063dSJacob Faibussowitsch   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
2239566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
2249566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, it, res));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)r, "residual"));
2269566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(r, NULL, "-res_vec_view"));
2279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
228c4762a1bSJed Brown   /* Calculate error */
2299566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2309566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2319566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &u));
2329566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &error));
2339566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, time, func, ctxs, INSERT_ALL_VALUES, error));
2349566063dSJacob Faibussowitsch   PetscCall(VecAXPY(error, -1.0, u));
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)error, "error"));
2369566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(error, NULL, "-err_vec_view"));
2379566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &error));
238c4762a1bSJed Brown   PetscFunctionReturn(0);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
241*9371c9d4SSatish Balay static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx) {
242c4762a1bSJed Brown   DM                   dm;
24330602db0SMatthew G. Knepley   PetscDS              ds;
24430602db0SMatthew G. Knepley   PetscSimplePointFunc func[1];
24530602db0SMatthew G. Knepley   void                *ctxs[1];
246c4762a1bSJed Brown   PetscReal            error;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   PetscFunctionBeginUser;
2499566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2509566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2519566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2529566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, crtime, func, ctxs, u, &error));
2539566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int)step, (double)crtime, (double)error));
254c4762a1bSJed Brown   PetscFunctionReturn(0);
255c4762a1bSJed Brown }
256c4762a1bSJed Brown 
257*9371c9d4SSatish Balay int main(int argc, char **argv) {
258c4762a1bSJed Brown   AppCtx    ctx;
259c4762a1bSJed Brown   DM        dm;
260c4762a1bSJed Brown   TS        ts;
261c4762a1bSJed Brown   Vec       u, r;
262c4762a1bSJed Brown   PetscReal t = 0.0;
263c4762a1bSJed Brown 
264327415f7SBarry Smith   PetscFunctionBeginUser;
2659566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2669566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
2679566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
2689566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &ctx));
2699566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &ctx));
270c4762a1bSJed Brown 
2719566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
2729566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "phi"));
2739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
274c4762a1bSJed Brown 
2759566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2769566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL));
2779566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
2789566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
2799566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
2809566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
2819566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2829566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
283c4762a1bSJed Brown 
28430602db0SMatthew G. Knepley   {
28530602db0SMatthew G. Knepley     PetscDS              ds;
28630602db0SMatthew G. Knepley     PetscSimplePointFunc func[1];
28730602db0SMatthew G. Knepley     void                *ctxs[1];
28830602db0SMatthew G. Knepley 
2899566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
2909566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2919566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, t, func, ctxs, INSERT_ALL_VALUES, u));
29230602db0SMatthew G. Knepley   }
293c4762a1bSJed Brown   {
294c4762a1bSJed Brown     SNES snes;
295c4762a1bSJed Brown     KSP  ksp;
296c4762a1bSJed Brown 
2979566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
2989566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
2999566063dSJacob Faibussowitsch     PetscCall(KSPMonitorSet(ksp, MonitorError, &ctx, NULL));
300c4762a1bSJed Brown   }
3019566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
3029566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
303c4762a1bSJed Brown 
3049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3069566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3079566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
309b122ec5aSJacob Faibussowitsch   return 0;
310c4762a1bSJed Brown }
311c4762a1bSJed Brown 
312c4762a1bSJed Brown /*TEST
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   # Full solves
315c4762a1bSJed Brown   test:
316c4762a1bSJed Brown     suffix: 2d_p1p1_r1
317c4762a1bSJed Brown     requires: triangle
318c4762a1bSJed Brown     args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -pc_type lu -snes_monitor_short -snes_converged_reason -ts_monitor
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   test:
321c4762a1bSJed Brown     suffix: 2d_p1p1_sor_r1
322c4762a1bSJed Brown     requires: triangle !single
323c4762a1bSJed Brown     args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ksp_rtol 1.0e-9 -pc_type sor -snes_monitor_short -snes_converged_reason -ksp_monitor_short -ts_monitor
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   test:
326c4762a1bSJed Brown     suffix: 2d_p1p1_mg_r1
327c4762a1bSJed Brown     requires: triangle !single
328c4762a1bSJed Brown     args: -dm_refine_hierarchy 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_dt 0.1 -ksp_type fgmres -ksp_rtol 1.0e-9 -pc_type mg -pc_mg_levels 2 -snes_monitor_short -snes_converged_reason -snes_view -ksp_monitor_true_residual -ts_monitor
329c4762a1bSJed Brown 
330c4762a1bSJed Brown TEST*/
331