find_DCCVn Subroutine

public subroutine find_DCCVn(qp, cells, Ifaces, Jfaces, Kfaces, dims)

Find gradient of the dot product between cell velocity and unit normal

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-2:dims%imx,-2:dims%jmx,-2:dims%kmx,-2:dims%n_var):: qp
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_dccvn~~CallsGraph proc~find_dccvn find_DCCVn proc~compute_gradient compute_gradient proc~find_dccvn->proc~compute_gradient proc~find_ccvn find_CCVn proc~find_dccvn->proc~find_ccvn

Called by

proc~~find_dccvn~~CalledByGraph proc~find_dccvn find_DCCVn proc~add_sst_source_lctm2015 add_sst_source_lctm2015 proc~add_sst_source_lctm2015->proc~find_dccvn 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~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


Source Code

    subroutine find_DCCVn(qp,cells, Ifaces, Jfaces, Kfaces,dims)
      !< Find gradient of the dot product between cell velocity and unit normal
      implicit none
      type(extent), intent(in) :: dims
      real(wp), dimension(-2:dims%imx,-2:dims%jmx,-2:dims%kmx,-2:dims%n_var), intent(in) :: qp
      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 find_CCVn(qp, dims)
      call compute_gradient(DCCVnX, dist, cells, Ifaces, Jfaces, Kfaces, 'x', dims)
      call compute_gradient(DCCVnY, dist, cells, Ifaces, Jfaces, Kfaces, 'y', dims)
      call compute_gradient(DCCVnZ, dist, cells, Ifaces, Jfaces, Kfaces, 'z', dims)
    end subroutine find_DCCVn