Compute viscous fluxes for additianal equations due to LCTM2015 transition model
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=wp), | intent(inout), | dimension(:, :, :, :) | :: | F | Flux array |
|
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 |
|
type(celltype), | intent(in), | dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2) | :: | cells | Input cell quantities: volume |
|
type(facetype), | intent(in), | dimension(-2:dims%imx+2+flags(1),-2:dims%jmx+2+flags(2),-2:dims%kmx+2+flags(3)) | :: | faces | Face quantities: area and unit normal |
|
integer, | intent(in), | dimension(3) | :: | flags | flags for direction switch |
|
type(extent), | intent(in) | :: | dims | Extent of the doamin:imx,jmx,kmx |
subroutine compute_viscous_fluxes_lctm2015(F, qp, cells, faces, flags,dims)
!< Compute viscous fluxes for additianal equations due to LCTM2015 transition model
implicit none
real(wp), dimension(:, :, :, :), intent(inout) :: F
!< Flux array
type(extent), intent(in) :: dims
!< Extent of the doamin:imx,jmx,kmx
integer, dimension(3), intent(in) :: flags
!<flags for direction switch
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
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+2+flags(1),-2:dims%jmx+2+flags(2),-2:dims%kmx+2+flags(3)), intent(in) :: faces
!< Face quantities: area and unit normal
! local variables
real(wp) :: rhoface
real(wp) :: normal_comp
real(wp) :: d_LR
real(wp) :: mu_f
real(wp) :: mut_f
real(wp) :: delx
real(wp) :: dely
real(wp) :: delz
real(wp) :: deltgm
real(wp) :: nx
real(wp) :: ny
real(wp) :: nz
real(wp) :: area
integer :: i, j, k
integer :: ii, jj, kk
!--- sa variable requirement ---!
real(wp) :: dtgmdx, dtgmdy, dtgmdz
ii = flags(1)
jj = flags(2)
kk = flags(3)
!---------------------------------------------------------------------
! Calculating the turbulent viscous fluxes at the faces
!--------------------------------------------------------------------
do k = 1, dims%kmx - 1 + kk
do j = 1, dims%jmx - 1 + jj
do i = 1, dims%imx - 1 + ii
!--- FACE Gradients ---!
! Gradients at face as average of gradients at cell centres
dtgmdx = 0.5 * (gradtgm_x(i-ii, j-jj, k-kk) + gradtgm_x(i, j, k))
dtgmdy = 0.5 * (gradtgm_y(i-ii, j-jj, k-kk) + gradtgm_y(i, j, k))
dtgmdz = 0.5 * (gradtgm_z(i-ii, j-jj, k-kk) + gradtgm_z(i, j, k))
!--- For ODD-EVEN coupling error ---!
! distance between cell center of adjacent cell for the i,j,k face
delx = cells(i, j, k)%centerx - cells(i-ii, j-jj, k-kk)%centerx
dely = cells(i, j, k)%centery - cells(i-ii, j-jj, k-kk)%centery
delz = cells(i, j, k)%centerz - cells(i-ii, j-jj, k-kk)%centerz
d_LR = sqrt(delx*delx + dely*dely + delz*delz)
! difference in state across face
deltgm = qp(i, j, k, 8) - qp(i-ii, j-jj, k-kk, 8)
!normal_comp = ( delta(phi) - (grad(phi).dot.delR) )/magnitudeR
!new grad(phi) = grad(phi) + correction(normal_comp.dot.delR/magnitudeR)
normal_comp = (deltgm - (dtgmdx*delx + dtgmdy*dely + dtgmdz*delz))/d_LR
dtgmdx = dtgmdx + (normal_comp * delx / d_LR)
dtgmdy = dtgmdy + (normal_comp * dely / d_LR)
dtgmdz = dtgmdz + (normal_comp * delz / d_LR)
!--- end of ODD-EVEN coupling correction ---!
rhoface = 0.5 * (qp(i-ii, j-jj, k-kk, 1) + qp(i, j, k, 1))
mu_f = 0.5*(mu(i-ii, j-jj, k-kk) + mu(i, j, k))
mut_f = 0.5*(mu_t(i-ii, j-jj, k-kk) + mu_t(i, j, k))
! calling some element from memory and keep them handy for calculation
nx = faces(i,j,k)%nx
ny = faces(i,j,k)%ny
nz = faces(i,j,k)%nz
area = faces(i,j,k)%A
! adding viscous fluxes to stored convective flux
F(i, j, k, dims%n_var) = F(i, j, k, dims%n_var) - (area*((mu_f + mut_f)*(dtgmdx*nx + dtgmdy*ny + dtgmdz*nz)))
end do
end do
end do
end subroutine compute_viscous_fluxes_lctm2015