compute_viscous_fluxes_sa Subroutine

private subroutine compute_viscous_fluxes_sa(F, qp, cells, faces, flags, dims)

Compute viscous fluxes for additianal equations due to SA turbulence model

Arguments

Type IntentOptional AttributesName
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 domain: imx,jmx,kmx


Called by

proc~~compute_viscous_fluxes_sa~~CalledByGraph proc~compute_viscous_fluxes_sa compute_viscous_fluxes_sa proc~compute_viscous_fluxes compute_viscous_fluxes proc~compute_viscous_fluxes->proc~compute_viscous_fluxes_sa proc~get_total_conservative_residue get_total_conservative_Residue proc~get_total_conservative_residue->proc~compute_viscous_fluxes proc~get_next_solution get_next_solution proc~get_next_solution->proc~get_total_conservative_residue proc~iterate_one_more_time_step iterate_one_more_time_step proc~iterate_one_more_time_step->proc~get_next_solution program~main main program~main->proc~iterate_one_more_time_step

Contents


Source Code

    subroutine compute_viscous_fluxes_sa(F, qp, cells, faces, flags, dims)
      !< Compute viscous fluxes for additianal equations due to SA turbulence model
      implicit none
      real(wp), dimension(:, :, :, :), intent(inout) :: F 
      !< Flux array
      type(extent), intent(in) :: dims
      !< Extent of the domain: 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) :: deltv
      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) :: dtvdx, dtvdy, dtvdz

      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
          dtvdx = 0.5 * (gradtv_x(i-ii, j-jj, k-kk) + gradtv_x(i, j, k))
          dtvdy = 0.5 * (gradtv_y(i-ii, j-jj, k-kk) + gradtv_y(i, j, k))
          dtvdz = 0.5 * (gradtv_z(i-ii, j-jj, k-kk) + gradtv_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
          deltv = qp(i, j, k, 6) - qp(i-ii, j-jj, k-kk, 6)

          !normal_comp   = ( delta(phi) - (grad(phi).dot.delR) )/magnitudeR
          !new grad(phi) =  grad(phi) + correction(normal_comp.dot.delR/magnitudeR)
          normal_comp = (deltv - (dtvdx*delx + dtvdy*dely + dtvdz*delz))/d_LR
          dtvdx       =  dtvdx + (normal_comp * delx / d_LR)
          dtvdy       =  dtvdy + (normal_comp * dely / d_LR)
          dtvdz       =  dtvdz + (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*(qp(i-ii, j-jj, k-kk, 6) + qp(i, j, k, 6))*rhoface

          ! 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, 6) = F(i, j, k, 6) - (area*((mu_f + mut_f)*(dtvdx*nx + dtvdy*ny + dtvdz*nz)))/sigma_sa
         
        end do
       end do
      end do
    end subroutine compute_viscous_fluxes_sa