get_total_conservative_Residue Subroutine

private subroutine get_total_conservative_Residue(qp, Temp, cells, residue, F, G, H, Ifaces, Jfaces, Kfaces, control, scheme, flow, bc, dims)

For each iteration it apply boundary conditions, use higher order method to reconstruct state at face, evalute fluxes at each face, calculate inviscid residual, and introuduce additional residual due to viscosity, turbulence and source terms.

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(inout), 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(in), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2):: Temp

Store Temperature variable at cell center

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

Input cell quantities: volume

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

Store residue at each cell-center

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

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

Input varaible which stores I faces' area and unit normal

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

Input varaible which stores J faces' area and unit normal

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

Input varaible which stores K faces' area and unit normal

type(controltype), intent(in) :: control

Control parameters

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

finite-volume Schemes

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~~get_total_conservative_residue~~CallsGraph proc~get_total_conservative_residue get_total_conservative_Residue proc~evaluate_all_gradients evaluate_all_gradients proc~get_total_conservative_residue->proc~evaluate_all_gradients proc~compute_face_interpolant compute_face_interpolant proc~get_total_conservative_residue->proc~compute_face_interpolant proc~calculate_viscosity calculate_viscosity proc~get_total_conservative_residue->proc~calculate_viscosity proc~reconstruct_boundary_state reconstruct_boundary_state proc~get_total_conservative_residue->proc~reconstruct_boundary_state proc~compute_fluxes~7 compute_fluxes proc~get_total_conservative_residue->proc~compute_fluxes~7 proc~add_source_term_residue add_source_term_residue proc~get_total_conservative_residue->proc~add_source_term_residue proc~compute_viscous_fluxes compute_viscous_fluxes proc~get_total_conservative_residue->proc~compute_viscous_fluxes proc~compute_residue compute_residue proc~get_total_conservative_residue->proc~compute_residue proc~apply_interface apply_interface proc~get_total_conservative_residue->proc~apply_interface proc~populate_ghost_primitive populate_ghost_primitive proc~get_total_conservative_residue->proc~populate_ghost_primitive

Called by

proc~~get_total_conservative_residue~~CalledByGraph proc~get_total_conservative_residue get_total_conservative_Residue 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 get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
        !< For each iteration it apply boundary conditions,
        !< use higher order method to reconstruct state at
        !< face, evalute fluxes at each face, calculate 
        !< inviscid residual, and introuduce additional 
        !< residual due to  viscosity, turbulence and source
        !< terms.
        implicit none
        type(controltype), intent(in) :: control
        !< Control parameters
        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
        type(boundarytype), intent(in) :: bc
        !< boundary conditions and fixed values
        real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2,1:dims%n_var), intent(inout):: qp
        !< Store primitive variable at cell center
        real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(in):: Temp
        !< Store Temperature variable at cell center
        real(wp), dimension(:, :, :, :), intent(inout)  :: residue
        !< Store residue at each cell-center
        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(celltype), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(in) :: cells
        !< Input cell quantities: volume
        type(facetype), dimension(-2:dims%imx+3,-2:dims%jmx+2,-2:dims%kmx+2), intent(in) :: Ifaces
        !< Input varaible which stores I faces' area and unit normal
        type(facetype), dimension(-2:dims%imx+2,-2:dims%jmx+3,-2:dims%kmx+2), intent(in) :: Jfaces
        !< Input varaible which stores J faces' area and unit normal
        type(facetype), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+3), intent(in) :: Kfaces
        !< Input varaible which stores K faces' area and unit normal

        call apply_interface(qp, control, bc, dims)
        call populate_ghost_primitive(qp, Ifaces, Jfaces, Kfaces, control, scheme, flow, bc, dims)
        call compute_face_interpolant(qp, cells, scheme, flow, dims)
        call reconstruct_boundary_state(qp, control, scheme, bc, dims)
        call compute_fluxes(F,G,H,Ifaces,Jfaces,Kfaces,scheme, flow, bc, dims)
        if (flow%mu_ref /= 0.0) then
          call evaluate_all_gradients(qp,Temp,cells,Ifaces,Jfaces,Kfaces,scheme,bc,dims)
          call calculate_viscosity(qp, scheme, flow, bc, dims)
          call compute_viscous_fluxes(F, G, H, qp, cells, Ifaces,Jfaces,Kfaces,scheme, flow, dims)
        end if
        call compute_residue(residue, F,G,H,dims)
        call add_source_term_residue(qp, residue, cells, Ifaces,Jfaces,Kfaces,scheme, flow, dims)

      end subroutine get_total_conservative_Residue