foreseer_riemann_pattern_object.f90 Source File

Define the abstract Riemann (states) pattern for FORESEER library.

This File Depends On

sourcefile~~foreseer_riemann_pattern_object.f90~~EfferentGraph sourcefile~foreseer_riemann_pattern_object.f90 foreseer_riemann_pattern_object.f90 sourcefile~foreseer_conservative_object.f90 foreseer_conservative_object.f90 sourcefile~foreseer_conservative_object.f90->sourcefile~foreseer_riemann_pattern_object.f90 sourcefile~foreseer_eos_object.f90 foreseer_eos_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_riemann_pattern_object.f90 sourcefile~foreseer_eos_object.f90->sourcefile~foreseer_conservative_object.f90
Help

Files Dependent On This One

sourcefile~~foreseer_riemann_pattern_object.f90~~AfferentGraph sourcefile~foreseer_riemann_pattern_object.f90 foreseer_riemann_pattern_object.f90 sourcefile~foreseer_riemann_pattern_compressible_object.f90 foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer_riemann_pattern_object.f90->sourcefile~foreseer_riemann_pattern_compressible_object.f90 sourcefile~foreseer.f90 foreseer.f90 sourcefile~foreseer_riemann_pattern_object.f90->sourcefile~foreseer.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_test_eos_compressible.f90 foreseer_test_eos_compressible.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_eos_compressible.f90 sourcefile~foreseer_test_primitive_compressible.f90 foreseer_test_primitive_compressible.F90 sourcefile~foreseer.f90->sourcefile~foreseer_test_primitive_compressible.f90 sourcefile~foreseer_test_riemann_solver_compressible_exact.f90 foreseer_test_riemann_solver_compressible_exact.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_exact.f90 sourcefile~foreseer_test_riemann_solver_compressible_hllc.f90 foreseer_test_riemann_solver_compressible_hllc.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_test_shock_tube.f90 foreseer_test_shock_tube.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_shock_tube.f90 sourcefile~foreseer_test_riemann_solver_compressible_llf.f90 foreseer_test_riemann_solver_compressible_llf.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_llf.f90 sourcefile~foreseer_test_riemann_solver_compressible_pvl.f90 foreseer_test_riemann_solver_compressible_pvl.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_test_riemann_solver_compressible_roe.f90 foreseer_test_riemann_solver_compressible_roe.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_riemann_solver_compressible_roe.f90 sourcefile~foreseer_test_compressible_transformations.f90 foreseer_test_compressible_transformations.f90 sourcefile~foreseer.f90->sourcefile~foreseer_test_compressible_transformations.f90 sourcefile~foreseer_test_conservative_compressible.f90 foreseer_test_conservative_compressible.F90 sourcefile~foreseer.f90->sourcefile~foreseer_test_conservative_compressible.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_llf.f90 foreseer_riemann_solver_compressible_llf.F90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_llf.f90 sourcefile~foreseer_riemann_solver_compressible_hllc.f90 foreseer_riemann_solver_compressible_hllc.F90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_hllc.f90 sourcefile~foreseer_riemann_solver_compressible_exact.f90 foreseer_riemann_solver_compressible_exact.F90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_exact.f90 sourcefile~foreseer_riemann_solver_compressible_pvl.f90 foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_pvl.f90 sourcefile~foreseer_riemann_solver_compressible_roe.f90 foreseer_riemann_solver_compressible_roe.F90 sourcefile~foreseer_riemann_pattern_compressible_pvl.f90->sourcefile~foreseer_riemann_solver_compressible_roe.f90 sourcefile~foreseer_riemann_solver_compressible_llf.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_hllc.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_exact.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_pvl.f90->sourcefile~foreseer.f90 sourcefile~foreseer_riemann_solver_compressible_roe.f90->sourcefile~foreseer.f90
Help


Source Code

!< Define the abstract Riemann (states) pattern for FORESEER library.

module foreseer_riemann_pattern_object
!< Define the abstract Riemann (states) pattern for FORESEER library.

use foreseer_conservative_object, only : conservative_object
use foreseer_eos_object, only : eos_object
use vecfor, only : vector

implicit none
private
public :: riemann_pattern_object

type, abstract :: riemann_pattern_object
   !< Riemann (states) pattern object class.
   !<
   !< This pattern is generated after the breaking of the initial discontinuity of the Riemann Problem.
   contains
      ! deferred methods
      procedure(compute_fluxes_interface), pass(self), deferred :: compute_fluxes   !< Compute fluxes at interface.
      procedure(compute_waves_interface),  pass(self), deferred :: compute_waves    !< Compute waves speed.
      procedure(description_interface),    pass(self), deferred :: description      !< Return pretty-printed object description.
      procedure(initialize_interface),     pass(self), deferred :: initialize       !< Initialize pattern with left/right states.
      procedure(assignment_interface),     pass(lhs),  deferred :: rpat_assign_rpat !< Operator `=`.
      ! operators
      generic :: assignment(=) => rpat_assign_rpat !< Overload `=`.
endtype riemann_pattern_object

abstract interface
   !< Abstract interfaces of deferred methods of [[riemann_pattern_object]].
   pure subroutine assignment_interface(lhs, rhs)
   !< Operator `=`.
   import :: riemann_pattern_object
   class(riemann_pattern_object), intent(inout) :: lhs !< Left hand side.
   class(riemann_pattern_object), intent(in)    :: rhs !< Right hand side.
   endsubroutine assignment_interface

   elemental subroutine compute_fluxes_interface(self, normal, fluxes)
   !< Compute fluxes at initial discontinuity interface.
   import :: conservative_object, riemann_pattern_object, vector
   class(riemann_pattern_object), intent(in)    :: self   !< Riemann (states) pattern solution.
   type(vector),                  intent(in)    :: normal !< Normal (versor) of face where fluxes are given.
   class(conservative_object),    intent(inout) :: fluxes !< Fluxes at initial discontinuity interface.
   endsubroutine compute_fluxes_interface

   pure subroutine compute_waves_interface(self)
   !< Compute fluxes at initial discontinuity interface.
   import :: riemann_pattern_object
   class(riemann_pattern_object), intent(inout) :: self !< Riemann (states) pattern solution.
   endsubroutine compute_waves_interface

   pure function description_interface(self, prefix) result(desc)
   !< Return a pretty-formatted object description.
   import :: riemann_pattern_object
   class(riemann_pattern_object), intent(in)           :: self   !< Riemann pattern.
   character(*),                  intent(in), optional :: prefix !< Prefixing string.
   character(len=:), allocatable                       :: desc   !< Description.
   endfunction description_interface

   elemental subroutine initialize_interface(self, eos_left, state_left, eos_right, state_right, normal)
   !< Initialize pattern with left and right states.
   import :: conservative_object, eos_object, riemann_pattern_object, vector
   class(riemann_pattern_object), intent(inout) :: self        !< Riemann (states) pattern solution.
   class(eos_object),             intent(in)    :: eos_left    !< Equation of state for left state.
   class(conservative_object),    intent(in)    :: state_left  !< Left Riemann state.
   class(eos_object),             intent(in)    :: eos_right   !< Equation of state for right state.
   class(conservative_object),    intent(in)    :: state_right !< Right Riemann state.
   type(vector),                  intent(in)    :: normal      !< Normal (versor) of face where fluxes are given.
   endsubroutine initialize_interface
endinterface
endmodule foreseer_riemann_pattern_object