|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Functions | |
| void | event_advection_term (void) |
| We make heavy use of the event inheritance mechanism. | |
| void | event_acceleration (void) |
| Event: acceleration (i++) | |
| 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. | |
Variables | |
| vector | Af |
| scalar | dp [] |
| vector | ab |
| After the advection and diffusion terms have been added to \(u\), we recover the update by adding the new face-interpolated value of the velocity field to the initial face velocity, and add the acceleration i.e. | |
Event: acceleration (i++)
We also add the centered gradient \(\mathbf{g^n}\) to the centered velocity field.
Step 2 above (i.e. \(u_f^{\star} = \overline{u^{\star}}\)) is performed by the acceleration event of the centered solver, but we need to reset the acceleration to zero.
Definition at line 135 of file double-projection.h.
References _i, a, ab, Af, correction(), dt, face_value, fm, u, vector::x, x, and zerof.
We make heavy use of the event inheritance mechanism.
All the events below are first defined in [centered.h]().
At the beginning of the timestep (i.e. before advection), we store the (interpolated) value of the initial velocity field in \(A_f\).
Event: advection_term (i++)
Definition at line 115 of file double-projection.h.
Step 4 is performed by the centered solver.
Step 5 is done at the end of the timestep. We first restore the fields modified above, then perform the second projection and compute the corresponding centered gradient \(g^{n+1}\).
Event: end_timestep (i++)
Definition at line 177 of file double-projection.h.
References a, ab, Af, alpha, centered_gradient(), dp, dt, g, mgp, mgstats::nrelax, p, project(), and swap.
The projection step 3 is also performed by the centered solver.
We want to store the resulting pressure in \(dp\) rather than \(p\), so we swap the two fields, before performing the projection.
Event: projection (i++)
Definition at line 163 of file double-projection.h.
References dp, p, scalar_clone, and swap.
| vector ab |
After the advection and diffusion terms have been added to \(u\), we recover the update by adding the new face-interpolated value of the velocity field to the initial face velocity, and add the acceleration i.e.
we perform step 1 above:
\[ A_f = \overline{A^{n+1/2}} + \Delta ta_f \]
Definition at line 132 of file double-projection.h.
Referenced by distance(), event_acceleration(), event_end_timestep(), and update_distance().
| vector Af |
This option for the centered Navier–Stokes solver is inspired by Almgren et al., 2000. These authors first recall that, while for the exact projection method, the definition of the pressure is unambiguous, in the case of the approximate projection method (used in [centered.h]()), several pressures can be defined.
Exploiting the different properties of these different definitions may be useful in some cases.
We first summarise the standard approximate projection scheme, as implemented in [centered.h]().
\[ \frac{u^{\star} - u^n}{\Delta t} = (u \cdot \nabla u)^{n + 1 / 2} + \alpha \nabla \cdot (\mu \nabla D)^{\star} \]
\[ u^{\star}_f = \overline{u^{\star}} + \Delta ta_f \]
\[ u^{n + 1}_f = P (u^{\star}_f, p^{n + 1}) \]
\[ g^{n + 1} = \overline{a_f - \alpha_f \nabla p^{n + 1}} \]
\[ u^{n + 1} = u^{\star} + \Delta tg^{n + 1} \]
with \(P(u,p)\) the projection operator and the overline designating either cell-center to cell-face averaging or reciprocally.
As its name indicates this scheme adds an extra projection step and defines two pressures: the standard one used to project the face-centered velocity field \(u_f\) (renamed \(\mathbf{\delta p}\) below), and \(p^{n+1}\) used to compute the cell-centered pressure gradient \(g\). This new pressure is obtained by projecting only the update* (i.e. its evolution in time) to the centered velocity field. Note that in the case of an exact projection these two projections are identical since the divergence of the velocity field at the start of the timestep is zero.
The scheme can be summarised as:
\[ \frac{u^{\star} - u^n}{\Delta t} = (u \cdot \nabla u)^{n + 1 / 2} + \alpha \nabla \cdot (\mu \nabla D)^{\star} \mathbf{+ g^n = A^{n+1/2} + g^n} \]
\[ \mathbf{A_f = \overline{A^{n+1/2}} + \Delta ta_f} \]
\[ u^{\star}_f = \overline{u^{\star}} \]
\[ u^{n + 1}_f = P (u^{\star}_f, \mathbf{\delta p}) \]
\[ u^{n + 1} = u^{\star} - \Delta t \overline{\alpha_f \nabla \delta p} \]
\[ \mathbf{P(A_f,p^n+1)} \]
\[ g^{n + 1} = \overline{a_f - \alpha_f \nabla p^{n + 1}} \]
where the additions to the previous scheme are highlighted in bold.
Why is this useful? The new pressure does not feel the history of divergence of the centered velocity field. This is useful in particular when this history includes the noise induced by adaptive mesh refinement.
The cost to pay is however significant since an extra (potentially expensive) projection is required.
We need a (temporary) field to store the update \(A_f\) and a (permanent) field to store the projection pressure \(\delta p\).
Definition at line 103 of file double-projection.h.
Referenced by event_acceleration(), event_advection_term(), and event_end_timestep().
| scalar dp[] |
Definition at line 104 of file double-projection.h.
Referenced by event_acceleration(), event_end_timestep(), event_pressure(), event_projection(), event_viscous_term(), is_local_prolongation(), and segment_flux().