Euler 1D (CAF enabled) PDEs system field.
This object does not use FOODIE.
The 1D Euler PDEs system considered is a non linear, hyperbolic (inviscid) system of conservation laws for compressible gas dynamics, that reads as $$ \begin{matrix} U_t = R(U) \Leftrightarrow U_t = F(U)_x \\ U = \begin{bmatrix} \rho \\ \rho u \\ \rho E \end{bmatrix}\;\;\; F(U) = \begin{bmatrix} \rho u \\ \rho u^2 + p \\ \rho u H \end{bmatrix} \end{matrix} $$ where \(\rho\) is the density, \(u\) is the velocity, \(p\) the pressure, \(E\) the total internal specific energy and \(H\) the total specific enthalpy. The PDEs system must completed with the proper initial and boundary conditions. Moreover, an ideal (thermally and calorically perfect) gas is considered $$ \begin{matrix} R = c_p - c_v \\ \gamma = \frac{c_p}{c_v}\\ e = c_v T \\ h = c_p T \end{matrix} $$ where R is the gas constant, \(c_p\,c_v\) are the specific heats at constant pressure and volume (respectively), e is the internal energy, h is the internal enthalpy and T is the temperature. The following addition equations of state hold: $$ \begin{matrix} T = \frac{p}{\rho R} \\ E = \rho e + \frac{1}{2} \rho u^2 \\ H = \rho h + \frac{1}{2} \rho u^2 \\ a = \sqrt{\frac{\gamma p}{\rho}} \end{matrix} $$
An extension of the above Euler system is considered allowing the modelling of a multi-fluid mixture of different gas (with different physical characteristics). The well known Standard Thermodynamic Model is used to model the gas mixture replacing the density with the density fraction of each specie composing the mixture. This led to the following system: $$ \begin{matrix} U_t = R(U) \Leftrightarrow U_t = F(U)_x \\ U = \begin{bmatrix} \rho_s \\ \rho u \\ \rho E \end{bmatrix}\;\;\; F(U) = \begin{bmatrix} \rho_s u \\ \rho u^2 + p \\ \rho u H \end{bmatrix}\;\;\; for\; s=1,2,...N_s \\ \rho = \sum_{s=1}^{N_s}\rho_s \\ c_p = \sum_{s=1}^{N_S} \frac{\rho_s}{\rho} c_{p,s} \quad c_v = \sum_{s=1}^{N_S} \frac{\rho_s}{\rho} c_{v,s} \end{matrix} $$ where \(N_s\) is the number of initial species composing the gas mixture.
The finite volume, Godunov's like approach is employed. The conservative variables (and the primitive ones) are co-located at the cell center. The cell and (inter)faces numeration is as follow.
cell (inter)faces | | v v |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------| | 1-Ng | 2-Ng | ..... | -1 | 0 | 1 | 2 | ..... | Ni | Ni+1 | Ni+1 | ..... |Ni+Ng-1| Ni+Ng | |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------| 0-Ng -1 0 1 2 Ni-1 Ni Ni+Ng
Where Ni are the finite volumes (cells) used for discretizing the domain and Ng are the ghost cells used for imposing the left and right boundary conditions (for a total of 2Ng cells).
Primitive variables are organized as an array of reals which the first index means:
Conservative variables are organized as an array (rank 2) of reals which the first index means:
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
character(len=:), | private, | allocatable | :: | BC_L | Left boundary condition type. |
||
character(len=:), | private, | allocatable | :: | BC_R | Right boundary condition type. |
||
real(kind=R_P), | private | :: | Dx | = | 0._R_P | Space step. |
|
integer(kind=I_P), | private | :: | Nc | = | 0 | Number of conservative variables, Ns+2. |
|
integer(kind=I_P), | private | :: | Ng | = | 0 | Number of ghost cells for boundary conditions handling. |
|
integer(kind=I_P), | private | :: | Ni | = | 0 | Space dimension. |
|
integer(kind=I_P), | private | :: | Np | = | 0 | Number of primitive variables, Ns+4. |
|
integer(kind=I_P), | private | :: | Ns | = | 0 | Number of initial species. |
|
real(kind=R_P), | private, | allocatable | :: | U(:,:) | Integrand (state) variables, whole physical domain [1:Nc,1:Ni]. |
||
real(kind=R_P), | private, | allocatable | :: | cp0(:) | Specific heat cp of initial species [1:Ns]. |
||
real(kind=R_P), | private, | allocatable | :: | cv0(:) | Specific heat cv of initial species [1:Ns]. |
||
integer(kind=I_P), | private | :: | me | = | 0 | ID of this_image(). |
|
integer(kind=I_P), | private | :: | ord | = | 0 | Space accuracy formal order. |
|
integer(kind=I_P), | private | :: | we | = | 0 | Number of CAF images used. |
|
type(weno_interpolator_upwind), | private | :: | weno | WENO interpolator. |
Euler = Euler.
Assign one Euler field to another.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(inout) | :: | lhs | Left hand side. |
||
class(euler_1D_caf_nf), | intent(in) | :: | rhs | Right hand side. |
Overloading = assignament.
Assign one Euler field to another.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(inout) | :: | lhs | Left hand side. |
||
class(euler_1D_caf_nf), | intent(in) | :: | rhs | Right hand side. |
Convert conservative variables to primitive ones.
Convert conservative variables to primitive variables.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
real(kind=R_P), | intent(in) | :: | conservative(:) | Conservative variables. |
Primitive variables.
Destroy field.
Destroy field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(inout) | :: | self | Euler field. |
Compute the current time step, by means of CFL condition.
Compute the current time step by means of CFL condition.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
integer(kind=I_P), | intent(in) | :: | Nmax | Maximun number of iterates. |
||
real(kind=R_P), | intent(in) | :: | Tmax | Maximum time (ignored if Nmax>0). |
||
real(kind=R_P), | intent(in) | :: | t | Time. |
||
real(kind=R_P), | intent(in) | :: | CFL | CFL value. |
Time step.
Impose boundary conditions.
Impose boundary conditions.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
real(kind=R_P), | intent(inout) | :: | primitive(1:self%Np,1-self%Ng:self%Ni+self%Ng) | Primitive variables [1:Np,1-Ng:Ni+Ng]. |
Init field.
Init field.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(inout) | :: | self | Euler field. |
||
integer(kind=I_P), | intent(in) | :: | Ni | Space dimension (local image). |
||
integer(kind=I_P), | intent(in) | :: | Ns | Number of initial species. |
||
real(kind=R_P), | intent(in) | :: | Dx | Space step. |
||
character(len=*), | intent(in) | :: | BC_L | Left boundary condition type. |
||
character(len=*), | intent(in) | :: | BC_R | Right boundary condition type. |
||
real(kind=R_P), | intent(in) | :: | initial_state(:,:) | Initial state of primitive variables. |
||
real(kind=R_P), | intent(in) | :: | cp0(:) | Initial specific heat, constant pressure. |
||
real(kind=R_P), | intent(in) | :: | cv0(:) | Initial specific heat, constant volume. |
||
integer(kind=I_P), | intent(in) | :: | me | ID of this_image(). |
||
integer(kind=I_P), | intent(in) | :: | we | Number of CAF images used. |
||
integer(kind=I_P), | intent(in), | optional | :: | ord | Space accuracy formal order. |
Extract Euler field.
Output the Euler field state (primitive variables).
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
logical, | intent(in), | optional | :: | conservative | Output conservative variables instead of primitive. |
Euler state vector.
Convert primitive variables to conservative ones.
Convert primitive variables to conservative variables.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
real(kind=R_P), | intent(in) | :: | primitive(:) | Primitive variables. |
Conservative variables.
Reconstruct interfaces states.
Reconstruct the interfaces states (into primitive variables formulation) by the requested order of accuracy.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
real(kind=R_P), | intent(in) | :: | primitive(1:self%Np,1-self%Ng:self%Ni+self%Ng) | Primitive variables. |
||
real(kind=R_P), | intent(inout) | :: | r_primitive(1:self%Np,1:2,0:self%Ni+1) | Reconstructed primitive variables. |
Solve the Riemann Problem at cell interfaces.
Solve the Riemann problem between the state $1$ and $4$ using the (local) Lax Friedrichs (Rusanov) solver.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
||
real(kind=R_P), | intent(in) | :: | p1 | Pressure of state 1. |
||
real(kind=R_P), | intent(in) | :: | r1 | Density of state 1. |
||
real(kind=R_P), | intent(in) | :: | u1 | Velocity of state 1. |
||
real(kind=R_P), | intent(in) | :: | g1 | Specific heats ratio of state 1. |
||
real(kind=R_P), | intent(in) | :: | p4 | Pressure of state 4. |
||
real(kind=R_P), | intent(in) | :: | r4 | Density of state 4. |
||
real(kind=R_P), | intent(in) | :: | u4 | Velocity of state 4. |
||
real(kind=R_P), | intent(in) | :: | g4 | Specific heats ratio of state 4. |
||
real(kind=R_P), | intent(out) | :: | F(1:self%Nc) | Resulting fluxes. |
Synchronize CAF images.
Synchronize CAF images.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
Time derivative, residuals function.
Time derivative of Euler field, the residuals function.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(euler_1D_caf_nf), | intent(in) | :: | self | Euler field. |
Euler field time derivative.
type :: euler_1D_caf_nf
!< Euler 1D (CAF enabled) PDEs system field.
!<
!< @note This object does not use FOODIE.
!<
!<### 1D Euler PDEs system
!< The 1D Euler PDEs system considered is a non linear, hyperbolic (inviscid) system of conservation laws for compressible gas
!< dynamics, that reads as
!<$$
!<\begin{matrix}
!<U_t = R(U) \Leftrightarrow U_t = F(U)_x \\
!<U = \begin{bmatrix}
!<\rho \\
!<\rho u \\
!<\rho E
!<\end{bmatrix}\;\;\;
!<F(U) = \begin{bmatrix}
!<\rho u \\
!<\rho u^2 + p \\
!<\rho u H
!<\end{bmatrix}
!<\end{matrix}
!<$$
!< where \(\rho\) is the density, \(u\) is the velocity, \(p\) the pressure, \(E\) the total internal specific energy and \(H\)
!< the total specific enthalpy. The PDEs system must completed with the proper initial and boundary conditions. Moreover, an ideal
!< (thermally and calorically perfect) gas is considered
!<$$
!<\begin{matrix}
!<R = c_p - c_v \\
!<\gamma = \frac{c_p}{c_v}\\
!<e = c_v T \\
!<h = c_p T
!<\end{matrix}
!<$$
!< where *R* is the gas constant, \(c_p\,c_v\) are the specific heats at constant pressure and volume (respectively), *e* is the
!< internal energy, *h* is the internal enthalpy and *T* is the temperature. The following addition equations of state hold:
!<$$
!<\begin{matrix}
!<T = \frac{p}{\rho R} \\
!<E = \rho e + \frac{1}{2} \rho u^2 \\
!<H = \rho h + \frac{1}{2} \rho u^2 \\
!<a = \sqrt{\frac{\gamma p}{\rho}}
!<\end{matrix}
!<$$
!<
!<### Multi-fluid Euler PDEs system
!< An extension of the above Euler system is considered allowing the modelling of a multi-fluid mixture of different gas (with
!< different physical characteristics). The well known Standard Thermodynamic Model is used to model the gas mixture replacing the
!< density with the density fraction of each specie composing the mixture. This led to the following system:
!<$$
!<\begin{matrix}
!<U_t = R(U) \Leftrightarrow U_t = F(U)_x \\
!<U = \begin{bmatrix}
!<\rho_s \\
!<\rho u \\
!<\rho E
!<\end{bmatrix}\;\;\;
!<F(U) = \begin{bmatrix}
!<\rho_s u \\
!<\rho u^2 + p \\
!<\rho u H
!<\end{bmatrix}\;\;\; for\; s=1,2,...N_s \\
!<\rho = \sum_{s=1}^{N_s}\rho_s \\
!<c_p = \sum_{s=1}^{N_S} \frac{\rho_s}{\rho} c_{p,s} \quad c_v = \sum_{s=1}^{N_S} \frac{\rho_s}{\rho} c_{v,s}
!<\end{matrix}
!<$$
!< where \(N_s\) is the number of initial species composing the gas mixture.
!<
!<#### Numerical grid organization
!< The finite volume, Godunov's like approach is employed. The conservative variables (and the primitive ones) are co-located at
!< the cell center. The cell and (inter)faces numeration is as follow.
!<```
!< cell (inter)faces
!< | |
!< v v
!< |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
!< | 1-Ng | 2-Ng | ..... | -1 | 0 | 1 | 2 | ..... | Ni | Ni+1 | Ni+1 | ..... |Ni+Ng-1| Ni+Ng |
!< |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
!< 0-Ng -1 0 1 2 Ni-1 Ni Ni+Ng
!<```
!< Where *Ni* are the finite volumes (cells) used for discretizing the domain and *Ng* are the ghost cells used for imposing the
!< left and right boundary conditions (for a total of *2Ng* cells).
!<
!<#### Primitive variables organization
!< Primitive variables are organized as an array of reals which the first index means:
!<
!< + 1 : density of species 1 (r1)
!< + 2 : density of species 2 (r2)
!< + ... :
!< + s : density of species s-th (rs)
!< + ... :
!< + Ns : density of species Ns (rNs)
!< + Ns+1 : velocity (u)
!< + Ns+2 : pressure (p)
!< + Ns+3 : density (r=sum(rs))
!< + Ns+4 : specific heats ratio (g)
!<
!<#### Conservative variables organization
!< Conservative variables are organized as an array (rank 2) of reals which the first index means:
!<
!< + 1 : mass conservation of species 1 (r1)
!< + 2 : mass conservation of species 2 (r2)
!< + ... :
!< + s : mass conservation of species s-th (rs)
!< + ... :
!< + Ns : mass conservation of species Ns (rNs)
!< + Ns+1 : momentum conservation (r*u)
!< + Ns+2 : energy conservation (r*E)
private
integer(I_P) :: ord=0 !< Space accuracy formal order.
integer(I_P) :: Ni=0 !< Space dimension.
integer(I_P) :: Ng=0 !< Number of ghost cells for boundary conditions handling.
integer(I_P) :: Ns=0 !< Number of initial species.
integer(I_P) :: Nc=0 !< Number of conservative variables, Ns+2.
integer(I_P) :: Np=0 !< Number of primitive variables, Ns+4.
real(R_P) :: Dx=0._R_P !< Space step.
type(weno_interpolator_upwind) :: weno !< WENO interpolator.
real(R_P), allocatable :: U(:,:) !< Integrand (state) variables, whole physical domain [1:Nc,1:Ni].
real(R_P), allocatable :: cp0(:) !< Specific heat cp of initial species [1:Ns].
real(R_P), allocatable :: cv0(:) !< Specific heat cv of initial species [1:Ns].
character(:), allocatable :: BC_L !< Left boundary condition type.
character(:), allocatable :: BC_R !< Right boundary condition type.
integer(I_P) :: me=0 !< ID of this_image().
integer(I_P) :: we=0 !< Number of CAF images used.
contains
! auxiliary methods
procedure, pass(self), public :: init !< Init field.
procedure, pass(self), public :: destroy !< Destroy field.
procedure, pass(self), public :: output !< Extract Euler field.
procedure, pass(self), public :: dt => compute_dt !< Compute the current time step, by means of CFL condition.
generic, public :: assignment(=) => assign_integrand !< Overloading = assignament.
! private methods
procedure, pass(self), private :: t => dEuler_dt !< Time derivative, residuals function.
procedure, pass(lhs), private :: assign_integrand => euler_assign_euler !< Euler = Euler.
procedure, pass(self), private :: primitive2conservative !< Convert primitive variables to conservative ones.
procedure, pass(self), private :: conservative2primitive !< Convert conservative variables to primitive ones.
procedure, pass(self), private :: synchronize !< Synchronize CAF images.
procedure, pass(self), private :: impose_boundary_conditions !< Impose boundary conditions.
procedure, pass(self), private :: reconstruct_interfaces_states !< Reconstruct interfaces states.
procedure, pass(self), private :: riemann_solver !< Solve the Riemann Problem at cell interfaces.
endtype euler_1D_caf_nf