Call PPM face-state reconstruction for each face with optional call for remove extrema based on input limter switch and call pressure based switching based on input pressure based switch
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(in), | 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(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var) | :: | x_qp_l | Store primitive state at the I-face |
|
real(kind=wp), | intent(inout), | dimension(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var) | :: | x_qp_r | Store primitive state at the I-face |
|
real(kind=wp), | intent(inout), | dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var) | :: | y_qp_l | Store primitive state at the J-face |
|
real(kind=wp), | intent(inout), | dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var) | :: | y_qp_r | Store primitive state at the J-face |
|
real(kind=wp), | intent(inout), | dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var) | :: | z_qp_l | Store primitive state at the K-face |
|
real(kind=wp), | intent(inout), | dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var) | :: | z_qp_r | Store primitive state at the K-face |
|
real(kind=wp), | intent(inout), | dimension(0:dims%imx,0:dims%jmx,0:dims%kmx) | :: | pdif | pressure difference |
|
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 |
subroutine compute_ppm_states(qp, x_qp_l, x_qp_r, y_qp_l, y_qp_r, z_qp_l, z_qp_r, pdif, scheme, flow, dims)
!< Call PPM face-state reconstruction for each face
!< with optional call for remove extrema based on
!< input limter switch and call pressure based switching
!< based on input pressure based switch
implicit none
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):: qp
!< Store primitive variable at cell center
real(wp), dimension(0:dims%imx+1,1:dims%jmx-1,1:dims%kmx-1,1:dims%n_var), intent(inout) :: x_qp_l, x_qp_r
!< Store primitive state at the I-face
real(wp), dimension(1:dims%imx-1,0:dims%jmx+1,1:dims%kmx-1,1:dims%n_var), intent(inout) :: y_qp_l, y_qp_r
!< Store primitive state at the J-face
real(wp), dimension(1:dims%imx-1,1:dims%jmx-1,0:dims%kmx+1,1:dims%n_var), intent(inout) :: z_qp_l, z_qp_r
!< Store primitive state at the K-face
real(wp), dimension(0:dims%imx,0:dims%jmx,0:dims%kmx), intent(inout) :: pdif
!< pressure difference
type(schemetype), intent(in) :: scheme
!< finite-volume Schemes
type(flowtype), intent(in) :: flow
!< Information about fluid flow: freestream-speed, ref-viscosity,etc.
integer, dimension(3) :: flags
!< flags for different directions
flags=(/1,0,0/)
call compute_face_estimates(qp, x_qp_l, x_qp_r, flags, dims)
if(scheme%ilimiter_switch==1)then
call remove_extrema(qp, x_qp_l, x_qp_r, flags, dims)
end if
if (scheme%iPB_switch==1)then
call pressure_based_switching(qp, x_qp_l, x_qp_r, pdif, flags, flow, dims)
end if
flags=(/0,1,0/)
call compute_face_estimates(qp, y_qp_l, y_qp_r, flags, dims)
if(scheme%jlimiter_switch==1)then
call remove_extrema(qp, y_qp_l, y_qp_r, flags, dims)
end if
if (scheme%jPB_switch==1)then
call pressure_based_switching(qp, y_qp_l, y_qp_r, pdif, flags, flow, dims)
end if
flags=(/0,0,1/)
call compute_face_estimates(qp, z_qp_l, z_qp_r, flags, dims)
if(scheme%klimiter_switch==1)then
call remove_extrema(qp, z_qp_l, z_qp_r, flags, dims)
end if
if (scheme%kPB_switch==1)then
call pressure_based_switching(qp, z_qp_l, z_qp_r, pdif, flags, flow, dims)
end if
end subroutine compute_ppm_states