evaluate_all_gradients Subroutine

public subroutine evaluate_all_gradients(qp, Temp, cells, Ifaces, Jfaces, Kfaces, scheme, bc, dims)

Call to all the required gradients and apply boundary condition for ghost cell gradients

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2, 1:dims%n_var), Target:: qp

Store primitive variable at cell center

real(kind=wp), intent(in), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2):: Temp

Store Temperature 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+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(schemetype), intent(in) :: scheme

finite-volume Schemes

type(boundarytype), intent(in) :: bc

boundary conditions and fixed values

type(extent), intent(in) :: dims

Extent of the domain:imx,jmx,kmx


Calls

proc~~evaluate_all_gradients~~CallsGraph proc~evaluate_all_gradients evaluate_all_gradients debugcall debugcall proc~evaluate_all_gradients->debugcall proc~compute_gradient_g compute_gradient_G proc~evaluate_all_gradients->proc~compute_gradient_g proc~apply_gradient_bc apply_gradient_bc proc~evaluate_all_gradients->proc~apply_gradient_bc proc~compute_gradient_g->debugcall proc~apply_gradient_bc->debugcall proc~apply_gradient_bc_face apply_gradient_bc_face proc~apply_gradient_bc->proc~apply_gradient_bc_face

Called by

proc~~evaluate_all_gradients~~CalledByGraph proc~evaluate_all_gradients evaluate_all_gradients proc~get_total_conservative_residue get_total_conservative_Residue proc~get_total_conservative_residue->proc~evaluate_all_gradients 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 evaluate_all_gradients(qp, Temp, cells, Ifaces, Jfaces, Kfaces, scheme, bc, dims)
      !< Call to all the required gradients and 
      !< apply boundary condition for ghost cell
      !< gradients

      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), Target :: qp
      !< Store primitive variable at cell center
      real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2), intent(in) :: Temp
      !< Store Temperature variable at cell center
      type(schemetype) , intent(in) :: scheme
      !< finite-volume Schemes
      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+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(boundarytype), intent(in) :: bc
      !< boundary conditions and fixed values
      real(wp), dimension(:, :, :), pointer :: x_speed      
       !< U pointer, point to slice of qp (:,:,:,2) 
      real(wp), dimension(:, :, :), pointer :: y_speed      
       !< V pointer, point to slice of qp (:,:,:,3) 
      real(wp), dimension(:, :, :), pointer :: z_speed      
       !< W pointer, point to slice of qp (:,:,:,4)
      real(wp), dimension(:, :, :), pointer :: tk   !< TKE/mass
      real(wp), dimension(:, :, :), pointer :: tw   !< Omega
      real(wp), dimension(:, :, :), pointer :: tv   !< SA visocity
      real(wp), dimension(:, :, :), pointer :: tkl  !< KL K-KL method
      real(wp), dimension(:, :, :), pointer :: tgm  !< Intermittency of LCTM2015

      DebugCall('evaluate_all_gradients')

      x_speed(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 2)
      y_speed(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 3)
      z_speed(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 4)

      !call compute_gradient_G(gradu_x, qp(:,:,:,2), cells, Ifaces, Jfaces, Kfaces, dims, 'x')
      call compute_gradient_G(gradu_x, x_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
      call compute_gradient_G(gradv_x, y_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
      call compute_gradient_G(gradw_x, z_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
      call compute_gradient_G(gradT_x, Temp,    cells, Ifaces, Jfaces, Kfaces, dims, 'x')
      !call compute_gradient_G(gradu_y, qp(:,:,:,2), cells, Ifaces, Jfaces, Kfaces, dims, 'y')
      call compute_gradient_G(gradu_y, x_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
      call compute_gradient_G(gradv_y, y_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
      call compute_gradient_G(gradw_y, z_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
      call compute_gradient_G(gradT_y, Temp   , cells, Ifaces, Jfaces, Kfaces, dims, 'y')
      if(dims%kmx>2) then
        call compute_gradient_G(gradu_z, x_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
        call compute_gradient_G(gradv_z, y_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
        call compute_gradient_G(gradw_z, z_speed, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
        !call compute_gradient_T(gradT_z         , cells, Ifaces, Jfaces, Kfaces, dims, 'z')
        call compute_gradient_G(gradT_z, Temp   , cells, Ifaces, Jfaces, Kfaces, dims, 'z')
      else
       gradqp_z=0.0
      end if

      select case (trim(scheme%turbulence))

        case ('none')
          !do nothing
          continue

        case ('sa', 'saBC')
          tv(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 6)
          call compute_gradient_G(gradtv_x, tv, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtv_y, tv, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          if(kmx>2)then
          call compute_gradient_G(gradtv_z, tv, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          end if

        case ('sst', 'sst2003')
          tk(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 6)
          tw(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 7)
          call compute_gradient_G(gradtk_x, tk, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtw_x, tw, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtk_y, tk, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          call compute_gradient_G(gradtw_y, tw, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          if(kmx>2)then
          call compute_gradient_G(gradtk_z, tk, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          call compute_gradient_G(gradtw_z, tw, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          end if

        case ('kkl')
          tk(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 6)
          tkl(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 7)
          call compute_gradient_G(gradtk_x , tk , cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtkl_x, tkl, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtk_y , tk , cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          call compute_gradient_G(gradtkl_y, tkl, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          if(kmx>2)then
          call compute_gradient_G(gradtk_z , tk , cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          call compute_gradient_G(gradtkl_z, tkl, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          end if

        case DEFAULT
          Fatal_error

      end select


      select case(trim(scheme%transition))
        case('lctm2015')
          tgm(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2) => qp(:, :, :, 8)
          call compute_gradient_G(gradtgm_x, tgm, cells, Ifaces, Jfaces, Kfaces, dims, 'x')
          call compute_gradient_G(gradtgm_y, tgm, cells, Ifaces, Jfaces, Kfaces, dims, 'y')
          if(kmx>2)then
          call compute_gradient_G(gradtgm_z, tgm, cells, Ifaces, Jfaces, Kfaces, dims, 'z')
          end if

        case('bc', 'none')
          !do nothing
          continue

        case DEFAULT
          Fatal_error

      end Select

      call apply_gradient_bc(qp, temp, cells, Ifaces, Jfaces, Kfaces, bc, dims)

    end subroutine evaluate_all_gradients