compute_face_interpolant Subroutine

public subroutine compute_face_interpolant(qp, cells, scheme, flow, dims)

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

type(celltype), intent(in), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2):: cells

Cell center quantities: volume, cellCenter

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_face_interpolant~~CallsGraph proc~compute_face_interpolant compute_face_interpolant proc~compute_ppm_states compute_ppm_states proc~compute_face_interpolant->proc~compute_ppm_states proc~extrapolate_cell_averages_to_faces extrapolate_cell_averages_to_faces proc~compute_face_interpolant->proc~extrapolate_cell_averages_to_faces proc~compute_muscl_states compute_muscl_states proc~compute_face_interpolant->proc~compute_muscl_states proc~compute_weno_nm_states compute_weno_NM_states proc~compute_face_interpolant->proc~compute_weno_nm_states proc~compute_weno_states compute_weno_states proc~compute_face_interpolant->proc~compute_weno_states proc~pressure_based_switching~2 pressure_based_switching proc~compute_ppm_states->proc~pressure_based_switching~2 proc~remove_extrema remove_extrema proc~compute_ppm_states->proc~remove_extrema proc~compute_face_estimates compute_face_estimates proc~compute_ppm_states->proc~compute_face_estimates debugcall debugcall proc~extrapolate_cell_averages_to_faces->debugcall 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 proc~compute_face_states~2 compute_face_states proc~compute_weno_nm_states->proc~compute_face_states~2 proc~compute_face_states compute_face_states proc~compute_weno_states->proc~compute_face_states proc~pressure_based_switching~2->debugcall proc~remove_extrema->debugcall proc~compute_face_estimates->debugcall proc~compute_face_state->debugcall proc~pressure_based_switching->debugcall

Called by

proc~~compute_face_interpolant~~CalledByGraph proc~compute_face_interpolant compute_face_interpolant 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_interpolant(qp, cells, scheme, flow, dims)
            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
            type(celltype), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(in) :: cells
            !< Cell center quantities: volume, cellCenter
            type(schemetype), intent(in) :: scheme
            !< finite-volume Schemes
            type(flowtype), intent(in) :: flow
            !< Information about fluid flow: freestream-speed, ref-viscosity,etc.
            select case (scheme%interpolant)
                case ("none")
                    call extrapolate_cell_averages_to_faces(qp, dims)
                case ("ppm")
                    call compute_ppm_states(qp,  x_qp_left, x_qp_right, y_qp_left, y_qp_right, z_qp_left, z_qp_right, pdif, scheme, flow, dims)
                case ("muscl")
                    call compute_muscl_states(qp, x_qp_left, x_qp_right, y_qp_left, y_qp_right, z_qp_left, z_qp_right, pdif, scheme, flow, dims)
                case ("weno")
                    call compute_weno_states(qp, x_qp_left, x_qp_right, y_qp_left, y_qp_right, z_qp_left, z_qp_right, dims)
                case ("weno_NM")
                    call compute_weno_NM_states(qp, x_qp_left, x_qp_right, y_qp_left, y_qp_right, z_qp_left, z_qp_right, cells, dims)
                case default
                    Fatal_error
            end select
        end subroutine compute_face_interpolant