compute_muscl_states Subroutine

public subroutine compute_muscl_states(qp, x_qp_l, x_qp_r, y_qp_l, y_qp_r, z_qp_l, z_qp_r, pdif, scheme, flow, dims)

Implement MUSCL scheme to get left and right states at each face. The computation is done through all cells and first level ghost cells

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(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var):: x_qp_l

Store primitive state at the I-face

real(kind=wp), intent(inout), dimension(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var):: x_qp_r

Store primitive state at the I-face

real(kind=wp), intent(inout), dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var):: y_qp_l

Store primitive state at the J-face

real(kind=wp), intent(inout), dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var):: y_qp_r

Store primitive state at the J-face

real(kind=wp), intent(inout), dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var):: z_qp_l

Store primitive state at the K-face

real(kind=wp), intent(inout), dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var):: z_qp_r

Store primitive state at the K-face

real(kind=wp), intent(inout), dimension(0:dims%imx,0:dims%jmx,0:dims%kmx):: pdif

pressure difference

type(schemetype), intent(in) :: scheme

finite-volume Schemes

type(flowtype), intent(in) :: flow

Information about fluid flow: freestream-speed, ref-viscosity,etc.

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

Extent of the domain:imx,jmx,kmx


Calls

proc~~compute_muscl_states~~CallsGraph proc~compute_muscl_states compute_muscl_states proc~compute_face_state compute_face_state proc~compute_muscl_states->proc~compute_face_state proc~pressure_based_switching pressure_based_switching proc~compute_muscl_states->proc~pressure_based_switching debugcall debugcall proc~compute_face_state->debugcall proc~pressure_based_switching->debugcall

Called by

proc~~compute_muscl_states~~CalledByGraph proc~compute_muscl_states compute_muscl_states proc~compute_face_interpolant compute_face_interpolant proc~compute_face_interpolant->proc~compute_muscl_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


Source Code

        subroutine compute_muscl_states(qp, x_qp_l, x_qp_r, y_qp_l, y_qp_r, z_qp_l, z_qp_r, pdif, scheme, flow, dims)
            !< Implement MUSCL scheme to get left and right states at
            !< each face. The computation is done through all cells
            !< and first level ghost cells
            !---------------------------------------------------------
            implicit none
            type(extent), intent(in) :: dims
            !< Extent of the domain:imx,jmx,kmx
            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(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var), intent(inout) :: x_qp_l, x_qp_r
            !< Store primitive state at the I-face 
            real(wp), dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var), intent(inout) :: y_qp_l, y_qp_r
            !< Store primitive state at the J-face 
            real(wp), dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var), intent(inout) :: z_qp_l, z_qp_r
            !< Store primitive state at the K-face 
            real(wp), dimension(0:dims%imx,0:dims%jmx,0:dims%kmx), intent(inout) :: pdif
            !< pressure difference
            type(schemetype), intent(in) :: scheme
            !< finite-volume Schemes
            type(flowtype), intent(in) ::flow
            !< Information about fluid flow: freestream-speed, ref-viscosity,etc.
            integer, dimension(3) :: flags
            !< Flags for direction

            
            flags=(/1,0,0/)
            call compute_face_state(qp, x_qp_l, x_qp_r, flags, scheme%ilimiter_switch, scheme%itlimiter_switch, dims)
            if(scheme%iPB_switch==1)then
              call pressure_based_switching(qp, x_qp_l, x_qp_r, pdif, flags, flow, dims)
            end if
            flags=(/0,1,0/)
            call compute_face_state(qp, y_qp_l, y_qp_r, flags, scheme%jlimiter_switch, scheme%jtlimiter_switch, dims)
            if(scheme%jPB_switch==1)then
              call pressure_based_switching(qp, y_qp_l, y_qp_r, pdif, flags, flow, dims)
            end if
            flags=(/0,0,1/)
            call compute_face_state(qp, z_qp_l, z_qp_r, flags, scheme%klimiter_switch, scheme%ktlimiter_switch, dims)
            if(scheme%kPB_switch==1)then
              call pressure_based_switching(qp, z_qp_l, z_qp_r, pdif, flags, flow, dims)
            end if

        end subroutine compute_muscl_states