compute_global_time_step Subroutine

private subroutine compute_global_time_step(qp, delta_t, cells, Ifaces, Jfaces, Kfaces, CFL, scheme, flow, dims)

Compute a common time step to be used at all cell centers

Global time stepping is generally used to get time accurate solutions; transients can be studied by employing this strategy.


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), target:: qp

Store primitive variable at cell center

real(kind=wp), intent(inout), dimension(1:dims%imx-1, 1:dims%jmx-1, 1:dims%kmx-1):: delta_t

Local time increment value at each cell center

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

Input cell quantities: volume

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

real(kind=wp), intent(in) :: CFL

CFL number

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_global_time_step~~CallsGraph proc~compute_global_time_step compute_global_time_step proc~compute_local_time_step compute_local_time_step proc~compute_global_time_step->proc~compute_local_time_step debugcall debugcall proc~compute_global_time_step->debugcall proc~compute_local_time_step->debugcall proc~add_viscous_time add_viscous_time proc~compute_local_time_step->proc~add_viscous_time proc~add_turbulent_time add_turbulent_time proc~compute_local_time_step->proc~add_turbulent_time proc~add_viscous_time->debugcall proc~add_turbulent_time->debugcall

Called by

proc~~compute_global_time_step~~CalledByGraph proc~compute_global_time_step compute_global_time_step proc~compute_time_step compute_time_step proc~compute_time_step->proc~compute_global_time_step proc~get_next_solution get_next_solution proc~get_next_solution->proc~compute_time_step 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_global_time_step(qp, delta_t, cells, Ifaces, Jfaces, Kfaces, CFL, scheme, flow, dims)
            !< Compute a common time step to be used at all cell centers
            !<
            !< Global time stepping is generally used to get time 
            !< accurate solutions; transients can be studied by 
            !< employing this strategy.
            !<-----------------------------------------------------------

            implicit none
            real(wp), intent(in) :: CFL
            !< CFL number
            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
            real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), intent(in), target :: qp
            !< Store primitive variable at cell center
            real(wp) , dimension(1:dims%imx-1, 1:dims%jmx-1, 1:dims%kmx-1), intent(inout) :: delta_t
            !< Local time increment value at each cell center
            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
            !< 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 
            
            DebugCall('compute_global_time_step')

            if (scheme%global_time_step > 0) then
                delta_t = scheme%global_time_step
            else
              call compute_local_time_step(qp, delta_t, cells, Ifaces, Jfaces, Kfaces, CFL, scheme, flow, dims)
                ! The global time step is the minimum of all the local time
                ! steps.
                delta_t = minval(delta_t)
            end if

        end subroutine compute_global_time_step