xref: /petsc/src/ts/tutorials/ex47.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 
34c4762a1bSJed Brown typedef enum {PRIMITIVE, INT_BY_PARTS} WeakFormType;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown typedef struct {
37c4762a1bSJed Brown   WeakFormType formType;
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown /* MMS1:
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   2D:
43c4762a1bSJed Brown   u   = <1, 1>
44c4762a1bSJed Brown   phi = x + y - 2t
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   phi_t + u . grad phi = -2 + <1, 1> . <1, 1> = 0
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   3D:
49c4762a1bSJed Brown   u   = <1, 1, 1>
50c4762a1bSJed Brown   phi = x + y + z - 3t
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   phi_t + u . grad phi = -3 + <1, 1, 1> . <1, 1, 1> = 0
53c4762a1bSJed Brown */
54c4762a1bSJed Brown 
55c4762a1bSJed Brown static PetscErrorCode analytic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
56c4762a1bSJed Brown {
57c4762a1bSJed Brown   PetscInt d;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   *u = -dim*time;
60c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d];
61c4762a1bSJed Brown   return 0;
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown static PetscErrorCode velocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   PetscInt d;
67c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = 1.0;
68c4762a1bSJed Brown   return 0;
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /* <psi, phi_t> + <psi, u . grad phi> */
72c4762a1bSJed Brown static void f0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
73c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
74c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
75c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
76c4762a1bSJed Brown {
77c4762a1bSJed Brown   PetscInt d;
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   f0[0] = u_t[0];
80c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += a[d] * u_x[d];
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown /* <psi, phi_t> */
84c4762a1bSJed Brown static void f0_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
85c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
86c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
87c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   f0[0] = u_t[0];
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92c4762a1bSJed Brown /* <grad psi, u phi> */
93c4762a1bSJed Brown static void f1_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
94c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
95c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
96c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   PetscInt d;
99c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = a[d]*u[0];
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /* <psi, phi_t> */
103c4762a1bSJed Brown static void g0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
104c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
105c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
106c4762a1bSJed Brown                         PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
107c4762a1bSJed Brown {
108c4762a1bSJed Brown   g0[0] = u_tShift*1.0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111c4762a1bSJed Brown /* <psi, u . grad phi> */
112c4762a1bSJed Brown static void g1_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
113c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
114c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
115c4762a1bSJed Brown                         PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
116c4762a1bSJed Brown {
117c4762a1bSJed Brown   PetscInt d;
118c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d] = a[d];
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown /* <grad psi, u phi> */
122c4762a1bSJed Brown static void g2_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
123c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
124c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
125c4762a1bSJed Brown                        PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   PetscInt d;
128c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d] = a[d];
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
132c4762a1bSJed Brown {
133c4762a1bSJed Brown   const char    *formTypes[2] = {"primitive", "int_by_parts"};
134c4762a1bSJed Brown   PetscInt       form;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBeginUser;
137c4762a1bSJed Brown   options->formType = PRIMITIVE;
138*d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX");
139c4762a1bSJed Brown   form = options->formType;
1409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL));
141c4762a1bSJed Brown   options->formType = (WeakFormType) form;
142*d0609cedSBarry Smith   PetscOptionsEnd();
143c4762a1bSJed Brown   PetscFunctionReturn(0);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
147c4762a1bSJed Brown {
148c4762a1bSJed Brown   PetscFunctionBeginUser;
1499566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1509566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1519566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1529566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
153c4762a1bSJed Brown   PetscFunctionReturn(0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
157c4762a1bSJed Brown {
15845480ffeSMatthew G. Knepley   PetscDS        ds;
15945480ffeSMatthew G. Knepley   DMLabel        label;
160c4762a1bSJed Brown   const PetscInt id = 1;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   PetscFunctionBeginUser;
1639566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
164c4762a1bSJed Brown   switch (ctx->formType) {
165c4762a1bSJed Brown   case PRIMITIVE:
1669566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_prim_phi, NULL));
1679566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL));
168c4762a1bSJed Brown     break;
169c4762a1bSJed Brown   case INT_BY_PARTS:
1709566063dSJacob Faibussowitsch     PetscCall(PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi));
1719566063dSJacob Faibussowitsch     PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL));
172c4762a1bSJed Brown     break;
173c4762a1bSJed Brown   }
1749566063dSJacob Faibussowitsch   PetscCall(PetscDSSetExactSolution(ds, 0, analytic_phi, ctx));
1759566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "marker", &label));
1769566063dSJacob Faibussowitsch   PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) analytic_phi, NULL, ctx, NULL));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user)
181c4762a1bSJed Brown {
18230602db0SMatthew G. Knepley   PetscSimplePointFunc funcs[1] = {velocity};
183c4762a1bSJed Brown   Vec                  v;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   PetscFunctionBeginUser;
1869566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmAux, &v));
1879566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v));
1889566063dSJacob Faibussowitsch   PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, v));
1899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
190c4762a1bSJed Brown   PetscFunctionReturn(0);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
193c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user)
194c4762a1bSJed Brown {
195c4762a1bSJed Brown   DM             dmAux, coordDM;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBegin;
198c4762a1bSJed Brown   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
1999566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &coordDM));
200c4762a1bSJed Brown   if (!feAux) PetscFunctionReturn(0);
2019566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmAux));
2029566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDM(dmAux, coordDM));
2039566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject) feAux));
2049566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmAux));
2059566063dSJacob Faibussowitsch   PetscCall(SetupVelocity(dm, dmAux, user));
2069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAux));
207c4762a1bSJed Brown   PetscFunctionReturn(0);
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
211c4762a1bSJed Brown {
212c4762a1bSJed Brown   DM             cdm = dm;
213c4762a1bSJed Brown   PetscFE        fe,   feAux;
214c4762a1bSJed Brown   MPI_Comm       comm;
21530602db0SMatthew G. Knepley   PetscInt       dim;
21630602db0SMatthew G. Knepley   PetscBool      simplex;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
2209566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
2229566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "phi_", -1, &fe));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) fe, "phi"));
2249566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", -1, &feAux));
2259566063dSJacob Faibussowitsch   PetscCall(PetscFECopyQuadrature(fe, feAux));
2269566063dSJacob Faibussowitsch   PetscCall(DMSetField(dm, 0, NULL, (PetscObject) fe));
2279566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dm));
2289566063dSJacob Faibussowitsch   PetscCall(SetupProblem(dm, ctx));
229c4762a1bSJed Brown   while (cdm) {
2309566063dSJacob Faibussowitsch     PetscCall(SetupAuxDM(cdm, feAux, ctx));
2319566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, cdm));
2329566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(cdm, &cdm));
233c4762a1bSJed Brown   }
2349566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&feAux));
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239798534f6SMatthew G. Knepley static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
240c4762a1bSJed Brown {
241c4762a1bSJed Brown   DM                   dm;
24230602db0SMatthew G. Knepley   PetscDS              ds;
24330602db0SMatthew G. Knepley   PetscSimplePointFunc func[1];
24430602db0SMatthew G. Knepley   void                *ctxs[1];
245c4762a1bSJed Brown   Vec                  u, r, error;
246c4762a1bSJed Brown   PetscReal            time = 0.5, res;
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   PetscFunctionBeginUser;
2499566063dSJacob Faibussowitsch   PetscCall(KSPGetDM(ksp, &dm));
2509566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, it, time));
251c4762a1bSJed Brown   /* Calculate residual */
2529566063dSJacob Faibussowitsch   PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
2539566063dSJacob Faibussowitsch   PetscCall(VecNorm(r, NORM_2, &res));
2549566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, it, res));
2559566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) r, "residual"));
2569566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(r, NULL, "-res_vec_view"));
2579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
258c4762a1bSJed Brown   /* Calculate error */
2599566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2609566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2619566063dSJacob Faibussowitsch   PetscCall(KSPBuildSolution(ksp, NULL, &u));
2629566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &error));
2639566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dm, time, func, ctxs, INSERT_ALL_VALUES, error));
2649566063dSJacob Faibussowitsch   PetscCall(VecAXPY(error, -1.0, u));
2659566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) error, "error"));
2669566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(error, NULL, "-err_vec_view"));
2679566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &error));
268c4762a1bSJed Brown   PetscFunctionReturn(0);
269c4762a1bSJed Brown }
270c4762a1bSJed Brown 
271c4762a1bSJed Brown static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
272c4762a1bSJed Brown {
273c4762a1bSJed Brown   DM                   dm;
27430602db0SMatthew G. Knepley   PetscDS              ds;
27530602db0SMatthew G. Knepley   PetscSimplePointFunc func[1];
27630602db0SMatthew G. Knepley   void                *ctxs[1];
277c4762a1bSJed Brown   PetscReal            error;
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   PetscFunctionBeginUser;
2809566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &dm));
2819566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &ds));
2829566063dSJacob Faibussowitsch   PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2839566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dm, crtime, func, ctxs, u, &error));
2849566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int) step, (double) crtime, (double) error));
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown int main(int argc, char **argv)
289c4762a1bSJed Brown {
290c4762a1bSJed Brown   AppCtx         ctx;
291c4762a1bSJed Brown   DM             dm;
292c4762a1bSJed Brown   TS             ts;
293c4762a1bSJed Brown   Vec            u, r;
294c4762a1bSJed Brown   PetscReal      t       = 0.0;
295c4762a1bSJed Brown 
2969566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2979566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
2989566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
2999566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(dm, &ctx));
3009566063dSJacob Faibussowitsch   PetscCall(SetupDiscretization(dm, &ctx));
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &u));
3039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "phi"));
3049566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &r));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3079566063dSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL));
3089566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, dm));
3099566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
3109566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
3119566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
3129566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3139566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
314c4762a1bSJed Brown 
31530602db0SMatthew G. Knepley   {
31630602db0SMatthew G. Knepley     PetscDS              ds;
31730602db0SMatthew G. Knepley     PetscSimplePointFunc func[1];
31830602db0SMatthew G. Knepley     void                *ctxs[1];
31930602db0SMatthew G. Knepley 
3209566063dSJacob Faibussowitsch     PetscCall(DMGetDS(dm, &ds));
3219566063dSJacob Faibussowitsch     PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
3229566063dSJacob Faibussowitsch     PetscCall(DMProjectFunction(dm, t, func, ctxs, INSERT_ALL_VALUES, u));
32330602db0SMatthew G. Knepley   }
324c4762a1bSJed Brown   {
325c4762a1bSJed Brown     SNES snes;
326c4762a1bSJed Brown     KSP  ksp;
327c4762a1bSJed Brown 
3289566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
3299566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
3309566063dSJacob Faibussowitsch     PetscCall(KSPMonitorSet(ksp, MonitorError, &ctx, NULL));
331c4762a1bSJed Brown   }
3329566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, u));
3339566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
334c4762a1bSJed Brown 
3359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
3369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
3379566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3389566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3399566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
340b122ec5aSJacob Faibussowitsch   return 0;
341c4762a1bSJed Brown }
342c4762a1bSJed Brown 
343c4762a1bSJed Brown /*TEST
344c4762a1bSJed Brown 
345c4762a1bSJed Brown   # Full solves
346c4762a1bSJed Brown   test:
347c4762a1bSJed Brown     suffix: 2d_p1p1_r1
348c4762a1bSJed Brown     requires: triangle
349c4762a1bSJed 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
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   test:
352c4762a1bSJed Brown     suffix: 2d_p1p1_sor_r1
353c4762a1bSJed Brown     requires: triangle !single
354c4762a1bSJed 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
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   test:
357c4762a1bSJed Brown     suffix: 2d_p1p1_mg_r1
358c4762a1bSJed Brown     requires: triangle !single
359c4762a1bSJed 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
360c4762a1bSJed Brown 
361c4762a1bSJed Brown TEST*/
362