compute_gradient Subroutine

private subroutine compute_gradient(grad, var, cells, Ifaces, Jfaces, Kfaces, dir, dims)

Generalized subroutine to calculate gradients

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(out), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2):: grad
real(kind=wp), intent(in), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2):: var
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

character(len=*), intent(in) :: dir
type(extent), intent(in) :: dims

Called by

proc~~compute_gradient~~CalledByGraph proc~compute_gradient compute_gradient proc~find_ccnormal find_CCnormal proc~find_ccnormal->proc~compute_gradient proc~find_dccvn find_DCCVn proc~find_dccvn->proc~compute_gradient proc~setupcc setupCC proc~setupcc->proc~find_ccnormal proc~add_sst_source_lctm2015 add_sst_source_lctm2015 proc~add_sst_source_lctm2015->proc~find_dccvn proc~setup_solver setup_solver proc~setup_solver->proc~setupcc proc~add_source_term_residue add_source_term_residue proc~add_source_term_residue->proc~add_sst_source_lctm2015 proc~get_total_conservative_residue get_total_conservative_Residue proc~get_total_conservative_residue->proc~add_source_term_residue proc~start_run start_run proc~start_run->proc~setup_solver program~main main program~main->proc~start_run proc~iterate_one_more_time_step iterate_one_more_time_step program~main->proc~iterate_one_more_time_step proc~get_next_solution get_next_solution proc~get_next_solution->proc~get_total_conservative_residue proc~iterate_one_more_time_step->proc~get_next_solution

Contents

Source Code


Source Code

    subroutine compute_gradient(grad, var, cells, Ifaces, Jfaces, Kfaces, dir, dims)
      !< Generalized subroutine to calculate gradients
      implicit none
      type(extent), intent(in) :: dims
      real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(out) :: grad
      real(wp), dimension(-2:dims%imx+2,-2:dims%jmx+2,-2:dims%kmx+2), intent(in) :: var
      character(len=*)                                          , intent(in) :: dir
      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
      
      integer :: i
      integer :: j
      integer :: k

      ! initialize
      grad = 0.0

      select case(dir)
        case('x')
          do k=0,dims%kmx
            do j=0,dims%jmx
              do i=0,dims%imx
                grad(i,j,k) =(-(var(i-1,j  ,k  )+var(i,j,k))*Ifaces(i,j,k)%nx*Jfaces(i,j,k)%A &
                              -(var(i  ,j-1,k  )+var(i,j,k))*Ifaces(i,j,k)%ny*Jfaces(i,j,k)%A &
                              -(var(i  ,j  ,k-1)+var(i,j,k))*Ifaces(i,j,k)%nz*Jfaces(i,j,k)%A &
                              +(var(i+1,j  ,k  )+var(i,j,k))*Ifaces(i+1,j  ,k  )%nx*Ifaces(i+1,j  ,k  )%A &
                              +(var(i  ,j+1,k  )+var(i,j,k))*Ifaces(i  ,j+1,k  )%ny*Ifaces(i  ,j+1,k  )%A &
                              +(var(i  ,j  ,k+1)+var(i,j,k))*Ifaces(i  ,j  ,k+1)%nz*Ifaces(i  ,j  ,k+1)%A &
                             )/(2*cells(i,j,k)%volume)
              end do
            end do
          end do
        case('y')
          do k=0,dims%kmx
            do j=0,dims%jmx
              do i=0,dims%imx
                grad(i,j,k) =(-(var(i-1,j  ,k  )+var(i,j,k))*Jfaces(i,j,k)%nx*Jfaces(i,j,k)%A &
                              -(var(i  ,j-1,k  )+var(i,j,k))*Jfaces(i,j,k)%ny*Jfaces(i,j,k)%A &
                              -(var(i  ,j  ,k-1)+var(i,j,k))*Jfaces(i,j,k)%nz*Jfaces(i,j,k)%A &
                              +(var(i+1,j  ,k  )+var(i,j,k))*Jfaces(i+1,j  ,k  )%nx*Jfaces(i+1,j  ,k  )%A &
                              +(var(i  ,j+1,k  )+var(i,j,k))*Jfaces(i  ,j+1,k  )%ny*Jfaces(i  ,j+1,k  )%A &
                              +(var(i  ,j  ,k+1)+var(i,j,k))*Jfaces(i  ,j  ,k+1)%nz*Jfaces(i  ,j  ,k+1)%A &
                             )/(2*cells(i,j,k)%volume)
              end do
            end do
          end do
        case('z')
          do k=0,dims%kmx
            do j=0,dims%jmx
              do i=0,dims%imx
                grad(i,j,k) =(-(var(i-1,j  ,k  )+var(i,j,k))*Kfaces(i,j,k)%nx*Kfaces(i,j,k)%A &
                              -(var(i  ,j-1,k  )+var(i,j,k))*Kfaces(i,j,k)%ny*Kfaces(i,j,k)%A &
                              -(var(i  ,j  ,k-1)+var(i,j,k))*Kfaces(i,j,k)%nz*Kfaces(i,j,k)%A &
                              +(var(i+1,j  ,k  )+var(i,j,k))*Kfaces(i+1,j  ,k  )%nx*Kfaces(i+1,j  ,k  )%A &
                              +(var(i  ,j+1,k  )+var(i,j,k))*Kfaces(i  ,j+1,k  )%ny*Kfaces(i  ,j+1,k  )%A &
                              +(var(i  ,j  ,k+1)+var(i,j,k))*Kfaces(i  ,j  ,k+1)%nz*Kfaces(i  ,j  ,k+1)%A &
                             )/(2*cells(i,j,k)%volume)
              end do
            end do
          end do
        case DEFAULT
          print*, "ERROR: gradient direction error"
          Fatal_error
      end select
      if(any(isnan(grad)))then
        Fatal_error
      end if

    end subroutine compute_gradient