Populate the state variables in the ghost cell with particular value based on the boundary conditio being applied at that face
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), target | :: | state | state variables |
|
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 | |||
type(schemetype), | intent(in) | :: | scheme | |||
type(flowtype), | intent(in) | :: | flow | |||
type(boundarytype), | intent(in) | :: | bc | |||
type(extent), | intent(in) | :: | dims |
subroutine populate_ghost_primitive(state, Ifaces, Jfaces, Kfaces, control, scheme, flow, bc, dims)
!< Populate the state variables in the ghost cell
!< with particular value based on the boundary conditio
!< being applied at that face
implicit none
type(extent), intent(in) :: dims
real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), intent(inout), target :: state
!< state variables
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
type(controltype), intent(in) :: control
type(schemetype), intent(in) :: scheme
type(flowtype), intent(in) :: flow
type(boundarytype), intent(in) :: bc
integer :: i
character(len=4) :: face
imx = dims%imx
jmx = dims%jmx
kmx = dims%kmx
n_var = dims%n_var
current_iter = control%current_iter
turbulence = trim(scheme%turbulence)
transition = trim(scheme%transition)
mu_ref = flow%mu_ref
gm = flow%gm
R_gas = flow%R_gas
T_ref = flow%T_ref
sutherland_temp = flow%sutherland_temp
x_speed_inf = flow%x_speed_inf
y_speed_inf = flow%y_speed_inf
z_speed_inf = flow%z_speed_inf
density_inf = flow%density_inf
pressure_inf = flow%pressure_inf
tk_inf = flow%tk_inf
tw_inf = flow%tw_inf
te_inf = flow%te_inf
tv_inf = flow%tv_inf
tgm_inf = flow%tgm_inf
tkl_inf = flow%tkl_inf
qp(-2:imx+2, -2:jmx+2, -2:kmx+2, 1:n_var) => state(:, :, :, :)
density(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 1)
x_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 2)
y_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 3)
z_speed(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 4)
pressure(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 5)
select case (trim(scheme%turbulence))
case ("none")
!include nothing
continue
case ("sst", "sst2003", "bsl", "des-sst", "kw")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
tw(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case ("kkl")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
tkl(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case ("sa", "saBC")
tv(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
case ("ke")
tk(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 6)
te(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, 7)
case ("les")
continue
! todo
case DEFAULT
Fatal_error
end select
! Transition modeling
select case(trim(scheme%transition))
case('lctm2015')
tgm(-2:imx+2, -2:jmx+2, -2:kmx+2) => state(:, :, :, n_var)
! tgm_inf => qp_inf(n_var)
case('bc', 'none')
!do nothing
continue
case DEFAULT
Fatal_error
end Select
do i = 1,6
face_num = i
face = bc%face_names(face_num)
select case(bc%id(face_num))
case(-1)
call supersonic_inlet(face, bc)
case(-2)
call supersonic_outlet(face, bc, dims)
case(-3)
call subsonic_inlet(face, bc, dims)
case(-4)
call subsonic_outlet(face, bc, dims)
case(-5)
call wall(face, bc, dims)
case(-6)
call slip_wall(face, Ifaces, Jfaces, Kfaces, bc, dims)
case(-7)
call pole(face, bc, dims)
case(-8)
call far_field(face, Ifaces, Jfaces, Kfaces, bc, dims)
case(-9)
call periodic_bc(face)
case(-11)
call total_pressure(face, Ifaces, Jfaces, Kfaces, bc, dims)
case Default
if(bc%id(i)>=0 .or. bc%id(i)==-10) then
continue !interface boundary
else
print*, " boundary condition not recognised -> id is :", bc%id(i)
end if
end select
end do
! qp(0,0,:,:) = 0.5*(qp(0,1,:,:)+qp(1,0,:,:))
! qp(0,jmx,:,:) = 0.5*(qp(0,jmx-1,:,:)+qp(1,jmx,:,:))
! qp(imx,0,:,:) = 0.5*(qp(imx,1,:,:)+qp(imx-1,0,:,:))
! qp(imx,jmx,:,:) = 0.5*(qp(imx,jmx-1,:,:)+qp(imx-1,jmx,:,:))
! qp(0,:,0,:) = 0.5*(qp(0,:,1,:)+qp(1,:,0,:))
! qp(0,:,kmx,:) = 0.5*(qp(0,:,kmx-1,:)+qp(1,:,kmx,:))
! qp(imx,:,0,:) = 0.5*(qp(imx,:,1,:)+qp(imx-1,:,0,:))
! qp(imx,:,kmx,:) = 0.5*(qp(imx,:,jmx-1,:)+qp(imx-1,:,jmx,:))
qp(:,0,0,:) = 0.33*(qp(:,1,1,:)+qp(:,0,1,:)+qp(:,1,0,:))
qp(:,0,kmx,:) = 0.33*(qp(:,1,kmx-1,:)+qp(:,0,kmx-1,:)+qp(:,1,kmx,:))
qp(:,jmx,0,:) = 0.33*(qp(:,jmx-1,1,:)+qp(:,jmx,1,:)+qp(:,jmx-1,0,:))
qp(:,jmx,kmx,:) = 0.33*(qp(:,jmx-1,kmx-1,:)+qp(:,jmx,kmx-1,:)+qp(:,jmx-1,kmx,:))
qp(imx,0,:,:) = 0.33*(qp(imx-1,1,:,:)+qp(imx-1,0,:,:)+qp(imx,1,:,:))
qp(0,0,:,:) = 0.33*(qp(1,1,:,:)+qp(1,0,:,:)+qp(0,1,:,:))
qp(0,jmx,:,:) = 0.33*(qp(1,jmx-1,:,:)+qp(1,jmx,:,:)+qp(0,jmx-1,:,:))
qp(imx,jmx,:,:) = 0.33*(qp(imx-1,jmx-1,:,:)+qp(imx-1,jmx,:,:)+qp(imx,jmx-1,:,:))
qp(0,0,0,:) = 0.33*(qp(1,0,0,:) + qp(0,1,0,:) + qp(0,0,1,:))
qp(imx,0,0,:) = 0.33*(qp(imx-1,0,0,:) + qp(imx,1,0,:) + qp(imx,0,1,:))
qp(0,jmx,0,:) = 0.33*(qp(1,jmx,0,:) + qp(0,jmx-1,0,:) + qp(0,jmx,1,:))
qp(0,0,kmx,:) = 0.33*(qp(1,0,kmx,:) + qp(0,1,kmx,:) + qp(0,0,kmx-1,:))
qp(imx,jmx,0,:) = 0.33*(qp(imx-1,jmx,0,:) + qp(imx,jmx-1,0,:) + qp(imx,jmx,1,:))
qp(imx,0,kmx,:) = 0.33*(qp(imx-1,0,kmx,:) + qp(imx,1,kmx,:) + qp(imx,0,kmx-1,:))
qp(0,jmx,kmx,:) = 0.33*(qp(1,jmx,kmx,:) + qp(0,jmx-1,kmx,:) + qp(0,jmx,kmx-1,:))
qp(imx,jmx,kmx,:) = 0.33*(qp(imx-1,jmx,kmx,:) + qp(imx,jmx-1,kmx,:) + qp(imx,jmx,kmx-1,:))
end subroutine populate_ghost_primitive