compute_fluxes Subroutine

public subroutine compute_fluxes(F, G, H, x_qp_l, x_qp_r, y_qp_l, y_qp_r, z_qp_l, z_qp_r, Ifaces, Jfaces, Kfaces, flow, bc, dims)

Call to compute fluxes throught faces in each direction

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(inout), dimension(:, :, :, :):: F

Store fluxes throught the I faces

real(kind=wp), intent(inout), dimension(:, :, :, :):: G

Store fluxes throught the J faces

real(kind=wp), intent(inout), dimension(:, :, :, :):: H

Store fluxes throught the K faces

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

type(facetype), intent(in), dimension(-2:dims%imx+3,-2:dims%jmx+2,-2:dims%kmx+2):: Ifaces

Store face quantites for I faces

type(facetype), intent(in), dimension(-2:dims%imx+2,-2:dims%jmx+3,-2:dims%kmx+2):: Jfaces

Store face quantites for J faces

type(facetype), intent(in), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+3):: Kfaces

Store face quantites for K faces

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

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

type(boundarytype), intent(in) :: bc

boundary conditions and fixed values

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

Extent of the domain:imx,jmx,kmx


Calls

proc~~compute_fluxes~5~~CallsGraph proc~compute_fluxes~5 compute_fluxes debugcall debugcall proc~compute_fluxes~5->debugcall proc~compute_flux~5 compute_flux proc~compute_fluxes~5->proc~compute_flux~5 proc~compute_flux~5->debugcall

Contents

Source Code


Source Code

        subroutine compute_fluxes(F,G,H, x_qp_l, x_qp_r, y_qp_l, y_qp_r, z_qp_l, z_qp_r, Ifaces, Jfaces, Kfaces, flow, bc, dims)
        !subroutine compute_fluxes(F,G,H, flow, dims)
          !< Call to compute fluxes throught faces in each direction
            
            implicit none
            type(extent), intent(in) :: dims
            !< Extent of the domain:imx,jmx,kmx
            type(flowtype), intent(in) :: flow
            !< Information about fluid flow: freestream-speed, ref-viscosity,etc.
            real(wp), dimension(:, :, :, :), intent(inout) :: F
            !< Store fluxes throught the I faces
            real(wp), dimension(:, :, :, :), intent(inout) :: G
            !< Store fluxes throught the J faces
            real(wp), dimension(:, :, :, :), intent(inout) :: H
            !< Store fluxes throught the K faces
            type(facetype), dimension(-2:dims%imx+3,-2:dims%jmx+2,-2:dims%kmx+2), intent(in) :: Ifaces
            !< Store face quantites for I faces 
            type(facetype), dimension(-2:dims%imx+2,-2:dims%jmx+3,-2:dims%kmx+2), intent(in) :: Jfaces
            !< Store face quantites for J faces 
            type(facetype), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+3), intent(in) :: Kfaces
            !< Store face quantites for K faces 
            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 
            type(boundarytype), intent(in) :: bc
            !< boundary conditions and fixed values
            integer, dimension(3) :: flags

            
            DebugCall('compute_fluxes')

            flags=(/1,0,0/)
            call compute_flux(F, x_qp_l, x_qp_r, Ifaces, flags, flow, bc, dims)
            if (any(isnan(F))) then
              Fatal_error
            end if    

            flags=(/0,1,0/)
            call compute_flux(G, y_qp_l, y_qp_r, Jfaces, flags, flow, bc, dims)
            if (any(isnan(G))) then 
              Fatal_error
            end if    
            
            if(dims%kmx==2) then
              H = 0.
            else
              flags=(/0,0,1/)
              call compute_flux(H, z_qp_l, z_qp_r, Kfaces, flags, flow, bc, dims)
            end if
            if (any(isnan(H))) then
              Fatal_error
            end if

        end subroutine compute_fluxes