compute_face_estimates Subroutine

private subroutine compute_face_estimates(qp, f_qp_left, f_qp_right, flags, dims)

Subroutine to calculate state at the face, generalized for

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var):: qp

Store primitive variable at cell center

real(kind=wp), intent(inout), dimension(1-flags(1):dims%imx-1+2*flags(1), 1-flags(2):dims%jmx-1+2*flags(2), 1-flags(3):dims%kmx-1+2*flags(3),1:dims%n_var):: f_qp_left

primitive state at faces

real(kind=wp), intent(inout), dimension(1-flags(1):dims%imx-1+2*flags(1), 1-flags(2):dims%jmx-1+2*flags(2), 1-flags(3):dims%kmx-1+2*flags(3),1:dims%n_var):: f_qp_right

primitive state at faces

integer, intent(in), dimension(3):: flags

Flags for direction switch

type(extent), intent(in) :: dims

Extent of the domain:imx,jmx,kmx


Calls

proc~~compute_face_estimates~~CallsGraph proc~compute_face_estimates compute_face_estimates debugcall debugcall proc~compute_face_estimates->debugcall

Called by

proc~~compute_face_estimates~~CalledByGraph proc~compute_face_estimates compute_face_estimates proc~compute_ppm_states compute_ppm_states proc~compute_ppm_states->proc~compute_face_estimates proc~compute_face_interpolant compute_face_interpolant proc~compute_face_interpolant->proc~compute_ppm_states proc~get_total_conservative_residue get_total_conservative_Residue proc~get_total_conservative_residue->proc~compute_face_interpolant proc~get_next_solution get_next_solution proc~get_next_solution->proc~get_total_conservative_residue proc~iterate_one_more_time_step iterate_one_more_time_step proc~iterate_one_more_time_step->proc~get_next_solution program~main main program~main->proc~iterate_one_more_time_step

Contents


Source Code

        subroutine compute_face_estimates(qp, f_qp_left, f_qp_right, flags, dims)
          !< Subroutine to calculate state at the face, generalized for

            implicit none
            type(extent), intent(in) :: dims
            !< Extent of the domain:imx,jmx,kmx
            integer, dimension(3), intent(in) :: flags
            !< Flags for direction switch
            real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), intent(in):: qp
            !< Store primitive variable at cell center
            real(wp), dimension(1-flags(1):dims%imx-1+2*flags(1), 1-flags(2):dims%jmx-1+2*flags(2),&
            1-flags(3):dims%kmx-1+2*flags(3),1:dims%n_var), intent(inout) :: f_qp_left, f_qp_right
            !< primitive state at faces
            integer :: i, j, k 
            integer :: i_f, j_f, k_f ! Flags to determine face direction

            DebugCall('compute_face_estimates')
            
            i_f = flags(1)
            j_f = flags(2)
            k_f = flags(3)

            ! Interior faces
            do k = (1 - k_f), dims%kmx - 1 + 2*k_f
             do j = (1 - j_f), dims%jmx - 1 + 2*j_f
              do i = (1 - i_f), dims%imx - 1 + 2*i_f
                f_qp_left(i, j, k, :) = (7. * (qp(i, j, k, :) + &
                    qp(i - i_f, j - j_f, k - k_f, :)) - (qp(i + i_f, j + j_f, k + k_f, :) + &
                    qp(i - 2*i_f, j - 2*j_f, k - 2*k_f, :))) / 12.
              end do
             end do
            end do
            f_qp_right= f_qp_left

        end subroutine compute_face_estimates