reconstruct_interfaces_characteristic Subroutine

private subroutine reconstruct_interfaces_characteristic(self, conservative, r_conservative)

Reconstruct interfaces states.

The reconstruction is done in pseudo characteristic variables.

Arguments

Type IntentOptional AttributesName
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.

Calls

proc~~reconstruct_interfaces_characteristic~~CallsGraph proc~reconstruct_interfaces_characteristic reconstruct_interfaces_characteristic proc~primitive_to_conservative_compressible primitive_to_conservative_compressible proc~reconstruct_interfaces_characteristic->proc~primitive_to_conservative_compressible wenoof_create wenoof_create proc~reconstruct_interfaces_characteristic->wenoof_create proc~conservative_to_primitive_compressible conservative_to_primitive_compressible proc~reconstruct_interfaces_characteristic->proc~conservative_to_primitive_compressible
Help


Source Code

   subroutine reconstruct_interfaces_characteristic(self, conservative, r_conservative)
   !< Reconstruct interfaces states.
   !<
   !< The reconstruction is done in pseudo characteristic variables.
   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.
   type(primitive_compressible)                   :: primitive(1-self%Ng:self%Ni+self%Ng) !< Primitive variables.
   type(primitive_compressible)                   :: r_primitive(1:2, 0:self%Ni+1)        !< Reconstructed primitive variables.
   type(primitive_compressible)                   :: Pm(1:2)                              !< Mean of primitive variables.
   real(R8P)                                      :: LPm(1:3, 1:3, 1:2)                   !< Mean left eigenvectors matrix.
   real(R8P)                                      :: RPm(1:3, 1:3, 1:2)                   !< Mean right eigenvectors matrix.
   real(R8P)                                      :: C(1:2, 1-self%Ng:-1+self%Ng, 1:3)    !< Pseudo characteristic variables.
   real(R8P)                                      :: CR(1:2, 1:3)                         !< Pseudo characteristic reconst.
   real(R8P)                                      :: buffer(1:3)                          !< Dummy buffer.
   integer(I4P)                                   :: i                                    !< Counter.
   integer(I4P)                                   :: j                                    !< Counter.
   integer(I4P)                                   :: f                                    !< Counter.
   integer(I4P)                                   :: v                                    !< Counter.
   class(interpolator_object), allocatable        :: interpolator                         !< WENO interpolator.

   select case(self%weno_order)
   case(1) ! 1st order piecewise constant reconstruction
      do i=0, self%Ni+1
         r_conservative(1, i) = conservative(i)
         r_conservative(2, i) = r_conservative(1, i)
      enddo
   case(3, 5, 7, 9, 11, 13, 15, 17) ! 3rd-17th order WENO reconstruction
      call wenoof_create(interpolator_type='reconstructor-JS', S=self%Ng, interpolator=interpolator)
      do i=1-self%Ng, self%Ni+self%Ng
         primitive(i) = conservative_to_primitive_compressible(conservative=conservative(i), eos=self%eos)
      enddo
      do i=0, self%Ni+1
         ! compute pseudo charteristic variables
         do f=1, 2
            Pm(f) = 0.5_R8P * (primitive(i+f-2) + primitive(i+f-1))
         enddo
         do f=1, 2
            LPm(:, :, f) = Pm(f)%left_eigenvectors(eos=self%eos)
            RPm(:, :, f) = Pm(f)%right_eigenvectors(eos=self%eos)
         enddo
         do j=i+1-self%Ng, i-1+self%Ng
            do f=1, 2
               do v=1, 3
                  C(f, j-i, v) = dot_product(LPm(v, :, f), [primitive(j)%density, primitive(j)%velocity%x, primitive(j)%pressure])
               enddo
            enddo
         enddo
         ! compute WENO reconstruction of pseudo charteristic variables
         do v=1, 3
            call interpolator%interpolate(stencil=C(:, :, v), interpolation=CR(:, v))
         enddo
         ! trasform back reconstructed pseudo charteristic variables to primitive ones
         do f=1, 2
            do v=1, 3
               buffer(v) = dot_product(RPm(v, :, f), CR(f, :))
            enddo
            r_primitive(f, i)%density  = buffer(1)
            r_primitive(f, i)%velocity = buffer(2) * ex
            r_primitive(f, i)%pressure = buffer(3)
         enddo
      enddo
      do i=0, self%Ni+1
         r_conservative(1, i) = primitive_to_conservative_compressible(primitive=r_primitive(1, i), eos=self%eos)
         r_conservative(2, i) = primitive_to_conservative_compressible(primitive=r_primitive(2, i), eos=self%eos)
      enddo
   endselect
   endsubroutine reconstruct_interfaces_characteristic


add add add_euler array array compute_derivate compute_dt compute_fluxes compute_fluxes compute_fluxes_from_primitive compute_post_rarefaction compute_post_shock compute_roe_state compute_states23_from_u23 compute_u23 compute_up23 compute_waves compute_waves_u23 compute_waves_up23 cons_assign_cons cons_divide_real cons_multiply_cons cons_multiply_real conservative_compressible conservative_compressible_instance conservative_compressible_pointer conservative_to_primitive_compressible cp cv delta density description description description description destroy destroy destroy dEuler_dt energy energy eos_assign_eos eos_compressible eos_compressible_instance eos_compressible_pointer eta euler_assign_euler euler_assign_real euler_local_error euler_multiply_euler euler_multiply_real g gm1 gp1 impose_boundary_conditions initialize initialize initialize initialize initialize initialize initialize initialize initialize initialize left_eigenvectors momentum negative negative output parse_command_line_interface positive positive pressure pressure prim_assign_prim prim_divide_real prim_multiply_prim prim_multiply_real primitive_compressible primitive_compressible_instance primitive_compressible_pointer primitive_to_conservative_compressible R real_multiply_cons real_multiply_euler real_multiply_prim reconstruct_interfaces_characteristic reconstruct_interfaces_conservative reconstruct_interfaces_primitive right_eigenvectors rpat_assign_rpat save_time_serie solve solve solve solve solve speed_of_sound sub sub sub_euler temperature velocity