xref: /petsc/src/ts/tutorials/ex47.c (revision 45480ffeba491aca5d3f3b8c78679069fb317d32)
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   /* Domain and mesh definition */
38c4762a1bSJed Brown   PetscInt          dim;               /* The topological mesh dimension */
39c4762a1bSJed Brown   PetscBool         simplex;           /* Use simplices or tensor product cells */
40c4762a1bSJed Brown   /* Problem definition */
41c4762a1bSJed Brown   WeakFormType      formType;
42c4762a1bSJed Brown   PetscErrorCode (**exactFuncs)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx);
43c4762a1bSJed Brown } AppCtx;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /* MMS1:
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   2D:
48c4762a1bSJed Brown   u   = <1, 1>
49c4762a1bSJed Brown   phi = x + y - 2t
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   phi_t + u . grad phi = -2 + <1, 1> . <1, 1> = 0
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   3D:
54c4762a1bSJed Brown   u   = <1, 1, 1>
55c4762a1bSJed Brown   phi = x + y + z - 3t
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   phi_t + u . grad phi = -3 + <1, 1, 1> . <1, 1, 1> = 0
58c4762a1bSJed Brown */
59c4762a1bSJed Brown 
60c4762a1bSJed Brown static PetscErrorCode analytic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   PetscInt d;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   *u = -dim*time;
65c4762a1bSJed Brown   for (d = 0; d < dim; ++d) *u += x[d];
66c4762a1bSJed Brown   return 0;
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown static PetscErrorCode velocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx)
70c4762a1bSJed Brown {
71c4762a1bSJed Brown   PetscInt d;
72c4762a1bSJed Brown   for (d = 0; d < dim; ++d) u[d] = 1.0;
73c4762a1bSJed Brown   return 0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown /* <psi, phi_t> + <psi, u . grad phi> */
77c4762a1bSJed Brown static void f0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
78c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
79c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
80c4762a1bSJed Brown                         PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
81c4762a1bSJed Brown {
82c4762a1bSJed Brown   PetscInt d;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   f0[0] = u_t[0];
85c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f0[0] += a[d] * u_x[d];
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /* <psi, phi_t> */
89c4762a1bSJed Brown static void f0_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
90c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
91c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
92c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
93c4762a1bSJed Brown {
94c4762a1bSJed Brown   f0[0] = u_t[0];
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown /* <grad psi, u phi> */
98c4762a1bSJed Brown static void f1_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
99c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
100c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
101c4762a1bSJed Brown                        PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   PetscInt d;
104c4762a1bSJed Brown   for (d = 0; d < dim; ++d) f1[d] = a[d]*u[0];
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /* <psi, phi_t> */
108c4762a1bSJed Brown static void g0_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
109c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
110c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
111c4762a1bSJed Brown                         PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
112c4762a1bSJed Brown {
113c4762a1bSJed Brown   g0[0] = u_tShift*1.0;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown /* <psi, u . grad phi> */
117c4762a1bSJed Brown static void g1_prim_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
118c4762a1bSJed Brown                         const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
119c4762a1bSJed Brown                         const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
120c4762a1bSJed Brown                         PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g1[])
121c4762a1bSJed Brown {
122c4762a1bSJed Brown   PetscInt d;
123c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g1[d] = a[d];
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /* <grad psi, u phi> */
127c4762a1bSJed Brown static void g2_ibp_phi(PetscInt dim, PetscInt Nf, PetscInt NfAux,
128c4762a1bSJed Brown                        const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
129c4762a1bSJed Brown                        const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
130c4762a1bSJed Brown                        PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g2[])
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   PetscInt d;
133c4762a1bSJed Brown   for (d = 0; d < dim; ++d) g2[d] = a[d];
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   const char    *formTypes[2] = {"primitive", "int_by_parts"};
139c4762a1bSJed Brown   PetscInt       form;
140c4762a1bSJed Brown   PetscErrorCode ierr;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBeginUser;
143c4762a1bSJed Brown   options->dim      = 2;
144c4762a1bSJed Brown   options->simplex  = PETSC_TRUE;
145c4762a1bSJed Brown   options->formType = PRIMITIVE;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX");CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = PetscOptionsInt("-dim", "The topological mesh dimension", "ex47.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = PetscOptionsBool("-simplex", "Simplicial (true) or tensor (false) mesh", "ex47.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
150c4762a1bSJed Brown   form = options->formType;
151c4762a1bSJed Brown   ierr = PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL);CHKERRQ(ierr);
152c4762a1bSJed Brown   options->formType = (WeakFormType) form;
153c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
154c4762a1bSJed Brown   PetscFunctionReturn(0);
155c4762a1bSJed Brown }
156c4762a1bSJed Brown 
157c4762a1bSJed Brown static PetscErrorCode CreateBCLabel(DM dm, const char name[])
158c4762a1bSJed Brown {
159408cafa0SMatthew G. Knepley   DM             plex;
160c4762a1bSJed Brown   DMLabel        label;
161c4762a1bSJed Brown   PetscErrorCode ierr;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   PetscFunctionBeginUser;
164c4762a1bSJed Brown   ierr = DMCreateLabel(dm, name);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = DMGetLabel(dm, name, &label);CHKERRQ(ierr);
166408cafa0SMatthew G. Knepley   ierr = DMConvert(dm, DMPLEX, &plex);CHKERRQ(ierr);
167c4762a1bSJed Brown   ierr = DMPlexMarkBoundaryFaces(dm, 1, label);CHKERRQ(ierr);
168408cafa0SMatthew G. Knepley   ierr = DMDestroy(&plex);CHKERRQ(ierr);
169c4762a1bSJed Brown   PetscFunctionReturn(0);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
173c4762a1bSJed Brown {
174c4762a1bSJed Brown   DM             pdm = NULL;
175c4762a1bSJed Brown   const PetscInt dim = ctx->dim;
176c4762a1bSJed Brown   PetscBool      hasLabel;
177c4762a1bSJed Brown   PetscErrorCode ierr;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   PetscFunctionBeginUser;
180c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
181c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
182c4762a1bSJed Brown   /* If no boundary marker exists, mark the whole boundary */
183c4762a1bSJed Brown   ierr = DMHasLabel(*dm, "marker", &hasLabel);CHKERRQ(ierr);
184c4762a1bSJed Brown   if (!hasLabel) {ierr = CreateBCLabel(*dm, "marker");CHKERRQ(ierr);}
185c4762a1bSJed Brown   /* Distribute mesh over processes */
186c4762a1bSJed Brown   ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
187c4762a1bSJed Brown   if (pdm) {
188c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
189c4762a1bSJed Brown     *dm  = pdm;
190c4762a1bSJed Brown   }
191c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
197c4762a1bSJed Brown {
198*45480ffeSMatthew G. Knepley   PetscDS        ds;
199*45480ffeSMatthew G. Knepley   DMLabel        label;
200c4762a1bSJed Brown   const PetscInt id = 1;
201c4762a1bSJed Brown   PetscErrorCode ierr;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   PetscFunctionBeginUser;
204*45480ffeSMatthew G. Knepley   ierr = DMGetDS(dm, &ds);CHKERRQ(ierr);
205c4762a1bSJed Brown   switch (ctx->formType) {
206c4762a1bSJed Brown   case PRIMITIVE:
207*45480ffeSMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_prim_phi, NULL);CHKERRQ(ierr);
208*45480ffeSMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL);CHKERRQ(ierr);
209c4762a1bSJed Brown     break;
210c4762a1bSJed Brown   case INT_BY_PARTS:
211*45480ffeSMatthew G. Knepley     ierr = PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi);CHKERRQ(ierr);
212*45480ffeSMatthew G. Knepley     ierr = PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL);CHKERRQ(ierr);
213c4762a1bSJed Brown     break;
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown   ctx->exactFuncs[0] = analytic_phi;
216*45480ffeSMatthew G. Knepley   ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr);
217*45480ffeSMatthew G. Knepley   ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void)) ctx->exactFuncs[0], NULL, ctx, NULL);CHKERRQ(ierr);
218c4762a1bSJed Brown   PetscFunctionReturn(0);
219c4762a1bSJed Brown }
220c4762a1bSJed Brown 
221c4762a1bSJed Brown static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user)
222c4762a1bSJed Brown {
223c4762a1bSJed Brown   PetscErrorCode (*funcs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) = {velocity};
224c4762a1bSJed Brown   Vec            v;
225c4762a1bSJed Brown   PetscErrorCode ierr;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   PetscFunctionBeginUser;
228c4762a1bSJed Brown   ierr = DMCreateLocalVector(dmAux, &v);CHKERRQ(ierr);
229c4762a1bSJed Brown   ierr = DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v);CHKERRQ(ierr);
2309a2a23afSMatthew G. Knepley   ierr = DMSetAuxiliaryVec(dm, NULL, 0, v);CHKERRQ(ierr);
231c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
232c4762a1bSJed Brown   PetscFunctionReturn(0);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user)
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   DM             dmAux, coordDM;
238c4762a1bSJed Brown   PetscErrorCode ierr;
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   PetscFunctionBegin;
241c4762a1bSJed Brown   /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
242c4762a1bSJed Brown   ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
243c4762a1bSJed Brown   if (!feAux) PetscFunctionReturn(0);
244c4762a1bSJed Brown   ierr = DMClone(dm, &dmAux);CHKERRQ(ierr);
245c4762a1bSJed Brown   ierr = DMSetCoordinateDM(dmAux, coordDM);CHKERRQ(ierr);
246c4762a1bSJed Brown   ierr = DMSetField(dmAux, 0, NULL, (PetscObject) feAux);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = DMCreateDS(dmAux);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = SetupVelocity(dm, dmAux, user);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = DMDestroy(&dmAux);CHKERRQ(ierr);
250c4762a1bSJed Brown   PetscFunctionReturn(0);
251c4762a1bSJed Brown }
252c4762a1bSJed Brown 
253c4762a1bSJed Brown static PetscErrorCode SetupDiscretization(DM dm, AppCtx* ctx)
254c4762a1bSJed Brown {
255c4762a1bSJed Brown   DM              cdm = dm;
256c4762a1bSJed Brown   const PetscInt  dim = ctx->dim;
257c4762a1bSJed Brown   PetscFE         fe,   feAux;
258c4762a1bSJed Brown   MPI_Comm        comm;
259c4762a1bSJed Brown   PetscErrorCode  ierr;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   PetscFunctionBeginUser;
262c4762a1bSJed Brown   /* Create finite element */
263c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, 1, ctx->simplex, "phi_", -1, &fe);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) fe, "phi");CHKERRQ(ierr);
266c4762a1bSJed Brown   /* Create velocity */
267c4762a1bSJed Brown   ierr = PetscFECreateDefault(comm, dim, dim, ctx->simplex, "vel_", -1, &feAux);CHKERRQ(ierr);
268c4762a1bSJed Brown   ierr = PetscFECopyQuadrature(fe, feAux);CHKERRQ(ierr);
269c4762a1bSJed Brown   /* Set discretization and boundary conditions for each mesh */
270c4762a1bSJed Brown   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr = DMCreateDS(dm);CHKERRQ(ierr);
272c4762a1bSJed Brown   ierr = SetupProblem(dm, ctx);CHKERRQ(ierr);
273c4762a1bSJed Brown   while (cdm) {
274c4762a1bSJed Brown     PetscBool hasLabel;
275c4762a1bSJed Brown 
276c4762a1bSJed Brown     ierr = SetupAuxDM(cdm, feAux, ctx);CHKERRQ(ierr);
277c4762a1bSJed Brown     ierr = DMHasLabel(cdm, "marker", &hasLabel);CHKERRQ(ierr);
278c4762a1bSJed Brown     if (!hasLabel) {ierr = CreateBCLabel(cdm, "marker");CHKERRQ(ierr);}
279408cafa0SMatthew G. Knepley     ierr = DMCopyDisc(dm, cdm);CHKERRQ(ierr);
280c4762a1bSJed Brown     ierr = DMGetCoarseDM(cdm, &cdm);CHKERRQ(ierr);
281c4762a1bSJed Brown   }
282c4762a1bSJed Brown   ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
283c4762a1bSJed Brown   ierr = PetscFEDestroy(&feAux);CHKERRQ(ierr);
284c4762a1bSJed Brown   PetscFunctionReturn(0);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287798534f6SMatthew G. Knepley static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
288c4762a1bSJed Brown {
289c4762a1bSJed Brown   AppCtx        *user = (AppCtx *) ctx;
290c4762a1bSJed Brown   DM             dm;
291c4762a1bSJed Brown   Vec            u, r, error;
292c4762a1bSJed Brown   PetscReal      time = 0.5, res;
293c4762a1bSJed Brown   PetscErrorCode ierr;
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   PetscFunctionBeginUser;
296c4762a1bSJed Brown   ierr = KSPGetDM(ksp, &dm);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = DMSetOutputSequenceNumber(dm, it, time);CHKERRQ(ierr);
298c4762a1bSJed Brown   /* Calculate residual */
299c4762a1bSJed Brown   ierr = KSPBuildResidual(ksp, NULL, NULL, &r);CHKERRQ(ierr);
300c4762a1bSJed Brown   ierr = VecNorm(r, NORM_2, &res);CHKERRQ(ierr);
301c4762a1bSJed Brown   ierr = DMSetOutputSequenceNumber(dm, it, res);CHKERRQ(ierr);
302c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) r, "residual");CHKERRQ(ierr);
303c4762a1bSJed Brown   ierr = VecViewFromOptions(r, NULL, "-res_vec_view");CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
305c4762a1bSJed Brown   /* Calculate error */
306c4762a1bSJed Brown   ierr = KSPBuildSolution(ksp, NULL, &u);CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = DMGetGlobalVector(dm, &error);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = DMProjectFunction(dm, time, user->exactFuncs, NULL, INSERT_ALL_VALUES, error);CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = VecAXPY(error, -1.0, u);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) error, "error");CHKERRQ(ierr);
311c4762a1bSJed Brown   ierr = VecViewFromOptions(error, NULL, "-err_vec_view");CHKERRQ(ierr);
312c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(dm, &error);CHKERRQ(ierr);
313c4762a1bSJed Brown   PetscFunctionReturn(0);
314c4762a1bSJed Brown }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, void *ctx)
317c4762a1bSJed Brown {
318c4762a1bSJed Brown   AppCtx        *user = (AppCtx *) ctx;
319c4762a1bSJed Brown   DM             dm;
320c4762a1bSJed Brown   PetscReal      error;
321c4762a1bSJed Brown   PetscErrorCode ierr;
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   PetscFunctionBeginUser;
324c4762a1bSJed Brown   ierr = TSGetDM(ts, &dm);CHKERRQ(ierr);
325c4762a1bSJed Brown   ierr = DMComputeL2Diff(dm, crtime, user->exactFuncs, NULL, u, &error);CHKERRQ(ierr);
326c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int) step, (double) crtime, (double) error);CHKERRQ(ierr);
327c4762a1bSJed Brown   PetscFunctionReturn(0);
328c4762a1bSJed Brown }
329c4762a1bSJed Brown 
330c4762a1bSJed Brown int main(int argc, char **argv)
331c4762a1bSJed Brown {
332c4762a1bSJed Brown   AppCtx         ctx;
333c4762a1bSJed Brown   DM             dm;
334c4762a1bSJed Brown   TS             ts;
335c4762a1bSJed Brown   Vec            u, r;
336c4762a1bSJed Brown   PetscReal      t       = 0.0;
337c4762a1bSJed Brown   PetscErrorCode ierr;
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
340c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
341c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &dm, &ctx);CHKERRQ(ierr);
342c4762a1bSJed Brown   ierr = DMSetApplicationContext(dm, &ctx);CHKERRQ(ierr);
343c4762a1bSJed Brown   ierr = PetscMalloc1(1, &ctx.exactFuncs);CHKERRQ(ierr);
344c4762a1bSJed Brown   ierr = SetupDiscretization(dm, &ctx);CHKERRQ(ierr);
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm, &u);CHKERRQ(ierr);
347c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) u, "phi");CHKERRQ(ierr);
348c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
351c4762a1bSJed Brown   ierr = TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL);CHKERRQ(ierr);
352c4762a1bSJed Brown   ierr = TSSetDM(ts, dm);CHKERRQ(ierr);
353c4762a1bSJed Brown   ierr = DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx);CHKERRQ(ierr);
354c4762a1bSJed Brown   ierr = DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx);CHKERRQ(ierr);
355c4762a1bSJed Brown   ierr = DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx);CHKERRQ(ierr);
356c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
357c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   ierr = DMProjectFunction(dm, t, ctx.exactFuncs, NULL, INSERT_ALL_VALUES, u);CHKERRQ(ierr);
360c4762a1bSJed Brown   {
361c4762a1bSJed Brown     SNES snes;
362c4762a1bSJed Brown     KSP  ksp;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown     ierr = TSGetSNES(ts, &snes);CHKERRQ(ierr);
365c4762a1bSJed Brown     ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
366798534f6SMatthew G. Knepley     ierr = KSPMonitorSet(ksp, MonitorError, &ctx, NULL);CHKERRQ(ierr);
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown   ierr = TSSolve(ts, u);CHKERRQ(ierr);
369c4762a1bSJed Brown   ierr = VecViewFromOptions(u, NULL, "-sol_vec_view");CHKERRQ(ierr);
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
372c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
373c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
374c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
375c4762a1bSJed Brown   ierr = PetscFree(ctx.exactFuncs);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = PetscFinalize();
377c4762a1bSJed Brown   return ierr;
378c4762a1bSJed Brown }
379c4762a1bSJed Brown 
380c4762a1bSJed Brown /*TEST
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   # Full solves
383c4762a1bSJed Brown   test:
384c4762a1bSJed Brown     suffix: 2d_p1p1_r1
385c4762a1bSJed Brown     requires: triangle
386c4762a1bSJed 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
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   test:
389c4762a1bSJed Brown     suffix: 2d_p1p1_sor_r1
390c4762a1bSJed Brown     requires: triangle !single
391c4762a1bSJed 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
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   test:
394c4762a1bSJed Brown     suffix: 2d_p1p1_mg_r1
395c4762a1bSJed Brown     requires: triangle !single
396c4762a1bSJed 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
397c4762a1bSJed Brown 
398c4762a1bSJed Brown TEST*/
399