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