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. |
subroutine initialize(self, Ni, Dx, BC_L, BC_R, initial_state, eos, weno_order, weno_variables, riemann_solver_scheme)
!< Initialize field.
class(euler_1d), intent(inout) :: self !< Euler field.
integer(I4P), intent(in) :: Ni !< Space dimension.
real(R8P), intent(in) :: Dx !< Space step.
character(*), intent(in) :: BC_L !< Left boundary condition type.
character(*), 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(I4P), intent(in), optional :: weno_order !< WENO reconstruction order.
character(*), intent(in), optional :: weno_variables !< Variables on which WENO reconstruction is done.
character(*), intent(in), optional :: riemann_solver_scheme !< Riemann solver scheme.
character(:), allocatable :: weno_variables_ !< WENO Variables, local variable.
character(:), allocatable :: riemann_solver_scheme_ !< Riemann solver scheme, local variable.
integer(I4P) :: i !< Space couner.
call self%destroy
self%weno_order = 1 ; if (present(weno_order)) self%weno_order = weno_order
self%Ni = Ni
self%Ng = (self%weno_order + 1) / 2
self%Dx = Dx
self%eos = eos
if (allocated(self%U)) deallocate(self%U) ; allocate(self%U(1-self%Ng:self%Ni+self%Ng))
do i=1, Ni
self%U(i) = primitive_to_conservative_compressible(primitive=initial_state(i), eos=eos)
enddo
self%BC_L = BC_L
self%BC_R = BC_R
if (self%weno_order>1) call wenoof_create(interpolator_type='reconstructor-JS', S=self%Ng, interpolator=self%interpolator)
weno_variables_ = 'characteristic'
if (present(weno_variables)) weno_variables_ = trim(adjustl(weno_variables))
select case(weno_variables_)
case('characteristic')
self%reconstruct_interfaces => reconstruct_interfaces_characteristic
case('conservative')
self%reconstruct_interfaces => reconstruct_interfaces_conservative
case('primitive')
self%reconstruct_interfaces => reconstruct_interfaces_primitive
case default
write(stderr, '(A)') 'error: WENO reconstruction variables set "'//weno_variables_//'" unknown!'
stop
endselect
riemann_solver_scheme_ = 'llf'
if (present(riemann_solver_scheme)) riemann_solver_scheme_ = trim(adjustl(riemann_solver_scheme))
select case(riemann_solver_scheme_)
case('exact')
allocate(riemann_solver_compressible_exact :: self%riemann_solver)
case('hllc')
allocate(riemann_solver_compressible_hllc :: self%riemann_solver)
case('llf')
allocate(riemann_solver_compressible_llf :: self%riemann_solver)
case('pvl')
allocate(riemann_solver_compressible_pvl :: self%riemann_solver)
case('roe')
allocate(riemann_solver_compressible_roe :: self%riemann_solver)
case default
write(stderr, '(A)') 'error: Riemann Solver scheme "'//riemann_solver_scheme_//'" unknown!'
stop
endselect
endsubroutine initialize