|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
Go to the source code of this file.
Data Structures | |
| struct | Thermal |
| We then define the corresponding relaxation and residual functions which will be used by the multigrid solver. More... | |
Macros | |
| #define | kappa(f) (clamp(f,0.,1.)*(kappa1 - kappa2) + kappa2) |
| #define | poisson(...) poisson_thermal(__VA_ARGS__) |
| Finally, we overload the poisson() function called by [two-phase.h]() with our new function. | |
Functions | |
| double | average_temperature (Point point, double p) |
| These functions are provided by the Equation Of State. | |
| double | thermal_expansion (Point point) |
| void | event_defaults (void) |
| If the thermal conductivities are not-zero, we need to allocate a new field. | |
| static void | relax_thermal (scalar *al, scalar *bl, int l, void *data) |
| static double | residual_thermal (scalar *al, scalar *bl, scalar *resl, void *data) |
| mgstats | poisson_thermal (scalar p, scalar rhs, const vector alpha, const scalar lambda, double tolerance=0., int nrelax=4, int minlevel=0, scalar *res=NULL) |
| This is the interface for the coupled solver. | |
| void | event_acceleration (void) |
| Event: acceleration (i++) | |
| void | event_end_timestep (void) |
| Event: end_timestep (i++) | |
Variables | |
| scalar | T [] |
| double | kappa1 = 0. |
| ... and the thermal conductivities and specific heat capacities. | |
| double | kappa2 = 0. |
| double | cp1 = 0. |
| double | cp2 = 0. |
| const vector | kappa0 [] = {0,0,0} |
| const vector | kappa = kappa0 |
| scalar | Ts [] |
| scalar | lambda1 [] |
| scalar | lambda2 [] |
| scalar | lambda4 [] |
| #define poisson | ( | ... | ) | poisson_thermal(__VA_ARGS__) |
These functions are provided by the Equation Of State.
See for example the Noble-Able Stiffened-Gas Equation Of State.
Definition at line 99 of file NASG.h.
References b1, b2, clamp(), cp1, cp2, cv1, cv2, f, frho1, frho2, p, PI1, PI2, and x.
Referenced by event_acceleration().
Event: acceleration (i++)
This is done before the projection, which will call the coupled Poisson–Helmholtz solver.
We compute \(\lambda_1\), \(\lambda_2\), \(\lambda_4\) and r.h.s. for the temperature which will be used by the coupled Poisson–Helmholtz solver. The \(\lambda\) parameter and r.h.s. for the pressure are computed by the two-phase compressible solver.
We also compute the face-averaged thermal conductivity.
Definition at line 254 of file thermal.h.
References _i, average_temperature(), beta, cm, cp1, cp2, dt, f, fm, frho1, frho2, kappa, kappa1, kappa2, lambda1, lambda2, lambda4, point, ps, thermal_expansion(), Ts, vector::x, and x.
Event: end_timestep (i++)
This adds the \(-\nabla\cdot\mathbf{q}_i\) term of the energy equation.
Definition at line 292 of file thermal.h.
References _i, clamp(), cm, dimension, dt, energy(), f, face_gradient_x, fE1, fE2, kappa, kappa1, kappa2, T, vector::x, and x.
| mgstats poisson_thermal | ( | scalar | p, |
| scalar | rhs, | ||
| const vector | alpha, | ||
| const scalar | lambda, | ||
| double | tolerance = 0., |
||
| int | nrelax = 4, |
||
| int | minlevel = 0, |
||
| scalar * | res = NULL |
||
| ) |
This is the interface for the coupled solver.
Definition at line 215 of file thermal.h.
References alpha, kappa, lambda, lambda1, lambda2, lambda4, max, mg_solve(), Thermal::minlevel, Thermal::nrelax, p, relax_thermal(), Thermal::res, residual_thermal(), restriction, s, T, Thermal::tolerance, TOLERANCE, tp, Ts, and x.
| double cp1 = 0. |
Definition at line 66 of file thermal.h.
Referenced by average_temperature(), and event_acceleration().
| double cp2 = 0. |
Definition at line 66 of file thermal.h.
Referenced by average_temperature(), and event_acceleration().
| double kappa1 = 0. |
... and the thermal conductivities and specific heat capacities.
Definition at line 65 of file thermal.h.
Referenced by event_acceleration(), event_defaults(), and event_end_timestep().
| double kappa2 = 0. |
Definition at line 65 of file thermal.h.
Referenced by event_acceleration(), event_defaults(), and event_end_timestep().
| scalar lambda1[] |
Definition at line 117 of file thermal.h.
Referenced by event_acceleration(), poisson_thermal(), relax_thermal(), and residual_thermal().
| scalar lambda2[] |
Definition at line 117 of file thermal.h.
Referenced by event_acceleration(), poisson_thermal(), relax_thermal(), and residual_thermal().
| scalar lambda4[] |
Definition at line 117 of file thermal.h.
Referenced by event_acceleration(), poisson_thermal(), relax_thermal(), and residual_thermal().
| scalar T[] |
This extension of the two-phase compressible solver is described in detail in Saade et al, 2023.
This solves the two-phase compressible Navier-Stokes equations including the total energy equation.
\[ \frac{\partial (f \rho_i)}{\partial t } + \nabla \cdot (f \rho_i \mathbf{u}) = 0 \]
\[ \frac{\partial (\rho_i \mathbf{u})}{\partial t } + \nabla \cdot ( \rho_i \mathbf{u} \mathbf{u}) = -\nabla p + \nabla \cdot \tau_i' \]
\[ \frac{\partial [\rho_i (e_i + \mathbf{u}^2/2)]}{\partial t } + \nabla \cdot [ \rho_i \mathbf{u} (e_i + \mathbf{u}^2/2)] = -\nabla \cdot (\mathbf{u} p_i) + \nabla \cdot \left( \tau'_i \mathbf{u} \right) {\color{blue} - \nabla \cdot \mathbf{q}_i} \]
an advection equation for the volume fraction \(f\)
\[ \frac{\partial f}{\partial t} + \mathbf{u} \cdot \nabla f = 0 \]
The term in blue in the energy equation is the heat flux which is added by this extension. The other terms are treated by the two-phase compressible solver. The heat flux is given by Fourier's relation as
\[ \nabla \cdot \mathbf{q}_i = - \kappa_i \nabla T_i \]
where \(\kappa_i\) is the thermal conductivity and \(T\) is the temperature field.
The temperature and pressure fields are then coupled through the two equations (8 and 10 in Saade et al, 2023)
\[ \rho_ic_{p,i}\frac{DT_i}{Dt} = \beta_iT_i\frac{Dp_i}{Dt} - \nabla\cdot\mathbf{q}_i \]
with \(c_p\) the specific heat capacity and \(\beta\) the thermal expansion coefficient.
\[ \left(\frac{\gamma_i}{\rho_ic_i^2} - \frac{\beta_i^2T_i}{\rho_ic_{p,i}}\right)\frac{Dp_i}{Dt} = -\frac{\beta_i}{\rho_ic_{p,i}}\nabla\cdot\mathbf{q}_i - \nabla\cdot\mathbf{u}_i \]
with \(\gamma\) the ratio of specific heats.
The system is closed by specifying an Equation Of State (EOS) relating \(p\), \(\rho\) and \(T\).
In addition to the primary variables of the compressible two-phase solver we have the temperature field...
Definition at line 60 of file thermal.h.
Referenced by diagonalization_2D(), energy(), equilibrium_tide_constituents(), event_end_timestep(), poisson_thermal(), relax_thermal(), residual_thermal(), and RPE().
| scalar Ts[] |
The robustness of the solver relies on solving equations (8) and (10) in a strongly-coupled manner with the multigrid solver. After temporal discretisation this leads to the following coupled Poisson–Helmholtz system for the temperature \(T^{n+1}\) and pressure \(p^{n+1}\) (see eqs. 25 and 26 in Saade et al, 2023)
\[ \nabla\cdot\left(\kappa\nabla T^{n+1}\right) + \lambda_1 T^{n+1} + \lambda_2 p^{n+1} = \lambda_T \]
\[ \nabla\cdot\left(\alpha\nabla p^{n+1}\right) + \lambda p^{n+1} + \lambda_4\nabla\cdot\left(\kappa\nabla T^{n+1}\right) = \lambda_p \]
We allocate fields to store the provisionary temperature \(T^n\) and the coefficients.
Definition at line 117 of file thermal.h.
Referenced by event_acceleration(), poisson_thermal(), and thermal_expansion().