find_CCnormal Subroutine

private subroutine find_CCnormal(cells, Ifaces, Jfaces, Kfaces, dims)

Find the cell-center unit normal

Arguments

Type IntentOptional AttributesName
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(extent), intent(in) :: dims

Calls

proc~~find_ccnormal~~CallsGraph proc~find_ccnormal find_CCnormal proc~compute_gradient compute_gradient proc~find_ccnormal->proc~compute_gradient

Called by

proc~~find_ccnormal~~CalledByGraph proc~find_ccnormal find_CCnormal proc~setupcc setupCC proc~setupcc->proc~find_ccnormal proc~setup_solver setup_solver proc~setup_solver->proc~setupcc proc~start_run start_run proc~start_run->proc~setup_solver program~main main program~main->proc~start_run

Contents

Source Code


Source Code

    subroutine find_CCnormal(cells, Ifaces, Jfaces, Kfaces,dims)
      !< Find the cell-center unit normal
      implicit none
      type(extent), intent(in) :: dims
      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
      call compute_gradient(CCnormalX, dist, cells, Ifaces, Jfaces, Kfaces, 'x', dims)
      call compute_gradient(CCnormalY, dist, cells, Ifaces, Jfaces, Kfaces, 'y', dims)
      call compute_gradient(CCnormalZ, dist, cells, Ifaces, Jfaces, Kfaces, 'z', dims)
      !using already allocated memeory for storing magnitude
      CCVn = sqrt(CCnormalX**2 + CCnormalY**2 + CCnormalZ**2)
      !CCVn hold the magnitude of CCnormal temporaraly and can be 
      !overwritten after next three lines of code.
      CCnormalX = CCnormalX/(CCVn + 1e-12)
      CCnormalY = CCnormalY/(CCVn + 1e-12)
      CCnormalZ = CCnormalZ/(CCVn + 1e-12)
    end subroutine find_CCnormal