Get solution at next time-step using scheme given in the input file.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
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 |
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