particle_set_Filter Subroutine

private impure subroutine particle_set_Filter(this, var, field)

Filters a quantity to the Eulerian grid.

Type Bound

particle_set

Arguments

Type IntentOptional Attributes Name
class(particle_set), intent(inout) :: this

Set of Lagrangian objects

character(len=*), intent(in) :: var

Variable to compute

type(eulerian_obj_r), intent(inout) :: field

Filtered quantity


Calls

proc~~particle_set_filter~~CallsGraph proc~particle_set_filter particle_set%particle_set_Filter extrapolate extrapolate proc~particle_set_filter->extrapolate particles particles proc~particle_set_filter->particles proc~eulerian_obj_addupghostcells eulerian_obj_base%eulerian_obj_AddUpGhostCells proc~particle_set_filter->proc~eulerian_obj_addupghostcells v v proc~particle_set_filter->v proc~eulerian_obj_addupghostcells_x eulerian_obj_base%eulerian_obj_AddUpGhostCells_x proc~eulerian_obj_addupghostcells->proc~eulerian_obj_addupghostcells_x proc~eulerian_obj_addupghostcells_y eulerian_obj_base%eulerian_obj_AddUpGhostCells_y proc~eulerian_obj_addupghostcells->proc~eulerian_obj_addupghostcells_y proc~eulerian_obj_addupghostcells_z eulerian_obj_base%eulerian_obj_AddUpGhostCells_z proc~eulerian_obj_addupghostcells->proc~eulerian_obj_addupghostcells_z proc~eulerian_obj_updateghostcells eulerian_obj_base%eulerian_obj_UpdateGhostCells proc~eulerian_obj_addupghostcells->proc~eulerian_obj_updateghostcells cell cell proc~eulerian_obj_addupghostcells_x->cell mpi_irecv mpi_irecv proc~eulerian_obj_addupghostcells_x->mpi_irecv mpi_isend mpi_isend proc~eulerian_obj_addupghostcells_x->mpi_isend mpi_waitall mpi_waitall proc~eulerian_obj_addupghostcells_x->mpi_waitall proc~eulerian_obj_addupghostcells_y->cell proc~eulerian_obj_addupghostcells_y->mpi_irecv proc~eulerian_obj_addupghostcells_y->mpi_isend proc~eulerian_obj_addupghostcells_y->mpi_waitall proc~eulerian_obj_addupghostcells_z->cell proc~eulerian_obj_addupghostcells_z->mpi_irecv proc~eulerian_obj_addupghostcells_z->mpi_isend proc~eulerian_obj_addupghostcells_z->mpi_waitall proc~eulerian_obj_updateghostcells_x eulerian_obj_base%eulerian_obj_UpdateGhostCells_x proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_x proc~eulerian_obj_updateghostcells_y eulerian_obj_base%eulerian_obj_UpdateGhostCells_y proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_y proc~eulerian_obj_updateghostcells_z eulerian_obj_base%eulerian_obj_UpdateGhostCells_z proc~eulerian_obj_updateghostcells->proc~eulerian_obj_updateghostcells_z proc~eulerian_obj_updateghostcells_x->cell proc~eulerian_obj_updateghostcells_x->mpi_irecv proc~eulerian_obj_updateghostcells_x->mpi_isend proc~eulerian_obj_updateghostcells_x->mpi_waitall proc~eulerian_obj_updateghostcells_y->cell proc~eulerian_obj_updateghostcells_y->mpi_irecv proc~eulerian_obj_updateghostcells_y->mpi_isend proc~eulerian_obj_updateghostcells_y->mpi_waitall proc~eulerian_obj_updateghostcells_z->cell proc~eulerian_obj_updateghostcells_z->mpi_irecv proc~eulerian_obj_updateghostcells_z->mpi_isend proc~eulerian_obj_updateghostcells_z->mpi_waitall

Source Code

    impure subroutine particle_set_Filter(this,var,field)
      !> Filters a quantity to the Eulerian grid.
      implicit none
      class(particle_set),  intent(inout) :: this                              !! Set of Lagrangian objects
      character(len=*),     intent(in)    :: var                               !! Variable to compute
      type(eulerian_obj_r), intent(inout) :: field                             !! Filtered quantity
      ! Work variables
      real(wp),allocatable :: bump(:,:,:)
      real(wp):: coef
      integer :: slo(3),shi(3)
      integer :: n,i,j,k
      real(wp),parameter :: Pi=4.0_wp*atan(1.0_wp)

      ! Initialize field
      field%cell=0.0_wp
      select type(particles=>this%p)
      class is (particle_obj)
        select case(trim(adjustl(var)))
        case ('volume fraction')
          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=particles(n)%c-this%stib/2
            shi=particles(n)%c+this%stib/2
            call particles(n)%extrapolate(this%l_filter,slo,shi,this%block,this%g1ex,bump)

            ! Scale the bump function
            ! -------------------------- !
            do k=slo(3),shi(3)
              do j=slo(2),shi(2)
                do i=slo(1),shi(1)
                  coef=Pi*particles(n)%d**3/6.0_wp
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do

        case ('V1')
          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=particles(n)%c-this%stib/2
            shi=particles(n)%c+this%stib/2
            call particles(n)%extrapolate(this%l_filter,slo,shi,this%block,this%g1ex,bump)

            ! Scale the bump function
            ! -------------------------- !
            do k=slo(3),shi(3)
              do j=slo(2),shi(2)
                do i=slo(1),shi(1)
                  coef=Pi*particles(n)%d**3/6.0_wp*particles(n)%v(1)
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do

        case ('V2')
          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=particles(n)%c-this%stib/2
            shi=particles(n)%c+this%stib/2
            call particles(n)%extrapolate(this%l_filter,slo,shi,this%block,this%g1ex,bump)

            ! Scale the bump function
            ! -------------------------- !
            do k=slo(3),shi(3)
              do j=slo(2),shi(2)
                do i=slo(1),shi(1)
                  coef=Pi*particles(n)%d**3/6.0_wp*particles(n)%v(2)
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do

        case ('V3')
          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=particles(n)%c-this%stib/2
            shi=particles(n)%c+this%stib/2
            call particles(n)%extrapolate(this%l_filter,slo,shi,this%block,this%g1ex,bump)

            ! Scale the bump function
            ! -------------------------- !
            do k=slo(3),shi(3)
              do j=slo(2),shi(2)
                do i=slo(1),shi(1)
                  coef=Pi*particles(n)%d**3/6.0_wp*particles(n)%v(3)
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do
        case default
          call this%parallel%Stop("Unkown quantity to filter")
        end select
      end select

      ! Free memory
      if (allocated(bump))      deallocate(bump)

      ! Update ghostcells
      call field%AddUpGhostCells

      return
    end subroutine particle_set_Filter