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