foreseer_test_riemann_solver_compressible_llf.f90 Source File

FORESEER test: Riemann solver compressible LLF class test.

This File Depends On

sourcefile~~foreseer_test_riemann_solver_compressible_llf.f90~~EfferentGraph sourcefile~foreseer_test_riemann_solver_compressible_llf.f90 foreseer_test_riemann_solver_compressible_llf.f90 sourcefile~foreseer.f90 foreseer.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_llf.f90 sourcefile~foreseer_compressible_transformations.f90 foreseer_compressible_transformations.f90 sourcefile~foreseer_compressible_transformations.f90->sourcefile~foreseer.f90 sourcefile~foreseer_eos_compressible.f90 foreseer_eos_compressible.f90 sourcefile~foreseer_eos_compressible.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_roe.f90 foreseer_riemann_solver_compressible_roe.F90 sourcefile~foreseer_eos_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_riemann_pattern_compressible_object.f90 foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_eos_compressible.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_riemann_solver_compressible_roe.f90->sourcefile~foreseer.f90 sourcefile~foreseer_primitive_object.f90 foreseer_primitive_object.f90 sourcefile~foreseer_primitive_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_primitive_compressible.f90 foreseer_primitive_compressible.f90 sourcefile~foreseer_primitive_object.f90->sourcefile~foreseer_primitive_compressible.f90 sourcefile~foreseer_riemann_pattern_compressible_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90 foreseer_riemann_pattern_compressible_pvl.f90 sourcefile~foreseer_riemann_pattern_compressible_object.f90->sourcefile~foreseer_riemann_pattern_compressible_pvl.f90 sourcefile~foreseer_conservative_compressible.f90 foreseer_conservative_compressible.f90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer.f90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_compressible_transformations.f90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_riemann_solver_compressible_pvl.f90 foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_riemann_solver_compressible_hllc.f90 foreseer_riemann_solver_compressible_hllc.F90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_riemann_solver_compressible_exact.f90 foreseer_riemann_solver_compressible_exact.F90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_riemann_solver_compressible_llf.f90 foreseer_riemann_solver_compressible_llf.F90 sourcefile~foreseer_conservative_compressible.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90 sourcefile~foreseer_riemann_solver_compressible_pvl.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_hllc.f90->sourcefile~foreseer.f90 sourcefile~foreseer_primitive_compressible.f90->sourcefile~foreseer.f90 sourcefile~foreseer_primitive_compressible.f90->sourcefile~foreseer_compressible_transformations.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90 sourcefile~foreseer_riemann_solver_object.f90 foreseer_riemann_solver_object.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_riemann_solver_object.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90 sourcefile~foreseer_riemann_pattern_object.f90 foreseer_riemann_pattern_object.f90 sourcefile~foreseer_riemann_pattern_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_pattern_object.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_conservative_object.f90 foreseer_conservative_object.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_conservative_compressible.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_pattern_compressible_pvl.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_object.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_pattern_object.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90 sourcefile~foreseer_riemann_solver_compressible_exact.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_llf.f90->sourcefile~foreseer.f90 sourcefile~foreseer_eos_object.f90 foreseer_eos_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_compressible_transformations.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_eos_compressible.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_primitive_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_conservative_compressible.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_primitive_compressible.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_pattern_compressible_pvl.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_pattern_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_conservative_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90
Help


Source Code

!< FORESEER test: Riemann solver compressible LLF class test.

program foreseer_test_riemann_solver_compressible_llf
!< FORESEER test: Riemann solver compressible LLF class test.

use foreseer, only : eos_compressible, conservative_compressible, riemann_solver_compressible_llf
use penf, only : R8P, str
use vecfor, only : ex

implicit none
type(eos_compressible)                :: eos                          !< The equation of state.
type(conservative_compressible)       :: state_left                   !< Left state.
type(conservative_compressible)       :: state_right                  !< Right state.
type(conservative_compressible)       :: fluxes                       !< Conservative fluxes.
type(riemann_solver_compressible_llf) :: riemann_solver               !< Riemann solver.
real(R8P)                             :: waves(1:5)                   !< Waves pattern.
real(R8P), parameter                  :: r_2=0.426319003105163574_R8P !< Exact value of density in state 2.
real(R8P), parameter                  :: r_3=0.265574008226394653_R8P !< Exact value of density in state 3.
real(R8P), parameter                  :: p23=0.303130000829696655_R8P !< Exact value of pressure in states 2 and 3.
real(R8P), parameter                  :: u23=0.927452981472015381_R8P !< Exact value of velocity in states 2 and 3.
logical                               :: are_tests_passed(1)          !< List of passed tests.

are_tests_passed = .true.

eos = eos_compressible(cp=1040.004_R8P, cv=742.86_R8P)
state_left  = conservative_compressible(density=1._R8P,    energy=1._R8P   *eos%energy(density=1._R8P,    pressure=1._R8P) )
state_right = conservative_compressible(density=0.125_R8P, energy=0.125_R8P*eos%energy(density=0.125_R8P, pressure=0.1_R8P))

print '(A)', 'Test solution with "u23" algorithm:'
call riemann_solver%initialize(config='u23')
call riemann_solver%solve(eos_left=eos, state_left=state_left, eos_right=eos, state_right=state_right, normal=ex, fluxes=fluxes)
print '(A)', 'Fluxes at interface:'
print '(A)', fluxes%description(prefix='  ')
call fluxes%compute_fluxes_from_primitive(eos=eos, p=p23, r=r_2, u=u23, normal=ex)
print '(A)', 'Exact fluxes at interface:'
print '(A)', fluxes%description(prefix='  ')
print '(A)', 'Exact intemediate states:'
print '(A)', '  r_2 = '//str(n=r_2)
print '(A)', '  r_3 = '//str(n=r_3)
print '(A)', '  p23 = '//str(n=p23)
print '(A)', '  u23 = '//str(n=u23)
print '(A)', 'Test solution with "up23" algorithm:'
call riemann_solver%initialize(config='up23')
call riemann_solver%solve(eos_left=eos, state_left=state_left, eos_right=eos, state_right=state_right, normal=ex, fluxes=fluxes)
print '(A)', 'Fluxes at interface:'
print '(A)', fluxes%description(prefix='  ')
call fluxes%compute_fluxes_from_primitive(eos=eos, p=p23, r=r_2, u=u23, normal=ex)
print '(A)', 'Exact fluxes at interface:'
print '(A)', fluxes%description(prefix='  ')
print '(A)', 'Exact intemediate states:'
print '(A)', '  r_2 = '//str(n=r_2)
print '(A)', '  r_3 = '//str(n=r_3)
print '(A)', '  p23 = '//str(n=p23)
print '(A)', '  u23 = '//str(n=u23)

print "(A,L1)", new_line('a')//'Are all tests passed? ', all(are_tests_passed)
endprogram foreseer_test_riemann_solver_compressible_llf