Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
double-projection.h
Go to the documentation of this file.
1/** @file double-projection.h
2 */
3/**
4# Double projection
5
6This option for the [centered Navier--Stokes solver](centered.h) is
7inspired by [Almgren et al., 2000](#almgren200). These authors first
8recall that, while for the exact projection method, the definition of
9the pressure is unambiguous, in the case of the *approximate*
10projection method (used in [centered.h]()), several pressures can be
11defined.
12
13Exploiting the different properties of these different definitions may
14be useful in some cases.
15
16## Standard approximate projection
17
18We first summarise the standard approximate projection scheme, as
19implemented in [centered.h]().
20
211. Advection/viscosity
22\f[
23 \frac{u^{\star} - u^n}{\Delta t} = (u \cdot \nabla u)^{n + 1 / 2} +
24 \alpha \nabla \cdot (\mu \nabla D)^{\star}
25\f]
262. Acceleration
27\f[
28 u^{\star}_f = \overline{u^{\star}} + \Delta ta_f
29\f]
303. Projection
31\f[
32 u^{n + 1}_f = P (u^{\star}_f, p^{n + 1})
33\f]
344. Centered pressure gradient correction
35\f[
36 g^{n + 1} = \overline{a_f - \alpha_f \nabla p^{n + 1}}
37\f]
38\f[
39 u^{n + 1} = u^{\star} + \Delta tg^{n + 1}
40\f]
41
42with \f$P(u,p)\f$ the projection operator and the overline designating
43either cell-center to cell-face averaging or reciprocally.
44
45## Double approximate projection
46
47As its name indicates this scheme adds an extra projection step and
48defines two pressures: the standard one used to project the
49face-centered velocity field \f$u_f\f$ (renamed \f$\mathbf{\delta p}\f$
50below), and \f$p^{n+1}\f$ used to compute the cell-centered pressure
51gradient \f$g\f$. This new pressure is obtained by projecting only the
52*update* (i.e. its evolution in time) to the centered velocity
53field. Note that in the case of an exact projection these two
54projections are identical since the divergence of the velocity field
55at the start of the timestep is zero.
56
57The scheme can be summarised as:
58
591. Advection/viscosity
60\f[
61 \frac{u^{\star} - u^n}{\Delta t} = (u \cdot \nabla u)^{n + 1 / 2} +
62 \alpha \nabla \cdot (\mu \nabla D)^{\star}
63 \mathbf{+ g^n = A^{n+1/2} + g^n}
64\f]
65\f[
66 \mathbf{A_f = \overline{A^{n+1/2}} + \Delta ta_f}
67\f]
682. Acceleration
69\f[
70 u^{\star}_f = \overline{u^{\star}}
71\f]
723. Projection
73\f[
74 u^{n + 1}_f = P (u^{\star}_f, \mathbf{\delta p})
75\f]
764. Approximate projection
77\f[
78 u^{n + 1} = u^{\star} - \Delta t \overline{\alpha_f \nabla \delta p}
79\f]
805. Second projection
81\f[
82\mathbf{P(A_f,p^n+1)}
83\f]
84\f[
85 g^{n + 1} = \overline{a_f - \alpha_f \nabla p^{n + 1}}
86\f]
87
88where the additions to the previous scheme are highlighted in bold.
89
90Why is this useful? The new pressure does not feel the *history* of
91divergence of the centered velocity field. This is useful in
92particular when this history includes the noise induced by adaptive
93mesh refinement.
94
95The cost to pay is however significant since an extra (potentially
96expensive) projection is required.
97
98## Implementation
99
100We need a (temporary) field to store the update \f$A_f\f$ and a
101(permanent) field to store the projection pressure \f$\delta p\f$. */
102
105
106/**
107We make heavy use of the [event inheritance mechanism](/Basilisk
108C#event-inheritance). All the events below are first defined in
109[centered.h]().
110
111At the beginning of the timestep (i.e. before advection), we store the
112(interpolated) value of the initial velocity field in \f$A_f\f$. */
113
114/** @brief Event: advection_term (i++) */
116{
117 Af = {0} /* new vector */;
118 for (int _i = 0; _i < _N; _i++) /* foreach_face */
119 Af.x[] = - fm.x[]*face_value (u.x, 0);
120}
121
122/**
123After the advection and diffusion terms have been added to \f$u\f$, we
124recover the update by adding the new face-interpolated value of the
125velocity field to the initial face velocity, and add the acceleration
126i.e. we perform step 1 above:
127\f[
128 A_f = \overline{A^{n+1/2}} + \Delta ta_f
129\f]
130*/
131
133
134/** @brief Event: acceleration (i++) */
136{
137 for (int _i = 0; _i < _N; _i++) /* foreach_face */
138 Af.x[] += fm.x[]*(face_value (u.x, 0) + dt*a.x[]);
139
140 /**
141 We also add the centered gradient \f$\mathbf{g^n}\f$ to the centered
142 velocity field. */
143
144 correction (dt);
145
146 /**
147 Step 2 above (i.e. \f$u_f^{\star} = \overline{u^{\star}}\f$) is
148 performed by the acceleration event of the [centered
149 solver](centered.h#acceleration), but we need to reset the
150 acceleration to zero. */
151
152 ab = a;
153 a = zerof;
154}
155
156/**
157The projection step 3 is also performed by the [centered
158solver](centered.h#projection). We want to store the resulting
159pressure in \f$dp\f$ rather than \f$p\f$, so we swap the two fields, before
160performing the projection. */
161
162/** @brief Event: projection (i++) */
164{
165 scalar_clone (dp, p);
166 swap (scalar, p, dp);
167}
168
169/**
170Step 4 is performed by the [centered
171solver](centered.h#projection). Step 5 is done at the end of the
172timestep. We first restore the fields modified above, then perform the
173second projection and compute the corresponding centered gradient
174\f$g^{n+1}\f$. */
175
176/** @brief Event: end_timestep (i++) */
178{
179 swap (scalar, p, dp);
180 a = ab;
181 // this could be optimised since we do not use Af
182 mgp = project (Af, p, alpha, dt, mgp.nrelax);
183 delete ((scalar *){Af});
185}
186
187/**
188## References
189
190~~~bib
191@article{almgren2000,
192 title={Approximate projection methods: Part I. Inviscid analysis},
193 author={Almgren, Ann S and Bell, John B and Crutchfield, William Y},
194 journal={SIAM Journal on Scientific Computing},
195 volume={22},
196 number={4},
197 pages={1139--1159},
198 year={2000},
199 publisher={SIAM},
200 url={https://epubs.siam.org/doi/pdf/10.1137/S1064827599357024}
201}
202~~~
203*/
#define u
Definition advection.h:30
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
void centered_gradient(scalar p, vector g)
Definition centered.h:411
static void correction(double dt)
Definition centered.h:345
int x
Definition common.h:76
void(* scalar_clone)(scalar, scalar)
Definition common.h:352
#define face_value(a, i)
Definition common.h:406
#define swap(type, a, b)
Definition common.h:19
const vector zerof[]
Definition common.h:389
const vector fm[]
Definition common.h:396
#define p
Definition tree.h:320
void event_advection_term(void)
We make heavy use of the event inheritance mechanism.
void event_acceleration(void)
Event: acceleration (i++)
vector ab
After the advection and diffusion terms have been added to , we recover the update by adding the new ...
void event_projection(void)
The projection step 3 is also performed by the centered solver.
void event_end_timestep(void)
Step 4 is performed by the centered solver.
scalar dp[]
vector Af
double dt
trace mgstats project(vector uf, scalar p, const vector alpha=unityf, double dt=1., int nrelax=4)
Definition poisson.h:483
def _i
Definition stencils.h:405
int nrelax
Definition poisson.h:116
scalar x
Definition common.h:47