find_wall_dist Subroutine

public subroutine find_wall_dist(nodes, dims)

Determine the minimum wall distance from the wall surface node points

Arguments

Type IntentOptional AttributesName
type(nodetype), intent(in), dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3):: nodes
type(extent), intent(in) :: dims

Calls

proc~~find_wall_dist~~CallsGraph proc~find_wall_dist find_wall_dist interface~alloc alloc proc~find_wall_dist->interface~alloc debugcall debugcall proc~find_wall_dist->debugcall proc~alloc_rank2_real alloc_rank2_real interface~alloc->proc~alloc_rank2_real proc~alloc_rank4_real alloc_rank4_real interface~alloc->proc~alloc_rank4_real proc~alloc_rank3_real alloc_rank3_real interface~alloc->proc~alloc_rank3_real proc~alloc_rank2_integer alloc_rank2_integer interface~alloc->proc~alloc_rank2_integer proc~alloc_rank1_integer alloc_rank1_integer interface~alloc->proc~alloc_rank1_integer proc~alloc_rank6_real alloc_rank6_real interface~alloc->proc~alloc_rank6_real proc~alloc_rank5_real alloc_rank5_real interface~alloc->proc~alloc_rank5_real proc~alloc_rank3_integer alloc_rank3_integer interface~alloc->proc~alloc_rank3_integer proc~alloc_rank1_real alloc_rank1_real interface~alloc->proc~alloc_rank1_real

Called by

proc~~find_wall_dist~~CalledByGraph proc~find_wall_dist find_wall_dist proc~setup_solver setup_solver proc~setup_solver->proc~find_wall_dist 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_wall_dist(nodes, dims)
      !< Determine the minimum wall distance from the wall surface node points

      implicit none
      type(extent), intent(in) :: dims
      type(nodetype), dimension(-2:dims%imx+3,-2:dims%jmx+3,-2:dims%kmx+3), intent(in) :: nodes

      integer :: i,j,k,n
      real(wp) :: current_dist
      real(wp), dimension(:,:,:), allocatable :: node_dist
      DebugCall('find_wall_dist')
      call alloc(node_dist,-2,imx+3,-2,jmx+3,-2,kmx+3)

      do k = -2,dims%kmx+3
        do j = -2,dims%jmx+3
          do i = -2,dims%imx+3
            node_dist(i,j,k) = 1.e+20
            do n = 1,n_surfnodes

            current_dist = sqrt((wall_x(n)-nodes(i,j,k)%x)**2&
                               +(wall_y(n)-nodes(i,j,k)%y)**2&
                               +(wall_z(n)-nodes(i,j,k)%z)**2&
                               ) 
            node_dist(i,j,k) = min(node_dist(i,j,k),current_dist)

            end do
          end do
        end do
      end do
      do k=-2,dims%kmx+2
        do j=-2,dims%jmx+2
          do i=-2,dims%imx+2
            dist(i,j,k) = 0.125*(node_dist(i  ,j  ,k  )&
                                +node_dist(i  ,j+1,k  )&
                                +node_dist(i  ,j+1,k+1)&
                                +node_dist(i  ,j  ,k+1)&
                                +node_dist(i+1,j  ,k+1)&
                                +node_dist(i+1,j  ,k  )&
                                +node_dist(i+1,j+1,k  )&
                                +node_dist(i+1,j+1,k+1)&
                                )
          end do
        end do
      end do
      deallocate(node_dist)
      DebugCall('find_wall_dist-> complete')

    end subroutine find_wall_dist