apply_gradient_bc Subroutine

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

Call same subroutine for all the face Apply/set value of all gradient in the ghost cells gradqp_G = (qp_I - qp_G)Area_Wunit_normal_G/(volume_G) volume_G = volume_I

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):: qp

Input variable of which graident is required

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

Intput Temperature variable

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(boundarytype), intent(in) :: bc
type(extent), intent(in) :: dims

Extent of the domain:imx,jmx,kmx


Calls

proc~~apply_gradient_bc~~CallsGraph proc~apply_gradient_bc apply_gradient_bc debugcall 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~~apply_gradient_bc~~CalledByGraph proc~apply_gradient_bc apply_gradient_bc proc~evaluate_all_gradients evaluate_all_gradients proc~evaluate_all_gradients->proc~apply_gradient_bc 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


Source Code

    subroutine apply_gradient_bc(qp, temp, cells, Ifaces, Jfaces, Kfaces, bc, dims)
      !< Call same subroutine for all the face
      !< Apply/set value of all gradient in the ghost cells
      !< gradqp_G = (qp_I - qp_G)*Area_W*unit_normal_G/(volume_G)
      !< volume_G = volume_I
      !-----------------------------------------------------------
      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) :: qp
      !< Input variable of which graident is required
      real(wp), dimension(-2:dims%imx+2, -2:dims%jmx+2, -2:dims%kmx+2), intent(in) :: Temp
      !< Intput Temperature variable
      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
      type(singlesub) :: domain

      DebugCall('apply_gradient_bc')


      domain%imin = 1
      domain%imax = 1
      domain%jmin = 1
      domain%jmax = dims%jmx-1
      domain%kmin = 1
      domain%kmax = dims%kmx-1
      domain%il   = 1; domain%jl = 0; domain%kl = 0
      domain%iu   = 0; domain%ju = 0; domain%ku = 0
      domain%sig  = 1
      call apply_gradient_bc_face(qp, temp, cells, Ifaces, dims, domain, bc%imin_id, bc%fixed_wall_temperature(1))

      domain%imin = dims%imx
      domain%imax = dims%imx
      domain%jmin = 1
      domain%jmax = dims%jmx-1
      domain%kmin = 1
      domain%kmax = dims%kmx-1
      domain%il   = 0; domain%jl = 0; domain%kl = 0
      domain%iu   = 1; domain%ju = 0; domain%ku = 0
      domain%sig  = -1
      call apply_gradient_bc_face(qp, temp, cells, Ifaces, dims, domain, bc%imax_id, bc%fixed_wall_temperature(2))

      domain%imin = 1
      domain%imax = dims%imx-1
      domain%jmin = 1
      domain%jmax = 1
      domain%kmin = 1
      domain%kmax = dims%kmx-1
      domain%il   = 0; domain%jl = 1; domain%kl = 0
      domain%iu   = 0; domain%ju = 0; domain%ku = 0
      domain%sig  = 1
      call apply_gradient_bc_face(qp, temp, cells, Jfaces, dims, domain, bc%jmin_id, bc%fixed_wall_temperature(3))

      domain%imin = 1
      domain%imax = dims%imx-1
      domain%jmin = dims%jmx
      domain%jmax = dims%jmx
      domain%kmin = 1
      domain%kmax = dims%kmx-1
      domain%il   = 0; domain%jl = 0; domain%kl = 0
      domain%iu   = 0; domain%ju = 1; domain%ku = 0
      domain%sig  = -1
      call apply_gradient_bc_face(qp, temp, cells, Jfaces, dims, domain, bc%jmax_id, bc%fixed_wall_temperature(4))

      domain%imin = 1
      domain%imax = dims%imx-1
      domain%jmin = 1
      domain%jmax = dims%jmx-1
      domain%kmin = 1
      domain%kmax = 1
      domain%il   = 0; domain%jl = 0; domain%kl = 1
      domain%iu   = 0; domain%ju = 0; domain%ku = 0
      domain%sig  = 1
      call apply_gradient_bc_face(qp, temp, cells, Kfaces, dims, domain, bc%kmin_id, bc%fixed_wall_temperature(5))

      domain%imin = 1
      domain%imax = dims%imx-1
      domain%jmin = 1
      domain%jmax = dims%jmx-1
      domain%kmin = dims%kmx
      domain%kmax = dims%kmx
      domain%il   = 0; domain%jl = 0; domain%kl = 0
      domain%iu   = 0; domain%ju = 0; domain%ku = 1
      domain%sig  = -1
      call apply_gradient_bc_face(qp, temp, cells, Kfaces, dims, domain, bc%kmax_id, bc%fixed_wall_temperature(6))
          
    end subroutine apply_gradient_bc