get_next_solution Subroutine

public subroutine get_next_solution(qp, Temp, residue, delta_t, cells, F, G, H, Ifaces, Jfaces, Kfaces, control, scheme, flow, bc, dims)

Get solution at next time-step using scheme given in the input file.

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(inout), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2):: Temp

Store Temperature variable at cell center

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

Store residue at each 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 Store residue 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

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

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(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_next_solution~~CallsGraph proc~get_next_solution get_next_solution proc~get_total_conservative_residue get_total_conservative_Residue proc~get_next_solution->proc~get_total_conservative_residue proc~update_with update_with proc~get_next_solution->proc~update_with proc~compute_time_step compute_time_step proc~get_next_solution->proc~compute_time_step proc~update_with_plusgs update_with_plusgs proc~get_next_solution->proc~update_with_plusgs proc~update_with_lusgs update_with_lusgs proc~get_next_solution->proc~update_with_lusgs proc~apply_interface apply_interface proc~get_total_conservative_residue->proc~apply_interface proc~add_source_term_residue add_source_term_residue proc~get_total_conservative_residue->proc~add_source_term_residue proc~compute_residue compute_residue proc~get_total_conservative_residue->proc~compute_residue proc~populate_ghost_primitive populate_ghost_primitive proc~get_total_conservative_residue->proc~populate_ghost_primitive proc~evaluate_all_gradients evaluate_all_gradients proc~get_total_conservative_residue->proc~evaluate_all_gradients 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~compute_face_interpolant compute_face_interpolant proc~get_total_conservative_residue->proc~compute_face_interpolant proc~compute_viscous_fluxes compute_viscous_fluxes proc~get_total_conservative_residue->proc~compute_viscous_fluxes debugcall debugcall proc~compute_time_step->debugcall proc~compute_local_time_step compute_local_time_step proc~compute_time_step->proc~compute_local_time_step proc~compute_global_time_step compute_global_time_step proc~compute_time_step->proc~compute_global_time_step proc~update_simulation_clock update_simulation_clock proc~compute_time_step->proc~update_simulation_clock proc~update_lctm2015 update_lctm2015 proc~update_with_plusgs->proc~update_lctm2015 proc~update_laminar_variables update_laminar_variables proc~update_with_plusgs->proc~update_laminar_variables proc~update_sa_variables update_SA_variables proc~update_with_plusgs->proc~update_sa_variables proc~update_sst_variables update_SST_variables proc~update_with_plusgs->proc~update_sst_variables proc~update_kkl_variables update_KKL_variables proc~update_with_lusgs->proc~update_kkl_variables proc~update_sa_variables~2 update_SA_variables proc~update_with_lusgs->proc~update_sa_variables~2 proc~update_lctm2015~2 update_lctm2015 proc~update_with_lusgs->proc~update_lctm2015~2 proc~update_with_lusgs->debugcall proc~update_laminar_variables~2 update_laminar_variables proc~update_with_lusgs->proc~update_laminar_variables~2 proc~update_sst_variables~2 update_SST_variables proc~update_with_lusgs->proc~update_sst_variables~2

Called by

proc~~get_next_solution~~CalledByGraph proc~get_next_solution get_next_solution 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 get_next_solution(qp, Temp, residue, delta_t, cells, F,G,H, Ifaces, Jfaces, Kfaces, control, scheme, flow, bc, dims)
        !< Get solution at next time-step using scheme
        !< given in the input file.
        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(boundarytype), intent(in) :: bc
        !< boundary conditions and fixed values
        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(inout):: qp
        !< Store primitive variable at cell center
        real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(inout):: Temp
        !< Store Temperature variable at cell center
        real(wp), dimension(:, :, :, :), intent(inout)  :: residue
        !< Store residue at each 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
        !< 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
        !< 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) :: CFL 
        CFL = control%CFL
        !finding the updated Temperature using ideal gas law
        !T=P/(R_gas*Rho)
        Temp = qp(:,:,:,5)/(flow%R_gas*qp(:,:,:,1) )
        select case (trim(scheme%time_step_accuracy))
            case ("none")
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1. ,1., .FALSE.) 
            case ("RK4")
              R_store=0.
              U_store = qp
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 0.5  , 1., .FALSE., R_store, U_store) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 0.5  , 2., .FALSE., R_store, U_store) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 2., .FALSE., R_store, U_store) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1./6., 1., .TRUE. , R_store, U_store) 
            case("RK2")
              R_store=0.
              U_store = qp
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 0.5  , 1., .FALSE., R_store, U_store) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 0.5  , 1., .TRUE., R_store, U_store) 
            case ("TVDRK3")
              U_store = qp
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 1.) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 1.) 
              qp = 0.75*U_store + 0.25*qp
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 1.) 
              qp = (1./3.)*U_store + (2./3.)*qp
            case ("TVDRK2")
              U_store = qp
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 1.) 
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call update_with(qp, residue, delta_t, cells, scheme, flow, "conservative", 1.0  , 1.) 
              qp = 0.5*U_store + 0.5*qp
            case ("implicit")
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with_lusgs(qp,residue, delta_t, cells,Ifaces,Jfaces,Kfaces, scheme, dims)
            case ("plusgs")
              call get_total_conservative_Residue(qp, Temp, cells, residue, F,G,H, Ifaces,Jfaces,Kfaces, control, scheme, flow, bc, dims)
              call compute_time_step(qp, delta_t, CFL, cells, Ifaces, Jfaces, Kfaces, scheme, flow, dims) ! has to be after get_..._Residue()
              call update_with_plusgs(qp,delta_t, cells,Ifaces,Jfaces,Kfaces, residue, scheme, dims)
            case default
              Fatal_error
        end select
      end subroutine get_next_solution