lagrangian_obj_Interpolate Function

private function lagrangian_obj_Interpolate(this, l_filter, slo, shi, block, g1in, f) result(inter)

Interpolates a field f defined on an Eulerian stencil to the location of a lagrangian object

Type Bound

lagrangian_obj

Arguments

Type IntentOptional Attributes Name
class(lagrangian_obj), intent(in) :: this

A Lagrangian object

real(kind=wp), intent(in) :: l_filter

Filter size

integer, intent(in) :: slo(3)

Stencil lower bound

integer, intent(in) :: shi(3)

Stencil higher bound

type(block_obj), intent(in), pointer :: block

A block object

procedure(kernel_1D), intent(in), pointer :: g1in

Filter kernel

real(kind=wp), intent(in) :: f(slo(1):shi(1),slo(2):shi(2),slo(3):shi(3))

Quantity to interpolate

Return Value real(kind=wp)


Source Code

    function lagrangian_obj_Interpolate(this,l_filter,slo,shi,block,g1in,f) result(inter)
      !> Interpolates a field f defined on an Eulerian stencil to the
      ! location of a lagrangian object
      implicit none
      class(lagrangian_obj), intent(in) :: this                                !! A Lagrangian object
      real(wp),              intent(in) :: l_filter                            !! Filter size
      integer,               intent(in) :: slo(3)                              !! Stencil lower bound
      integer,               intent(in) :: shi(3)                              !! Stencil higher bound
      type(block_obj),       intent(in), &
                               pointer  :: block                               !! A block object
      procedure(kernel_1D),  intent(in), &
                               pointer  :: g1in                                !! Filter kernel
      real(wp),              intent(in) :: f(slo(1):shi(1),slo(2):shi(2),slo(3):shi(3))                     
                                                                               !! Quantity to interpolate                                          
      ! Work variables
      real(wp):: inter
      integer :: i,j,k
      real(wp):: bump,norm
      real(wp):: dVol

      inter=0.0_wp
      norm =0.0_wp
      do k=slo(3),shi(3)
        do j=slo(2),shi(2)
          do i=slo(1),shi(1)
            ! Cell volume
            dVol = (block%x(i+1)-block%x(i))*(block%y(j+1)-block%y(j))*(block%z(k+1)-block%z(k))

            bump = g1in((block%xm(i)-this%p(1))/l_filter)/l_filter &
                 * g1in((block%ym(j)-this%p(2))/l_filter)/l_filter &
                 * g1in((block%zm(k)-this%p(3))/l_filter)/l_filter

            inter = inter +f(i,j,k)*bump*dVol
            norm  = norm + bump*dVol
          end do
        end do
      end do

      inter=inter/norm

      return
    end function lagrangian_obj_Interpolate