Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
poisson.h
Go to the documentation of this file.
1/** @file poisson.h
2 */
3/**
4# Multigrid Poisson--Helmholtz solvers
5
6We want to solve Poisson--Helmholtz equations of the general form
7\f[
8L(a) = \nabla\cdot (\alpha\nabla a) + \lambda a = b
9\f]
10This can be done efficiently using a multigrid solver.
11
12An important aspect of Poisson--Helmholtz equations is that the
13operator \f$L()\f$ is linear. This property can be used to build better
14estimates of a solution by successive *corrections* to an initial
15guess. If we define an approximate solution \f$\tilde{a}\f$ as
16\f[
17\tilde{a} + da = a
18\f]
19where \f$a\f$ is the exact (unknown) solution, using the linearity of the
20operator we find that \f$da\f$ verifies
21\f[
22L(da) = b - L(\tilde{a})
23\f]
24where the right-hand-side is often called the *residual* of the
25approximate solution \f$\tilde{a}\f$.
26
27## Multigrid cycle
28
29Here we implement the multigrid cycle proper. Given an initial guess
30*a*, a residual *res*, a correction field *da* and a relaxation
31function *relax*, we will provide an improved guess at the end of the
32cycle. */
33
35void mg_cycle (scalar * a, scalar * res, scalar * da,
36 void (* relax) (scalar * da, scalar * res,
37 int depth, void * data),
38 void * data,
39 int nrelax, int minlevel, int maxlevel)
40{
41
42 /**
43 We first define the residual on all levels. */
44
45 restriction (res);
46
47 /**
48 We then proceed from the coarsest grid (*minlevel*) down to the
49 finest grid. */
50
51 minlevel = min (minlevel, maxlevel);
52 for (int l = minlevel; l <= maxlevel; l++) {
53
54 /**
55 On the coarsest grid, we take zero as initial guess. */
56
57 if (l == minlevel)
58 for (int _i = 0; _i < l; _i++)
59 for (int _s = 0; _s < 1; _s++) /* scalar in da */
60 /* foreach_blockf */
61 s[] = 0.;
62
63 /**
64 On all other grids, we take as initial guess the approximate solution
65 on the coarser grid bilinearly interpolated onto the current grid. */
66
67 else {
68 boundary_level (da, l - 1);
69 for (int _l = 0; _l < l; _l++)
70 for (int _s = 0; _s < 1; _s++) /* scalar in da */
71 /* foreach_blockf */
72 s[] = bilinear (point, s);
73 }
74
75 /**
76 We then apply homogeneous boundary conditions and do several
77 iterations of the relaxation function to refine the initial guess. */
78
79 for (int i = 0; i < nrelax; i++) {
81 relax (da, res, l, data);
82 }
83 }
84
85 /**
86 And finally we apply the resulting correction to *a*. */
87
88 for (int _i = 0; _i < _N; _i++) /* foreach */ {
89 scalar s, ds;
90 for (s, ds in a, da)
91 /* foreach_blockf */
92 s[] += ds[];
93 }
94}
95
96/**
97## Multigrid solver
98
99The multigrid solver itself uses successive calls to the multigrid
100cycle to refine an initial guess until a specified tolerance is
101reached.
102
103The maximum number of iterations is controlled by *NITERMAX* and the
104tolerance by *TOLERANCE* with the default values below. */
105
106int NITERMAX = 100, NITERMIN = 1;
107double TOLERANCE = 1e-3 [*];
108
109/**
110Information about the convergence of the solver is returned in a structure. */
111
112typedef struct {
113 int i; // number of iterations
114 double resb, resa; // maximum residual before and after the iterations
115 double sum; // sum of r.h.s.
116 int nrelax; // number of relaxations
117 int minlevel; // minimum level of the multigrid hierarchy
118} mgstats;
119
120/**
121The user needs to provide a function which computes the residual field
122(and returns its maximum) as well as the relaxation function. The
123user-defined pointer *data* can be used to pass arguments to these
124functions. The optional number of relaxations is *nrelax* and *res* is
125an optional list of fields used to store the residuals. The minimum
126level of the hierarchy can be set (default is zero i.e. the root
127cell). */
128
129trace
131 double (* residual) (scalar * a, scalar * b, scalar * res,
132 void * data),
133 void (* relax) (scalar * da, scalar * res, int depth,
134 void * data),
135 void * data = NULL,
136 int nrelax = 4,
137 scalar * res = NULL,
138 int minlevel = 0,
139 double tolerance = TOLERANCE)
140{
141
142 /**
143 We allocate a new correction and residual field for each of the scalars
144 in *a*. */
145
146 scalar * da = list_clone (a), * pres = res;
147 if (!res)
148 res = list_clone (b);
149
150 /**
151 The boundary conditions for the correction fields are the
152 *homogeneous* equivalent of the boundary conditions applied to
153 *a*. */
154
155 for (int b = 0; b < nboundary; b++)
156 for (int _s = 0; _s < 1; _s++) /* scalar in da */
157 s.boundary[b] = s.boundary_homogeneous[b];
158
159 /**
160 We initialise the structure storing convergence statistics. */
161
162 mgstats s = {0};
163 double sum = 0.;
164 scalar rhs = b[0];
165 foreach (reduction(+:sum))
166 sum += rhs[];
167 s.sum = sum;
168 s.nrelax = nrelax > 0 ? nrelax : 4;
169
170 /**
171 Here we compute the initial residual field and its maximum. */
172
173 double resb;
174 resb = s.resb = s.resa = (* residual) (a, b, res, data);
175
176 /**
177 We then iterate until convergence or until *NITERMAX* is reached. Note
178 also that we force the solver to apply at least one cycle, even if the
179 initial residual is lower than *TOLERANCE*. */
180
181 for (s.i = 0;
182 s.i < NITERMAX && (s.i < NITERMIN || s.resa > tolerance);
183 s.i++) {
184 mg_cycle (a, res, da, relax, data,
185 s.nrelax,
186 minlevel,
187 grid->maxdepth);
188 s.resa = (* residual) (a, b, res, data);
189
190 /**
191 We tune the number of relaxations so that the residual is reduced
192 by between 2 and 20 for each cycle. This is particularly useful
193 for stiff systems which may require a larger number of relaxations
194 on the finest grid. */
195
196#if 1
197 if (s.resa > tolerance) {
198 if (resb/s.resa < 1.2 && s.nrelax < 100)
199 s.nrelax++;
200 else if (resb/s.resa > 10 && s.nrelax > 2)
201 s.nrelax--;
202 }
203#else
204 if (s.resa == resb) /* convergence has stopped!! */
205 break;
206 if (s.resa > resb/1.1 && p.minlevel < grid->maxdepth)
207 p.minlevel++;
208#endif
209
210 resb = s.resa;
211 }
212 s.minlevel = minlevel;
213
214 /**
215 If we have not satisfied the tolerance, we warn the user. */
216
217 if (s.resa > tolerance) {
218 scalar v = a[0]; // fixme: should not be necessary
219 fprintf (ferr,
220 "src/poisson.h:%d: warning: convergence for %s not reached after %d iterations\n"
221 " res: %g sum: %g nrelax: %d tolerance: %g\n", LINENO, v.name,
222 s.i, s.resa, s.sum, s.nrelax, tolerance), fflush (ferr);
223 }
224
225 /**
226 We deallocate the residual and correction fields and free the lists. */
227
228 if (!pres)
229 delete (res), free (res);
230 delete (da), free (da);
231
232 return s;
233}
234
235/**
236## Application to the Poisson--Helmholtz equation
237
238We now apply the generic multigrid solver to the Poisson--Helmholtz equation
239\f[
240\nabla\cdot (\alpha\nabla a) + \lambda a = b
241\f]
242We first setup the data structure required to pass the extra
243parameters \f$\alpha\f$ and \f$\lambda\f$. We define \f$\alpha\f$ as a face
244vector field because we need values at the face locations
245corresponding to the face gradients of field \f$a\f$.
246
247*alpha* and *lambda* are declared as *(const)* to indicate that the
248function works also when *alpha* and *lambda* are constant vector
249(resp. scalar) fields. If *tolerance* is set, it supersedes the
250default *TOLERANCE* of the multigrid solver, *nrelax* controls the
251initial number of relaxations (default is one), *minlevel* controls
252the minimum level of the hierarchy (default is one) and *res* is an
253optional list of fields used to store the final residual (which can be
254useful to monitor convergence). */
255
256struct Poisson {
260 double tolerance;
263#if EMBED
264 double (* embed_flux) (Point, scalar, vector, double *);
265#endif
266};
267
268/**
269We can now write the relaxation function. We first recover the extra
270parameters from the data pointer. */
271
272static void relax (scalar * al, scalar * bl, int l, void * data)
273{
274 scalar a = al[0], b = bl[0];
275 struct Poisson * p = (struct Poisson *) data;
276 const vector alpha = p->alpha;
277 const scalar lambda = p->lambda;
278
279 /**
280 We use either Jacobi (under)relaxation, Gauss-Seidel or we directly
281 reuse values as soon as they are updated. For Jacobi, we need to
282 allocate space for the new field *c*. Jacobi is useful mostly as it
283 gives results which are independent of the order in which the cells
284 are traversed. This is not the case for the simple traversal, which
285 means for example that results will depend on whether a tree or a
286 multigrid is used (because cells will be traversed in a different
287 order). The same comment applies to OpenMP or MPI parallelism. In
288 practice however Jacobi convergence tends to be slower than simple
289 reuse. */
290
291#if JACOBI
292 scalar c[];
293#else
294 scalar c = a;
295#endif
296
297 /**
298 On GPUs, we use red/black Gauss-Seidel relaxation, which requires
299 two loops (for odd/even indices). Note also that, unlike the other
300 option, red/black relaxation should be deterministic. */
301
302#if GAUSS_SEIDEL || _GPU
303 for (int parity = 0; parity < 2; parity++)
304 for (int _i = 0; _i < l, nowarning; _i++)
305 if (level == 0 || ((point.i + parity) % 2) != (point.j % 2))
306#else
307 for (int _i = 0; _i < l, nowarning; _i++)
308#endif
309 {
310
311 /**
312 We use the face values of \f$\alpha\f$ to weight the gradients of the
313 5-points Laplacian operator. We get the relaxation function. */
314
315 double n = - sq(Delta)*b[], d = - lambda[]*sq(Delta);
316 for (int _d = 0; _d < dimension; _d++) {
317 n += alpha.x[1]*a[1] + alpha.x[]*a[-1];
318 d += alpha.x[1] + alpha.x[];
319 }
320#if EMBED
321 if (p->embed_flux) {
322 double c, e = p->embed_flux (point, a, alpha, &c);
323 n -= c*sq(Delta);
324 d += e*sq(Delta);
325 }
326 if (!d)
327 c[] = 0., b[] = 0.;
328 else
329#endif // EMBED
330 c[] = n/d;
331 }
332
333 /**
334 For weighted Jacobi we under-relax with a weight of 2/3. */
335
336#if JACOBI
337 for (int _i = 0; _i < l; _i++)
338 a[] = (a[] + 2.*c[])/3.;
339#endif
340
341#if TRASH
342 scalar a1[];
343 for (int _i = 0; _i < l; _i++)
344 a1[] = a[];
345 /* trash: a */;
346 for (int _i = 0; _i < l; _i++)
347 a[] = a1[];
348#endif
349}
350
351/**
352The equivalent residual function is obtained in a similar way in the
353case of a Cartesian grid, however the case of the tree mesh
354requires more careful consideration... */
355
356static double residual (scalar * al, scalar * bl, scalar * resl, void * data)
357{
358 scalar a = al[0], b = bl[0], res = resl[0];
359 struct Poisson * p = (struct Poisson *) data;
360 const vector alpha = p->alpha;
361 const scalar lambda = p->lambda;
362 double maxres = 0.;
363#if TREE
364 /* conservative coarse/fine discretisation (2nd order) */
365 vector g[];
366 for (int _i = 0; _i < _N; _i++) /* foreach_face */
367 g.x[] = alpha.x[]*face_gradient_x (a, 0);
368 foreach (reduction(max:maxres), nowarning) {
369 res[] = b[] - lambda[]*a[];
370 for (int _d = 0; _d < dimension; _d++)
371 res[] -= (g.x[1] - g.x[])/Delta;
372#if EMBED
373 if (p->embed_flux) {
374 double c, e = p->embed_flux (point, a, alpha, &c);
375 res[] += c - e*a[];
376 }
377#endif // EMBED
378 if (fabs (res[]) > maxres)
379 maxres = fabs (res[]);
380 }
381#else // !TREE
382 /* "naive" discretisation (only 1st order on trees) */
383 foreach (reduction(max:maxres), nowarning) {
384 res[] = b[] - lambda[]*a[];
385 for (int _d = 0; _d < dimension; _d++)
386 res[] += (alpha.x[0]*face_gradient_x (a, 0) -
387 alpha.x[1]*face_gradient_x (a, 1))/Delta;
388#if EMBED
389 if (p->embed_flux) {
390 double c, e = p->embed_flux (point, a, alpha, &c);
391 res[] += c - e*a[];
392 }
393#endif // EMBED
394 if (fabs (res[]) > maxres)
395 maxres = fabs (res[]);
396 }
397#endif // !TREE
398 return maxres;
399}
400
401/**
402## User interface
403
404Finally we provide a generic user interface for a Poisson--Helmholtz
405equation of the form
406\f[
407\nabla\cdot (\alpha\nabla a) + \lambda a = b
408\f] */
409
410trace
412 const vector alpha = {{-1}},
413 const scalar lambda = {-1},
414 double tolerance = 0.,
415 int nrelax = 4,
416 int minlevel = 0,
417 scalar * res = NULL,
418 double (* flux) (Point, scalar, vector, double *) = NULL)
419{
420
421 /**
422 If \f$\alpha\f$ or \f$\lambda\f$ are not set, we replace them with constant
423 unity vector (resp. zero scalar) fields. Note that the user is free to
424 provide \f$\alpha\f$ and \f$\beta\f$ as constant fields. */
425
426 if (alpha.x.i < 0)
427 alpha[] = {1.,1.,1.};
428 if (lambda.i < 0)
429 lambda[] = 0.;
430
431
432
433 /**
434 We need \f$\alpha\f$ and \f$\lambda\f$ on all levels of the grid. */
435
437
438 /**
439 If *tolerance* is set it supersedes the default of the multigrid
440 solver. */
441
442 double defaultol = TOLERANCE;
443 if (tolerance)
445
446 struct Poisson p = {a, b, alpha, lambda, tolerance, nrelax, minlevel, res };
447#if EMBED
448 if (!flux && a.boundary[embed] != symmetry)
449 p.embed_flux = embed_flux;
450 else
451 p.embed_flux = flux;
452#endif // EMBED
453 mgstats s = mg_solve ({a}, {b}, residual, relax, &p,
454 nrelax, res, max(1, minlevel));
455
456 /**
457 We restore the default. */
458
459 if (tolerance)
461
462 return s;
463}
464
465/**
466## Projection of a velocity field
467
468The function below "projects" the velocity field *u* onto the space of
469divergence-free velocity fields i.e.
470\f[
471\mathbf{u}_f^{n+1} \leftarrow \mathbf{u}_f - \Delta t\alpha\nabla p
472\f]
473so that
474\f[
475\nabla\cdot\mathbf{u}_f^{n+1} = 0
476\f]
477This gives the Poisson equation for the pressure
478\f[
479\nabla\cdot(\alpha\nabla p) = \frac{\nabla\cdot\mathbf{u}_f}{\Delta t}
480\f] */
481
482trace
484 const vector alpha = unityf,
485 double dt = 1.,
486 int nrelax = 4)
487{
488
489 /**
490 We allocate a local scalar field and compute the divergence of
491 \f$\mathbf{u}_f\f$. The divergence is scaled by *dt* so that the
492 pressure has the correct dimension. */
493
494 scalar div[];
495 for (int _i = 0; _i < _N; _i++) /* foreach */ {
496 div[] = 0.;
497 for (int _d = 0; _d < dimension; _d++)
498 div[] += uf.x[1] - uf.x[];
499 div[] /= dt*Delta;
500 }
501
502 /**
503 We solve the Poisson problem. The tolerance (set with *TOLERANCE*) is
504 the maximum relative change in volume of a cell (due to the divergence
505 of the flow) during one timestep i.e. the non-dimensional quantity
506 \f[
507 |\nabla\cdot\mathbf{u}_f|\Delta t
508 \f]
509 Given the scaling of the divergence above, this gives */
510
513
514 /**
515 And compute \f$\mathbf{u}_f^{n+1}\f$ using \f$\mathbf{u}_f\f$ and \f$p\f$. */
516
517 for (int _i = 0; _i < _N; _i++) /* foreach_face */
518 uf.x[] -= dt*alpha.x[]*face_gradient_x (p, 0);
519
520 return mgp;
521}
vector uf[]
We allocate the (face) velocity field.
Definition advection.h:29
vector g[]
We store the combined pressure gradient and acceleration field in g*.
Definition all-mach.h:65
const vector alpha
Definition all-mach.h:47
const vector a
Definition all-mach.h:59
mgstats mgp
Definition all-mach.h:66
int min
Definition balance.h:192
#define dimension
Definition bitree.h:3
scalar * list_clone(scalar *l)
void(* boundary_level)(scalar *, int l)
static double symmetry(Point point, Point neighbor, scalar s, bool *data)
define l
free(list1)
int x
Definition common.h:76
#define face_gradient_x(a, i)
Definition common.h:403
static number sq(number x)
Definition common.h:11
Grid * grid
Definition common.h:32
#define vector(x)
Definition common.h:354
int nboundary
Definition common.h:127
const vector unityf[]
Definition common.h:390
#define p
Definition tree.h:320
#define LINENO
Definition config.h:105
Point point
Definition conserving.h:86
else return n
Definition curvature.h:101
double dt
scalar s
Definition embed-tree.h:56
double embed_flux(Point point, scalar s, vector mu, double *val)
Definition embed.h:657
double v[2]
Definition embed.h:383
double d[2]
Definition embed.h:383
scalar int i
Definition embed.h:74
bid embed
Definition embed.h:463
#define depth()
Definition cartesian.h:14
double * sum
#define data(k, l)
Definition linear.h:26
size_t max
Definition mtrace.h:8
void(* restriction)(Point, scalar)
static double bilinear(Point point, scalar s)
size *double * b
int int int level
int NITERMAX
Definition poisson.h:106
int NITERMIN
Definition poisson.h:106
trace mgstats mg_solve(scalar *a, scalar *b, double(*residual)(scalar *a, scalar *b, scalar *res, void *data), void(*relax)(scalar *da, scalar *res, int depth, void *data), void *data=NULL, int nrelax=4, scalar *res=NULL, int minlevel=0, double tolerance=TOLERANCE)
The user needs to provide a function which computes the residual field (and returns its maximum) as w...
Definition poisson.h:130
trace void mg_cycle(scalar *a, scalar *res, scalar *da, void(*relax)(scalar *da, scalar *res, int depth, void *data), void *data, int nrelax, int minlevel, int maxlevel)
Definition poisson.h:35
static void relax(scalar *al, scalar *bl, int l, void *data)
We can now write the relaxation function.
Definition poisson.h:272
static double residual(scalar *al, scalar *bl, scalar *resl, void *data)
The equivalent residual function is obtained in a similar way in the case of a Cartesian grid,...
Definition poisson.h:356
trace mgstats project(vector uf, scalar p, const vector alpha=unityf, double dt=1., int nrelax=4)
Definition poisson.h:483
double TOLERANCE
Definition poisson.h:107
def _i
Definition stencils.h:405
int maxdepth
Definition common.h:30
Definition linear.h:21
int i
Definition linear.h:23
int nrelax
Definition poisson.h:261
const scalar lambda
Definition poisson.h:259
scalar * res
Definition poisson.h:262
double tolerance
Definition poisson.h:260
int minlevel
Definition poisson.h:261
scalar a
Definition poisson.h:257
const vector alpha
Definition poisson.h:258
scalar b
Definition poisson.h:257
Information about the convergence of the solver is returned in a structure.
Definition poisson.h:112
int nrelax
Definition poisson.h:116
int minlevel
Definition poisson.h:117
double sum
Definition poisson.h:115
double resa
Definition poisson.h:114
int i
Definition poisson.h:113
int i
Definition common.h:44
scalar x
Definition common.h:47
#define poisson(...)
Finally, we overload the poisson() function called by [two-phase.h]() with our new function.
Definition thermal.h:243
#define lambda
scalar flux[]
Definition vof.h:166
scalar c
Definition vof.h:57
src vof fflush(ferr)