Euler 1D PDEs system field.
It is a FOODIE integrand class concrete extension.
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} $$
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).
| Type | Visibility | Attributes | Name | Initial | |||
|---|---|---|---|---|---|---|---|
| character(len=:), | public, | allocatable | :: | BC_L | Left boundary condition type. |
||
| character(len=:), | public, | allocatable | :: | BC_R | Right boundary condition type. |
||
| real(kind=R8P), | public | :: | Dx | = | 0._R8P | Space step. |
|
| integer(kind=I4P), | public | :: | Ng | = | 0 | Ghost cells number. |
|
| integer(kind=I4P), | public | :: | Ni | = | 0 | Space dimension. |
|
| type(conservative_compressible), | public, | allocatable | :: | U(:) | Integrand (state) variables. |
||
| type(eos_compressible), | public | :: | eos | Equation of state. |
|||
| class(interpolator_object), | public, | allocatable | :: | interpolator | WENO interpolator. |
||
| procedure(reconstruct_interfaces_), | public, | pointer | :: | reconstruct_interfaces | => | reconstruct_interfaces_characteristic | Reconstruct interface states. |
| class(riemann_solver_object), | public, | allocatable | :: | riemann_solver | Riemann solver. |
||
| integer(kind=I4P), | public | :: | weno_order | = | 0 | WENO reconstruction order. |
Operator +.
Operator =.
Assign one Euler field to another.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(inout) | :: | lhs | Left hand side. |
||
| class(integrand), | intent(in) | :: | rhs | Right hand side. |
Operator euler = real.
Assign one real to an Euler field.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(inout) | :: | lhs | Left hand side. |
||
| real(kind=R8P), | intent(in) | :: | rhs | Right hand side. |
Destroy 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), | intent(in) | :: | self | Euler field. |
||
| integer(kind=I4P), | intent(in) | :: | steps_max | Maximun number of time steps. |
||
| real(kind=R8P), | intent(in) | :: | t_max | Maximum integration time. |
||
| real(kind=R8P), | intent(in) | :: | t | Time. |
||
| real(kind=R8P), | intent(in) | :: | CFL | CFL value. |
Time step.
Impose boundary conditions.
Impose boundary conditions.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| type(conservative_compressible), | intent(inout) | :: | U(1-self%Ng:) | Conservative variables. |
Initialize field.
Initialize field.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(inout) | :: | self | Euler field. |
||
| integer(kind=I4P), | intent(in) | :: | Ni | Space dimension. |
||
| real(kind=R8P), | 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. |
||
| type(primitive_compressible), | intent(in) | :: | initial_state(1:) | Initial state of primitive variables. |
||
| type(eos_compressible), | intent(in) | :: | eos | Equation of state. |
||
| integer(kind=I4P), | intent(in), | optional | :: | weno_order | WENO reconstruction order. |
|
| character(len=*), | intent(in), | optional | :: | weno_variables | Variables on which WENO reconstruction is done. |
|
| character(len=*), | intent(in), | optional | :: | riemann_solver_scheme | Riemann solver scheme. |
Operator *.
Multiply an Euler field by another one.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | lhs | Left hand side. |
||
| class(integrand), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Operator euler * real.
Multiply an Euler field by a real scalar.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | lhs | Left hand side. |
||
| real(kind=R8P), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Operator ||euler-euler||.
Estimate local truncation error between 2 euler approximations.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | lhs | Left hand side. |
||
| class(integrand), | intent(in) | :: | rhs | Right hand side. |
Error estimation.
Extract Euler field.
Output the Euler field state.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| logical, | intent(in), | optional | :: | is_primitive | Output in primitive variables. |
Euler state vector.
Operator real * euler.
Multiply a real scalar by an Euler field.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| real(kind=R8P), | intent(in) | :: | lhs | Left hand side. |
||
| class(euler_1d), | intent(in) | :: | rhs | Right hand side. |
Operator result.
Reconstruct (charc.) interface states.
Reconstruct interfaces states.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| type(conservative_compressible), | intent(in) | :: | conservative(1-self%Ng:) | Conservative variables. |
||
| type(conservative_compressible), | intent(inout) | :: | r_conservative(1:,0:) | Reconstructed conservative vars. |
Reconstruct (cons.) interface states.
Reconstruct interfaces states.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| type(conservative_compressible), | intent(in) | :: | conservative(1-self%Ng:) | Conservative variables. |
||
| type(conservative_compressible), | intent(inout) | :: | r_conservative(1:,0:) | Reconstructed conservative vars. |
Reconstruct (prim.) interface states.
Reconstruct interfaces states.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| type(conservative_compressible), | intent(in) | :: | conservative(1-self%Ng:) | Conservative variables. |
||
| type(conservative_compressible), | intent(inout) | :: | r_conservative(1:,0:) | Reconstructed conservative vars. |
Operator -.
Time derivative, residuals function.
Time derivative of Euler field, the residuals function.
| Type | Intent | Optional | Attributes | Name | ||
|---|---|---|---|---|---|---|
| class(euler_1d), | intent(in) | :: | self | Euler field. |
||
| real(kind=R8P), | intent(in), | optional | :: | t | Time. |
Euler field time derivative.
type, extends(integrand) :: euler_1d
!< Euler 1D PDEs system field.
!<
!< It is a FOODIE integrand class concrete extension.
!<
!<### 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}
!<$$
!<
!<#### 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).
integer(I4P) :: weno_order=0 !< WENO reconstruction order.
integer(I4P) :: Ni=0 !< Space dimension.
integer(I4P) :: Ng=0 !< Ghost cells number.
real(R8P) :: Dx=0._R8P !< Space step.
type(eos_compressible) :: eos !< Equation of state.
type(conservative_compressible), allocatable :: U(:) !< Integrand (state) variables.
character(:), allocatable :: BC_L !< Left boundary condition type.
character(:), allocatable :: BC_R !< Right boundary condition type.
class(interpolator_object), allocatable :: interpolator !< WENO interpolator.
procedure(reconstruct_interfaces_), pointer :: reconstruct_interfaces=>&
reconstruct_interfaces_characteristic !< Reconstruct interface states.
class(riemann_solver_object), allocatable :: riemann_solver !< Riemann solver.
contains
! auxiliary methods
procedure, pass(self) :: initialize !< Initialize field.
procedure, pass(self) :: destroy !< Destroy field.
procedure, pass(self) :: output !< Extract Euler field.
procedure, pass(self) :: dt => compute_dt !< Compute the current time step, by means of CFL condition.
! ADT integrand deferred methods
procedure, pass(self) :: t => dEuler_dt !< Time derivative, residuals function.
procedure, pass(lhs) :: local_error => euler_local_error !< Operator `||euler-euler||`.
procedure, pass(lhs) :: integrand_multiply_integrand => euler_multiply_euler !< Operator `*`.
procedure, pass(lhs) :: integrand_multiply_real => euler_multiply_real !< Operator `euler * real`.
procedure, pass(rhs) :: real_multiply_integrand => real_multiply_euler !< Operator `real * euler`.
procedure, pass(lhs) :: add => add_euler !< Operator `+`.
procedure, pass(lhs) :: sub => sub_euler !< Operator `-`.
procedure, pass(lhs) :: assign_integrand => euler_assign_euler !< Operator `=`.
procedure, pass(lhs) :: assign_real => euler_assign_real !< Operator `euler = real`.
! private methods
procedure, pass(self), private :: impose_boundary_conditions !< Impose boundary conditions.
procedure, pass(self), private :: reconstruct_interfaces_characteristic !< Reconstruct (charc.) interface states.
procedure, pass(self), private :: reconstruct_interfaces_conservative !< Reconstruct (cons.) interface states.
procedure, pass(self), private :: reconstruct_interfaces_primitive !< Reconstruct (prim.) interface states.
endtype euler_1d