Reconstruct state at the JMIN boundary face with MUSCL scheme
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | dimension(-2:imx+2, -2:jmx+2, -2:kmx+2, 1:n_var) | :: | qp | ||
type(schemetype), | intent(in) | :: | scheme | |||
type(boundarytype), | intent(in) | :: | bc |
subroutine reconstruct_jmin(qp, scheme, bc)
!< Reconstruct state at the JMIN boundary face with MUSCL scheme
implicit none
real(wp), dimension(-2:imx+2, -2:jmx+2, -2:kmx+2, 1:n_var), intent(in) :: qp
type(schemetype), intent(in) :: scheme
type(boundarytype), intent(in) :: bc
integer :: i, j, k, l
real(wp) :: psi1, psi2, fd, bd, r
real(wp) :: kappa, phi
phi = 1.0
kappa = 1./3.
switch_L=scheme%jlimiter_switch
if (ppm_flag==1) then
do l = 1, n_var
if(l>=6) switch_L=scheme%jtlimiter_switch
do k = 1, kmx - 1
do j = 1, 1
do i = 1, imx - 1
fd = qp(i, j+1, k, l) - qp(i, j, k, l)
bd = qp(i, j, k, l) - qp(i, j-1, k, l)
r = fd / bd
psi1 = max(0., min(2*r, (2 + r)/3., 2.))
psi1 = (1 - (1 - psi1)*switch_L )
r = bd / fd
psi2 = max(0., min(2*r, (2 + r)/3., 2.))
psi2 = (1 - (1 - psi2)*switch_L )
! right face of first j cell
y_qp_left(i, j+1, k, l) = qp(i, j, k, l) + 0.25*phi* &
(((1-kappa) * psi1 * bd) + ((1+kappa) * psi2 * fd))
! left face of first j cell
y_qp_right(i, j, k, l) = qp(i, j, k, l) - 0.25*phi* &
(((1+kappa) * psi1 * bd) + ((1-kappa) * psi2 * fd))
end do
end do
end do
end do
end if
if(bc%jmin_id==-8 .or. bc%jmin_id==-9)then
y_qp_left(1:imx-1,1,1:kmx-1,1:n_var) = qp(1:imx-1,0,1:kmx-1,1:n_var)
y_qp_right(1:imx-1,1,1:kmx-1,1:n_var) = qp(1:imx-1,0,1:kmx-1,1:n_var)
else
y_qp_left(1:imx-1,1,1:kmx-1,1:n_var) = 0.5*(qp(1:imx-1,0,1:kmx-1,1:n_var)&
+qp(1:imx-1,1,1:kmx-1,1:n_var))
end if
end subroutine reconstruct_jmin