marker_set_Filter Subroutine

private impure subroutine marker_set_Filter(this, var, field)

Computes a filtered quantity on the Eulerian grid.

Type Bound

marker_set

Arguments

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

A collection of tessellation elements

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

Variable to compute

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

Filtered quantity


Calls

proc~~marker_set_filter~~CallsGraph proc~marker_set_filter marker_set%marker_set_Filter extrapolate extrapolate proc~marker_set_filter->extrapolate f f proc~marker_set_filter->f markers markers proc~marker_set_filter->markers n n proc~marker_set_filter->n proc~eulerian_obj_addupghostcells eulerian_obj_base%eulerian_obj_AddUpGhostCells proc~marker_set_filter->proc~eulerian_obj_addupghostcells v v proc~marker_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

Called by

proc~~marker_set_filter~~CalledByGraph proc~marker_set_filter marker_set%marker_set_Filter proc~marker_set_computesolidvolfrac marker_set%marker_set_ComputeSolidVolFrac proc~marker_set_computesolidvolfrac->proc~marker_set_filter proc~marker_set_getibforcing marker_set%marker_set_GetIBForcing proc~marker_set_getibforcing->proc~marker_set_filter proc~marker_set_updatenormals marker_set%marker_set_UpdateNormals proc~marker_set_updatenormals->proc~marker_set_filter proc~marker_set_updatesdf marker_set%marker_set_UpdateSDF proc~marker_set_updatesdf->proc~marker_set_filter proc~respart_set_filter ResPart_set%ResPart_set_Filter proc~respart_set_filter->proc~marker_set_filter proc~solid_set_filter solid_set%solid_set_Filter proc~solid_set_filter->proc~marker_set_filter proc~respart_set_getibforcing ResPart_set%ResPart_set_GetIBForcing proc~respart_set_getibforcing->proc~marker_set_getibforcing proc~respart_set_updatenormals ResPart_set%ResPart_set_UpdateNormals proc~respart_set_updatenormals->proc~marker_set_updatenormals proc~respart_set_updatesdf ResPart_set%ResPart_set_UpdateSDF proc~respart_set_updatesdf->proc~marker_set_updatesdf proc~cdifs_obj_advancesolutionib cdifs_obj_AdvanceSolutionIB proc~cdifs_obj_advancesolutionib->proc~respart_set_getibforcing proc~cdifs_obj_advancesolutionib->proc~respart_set_updatesdf proc~cdifs_obj_computesolidvf cdifs_obj%cdifs_obj_ComputeSolidVF proc~cdifs_obj_computesolidvf->proc~respart_set_updatenormals proc~cdifs_obj_preparesolver cdifs_obj_PrepareSolver proc~cdifs_obj_preparesolver->proc~respart_set_updatesdf proc~cdifs_obj_preparesolver->proc~cdifs_obj_computesolidvf proc~grans_obj_computesolidvf grans_obj%grans_obj_ComputeSolidVF proc~grans_obj_computesolidvf->proc~respart_set_updatenormals proc~grans_obj_preparesolver grans_obj_PrepareSolver proc~grans_obj_preparesolver->proc~respart_set_updatesdf proc~grans_obj_preparesolver->proc~grans_obj_computesolidvf interface~grans_obj_writeoutputdata grans_obj%grans_obj_WriteOutputData proc~grans_obj_preparesolver->interface~grans_obj_writeoutputdata proc~grans_obj_writeoutputdata grans_obj_WriteOutputData proc~grans_obj_writeoutputdata->proc~respart_set_updatesdf proc~grans_obj_writeoutputdata->proc~grans_obj_computesolidvf interface~cdifs_obj_preparesolver cdifs_obj%cdifs_obj_PrepareSolver interface~cdifs_obj_preparesolver->proc~cdifs_obj_preparesolver interface~grans_obj_preparesolver grans_obj%grans_obj_PrepareSolver interface~grans_obj_preparesolver->proc~grans_obj_preparesolver interface~grans_obj_writeoutputdata->proc~grans_obj_writeoutputdata proc~cdifs_obj_advancesolution cdifs_obj_AdvanceSolution proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionib proc~cdifs_obj_advancesolutionrp cdifs_obj_AdvanceSolutionRP proc~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolutionrp proc~cdifs_obj_advancesolutionrp->proc~cdifs_obj_computesolidvf interface~cdifs_obj_advancesolution cdifs_obj%cdifs_obj_AdvanceSolution interface~cdifs_obj_advancesolution->proc~cdifs_obj_advancesolution

Source Code

    impure subroutine marker_set_Filter(this,var,field)
      !> Computes a filtered quantity on the Eulerian grid.
      implicit none
      class(marker_set),    intent(inout) :: this                              !! A collection of tessellation elements
      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

      ! Initialize field
      field%cell=0.0_wp
      select type(markers=>this%p)
      type is (marker_obj)
        select case(trim(adjustl(var)))
        case ('SA')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA
                  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=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(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=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(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=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(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 ('F1')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%f(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 ('F2')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%f(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 ('F3')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%f(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 ('N1')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%n(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 ('N2')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%n(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 ('N3')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*markers(n)%n(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 ('V.N')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*sum(markers(n)%v(:)*markers(n)%n(:))
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do
        case ('VISC11')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(1)*markers(n)%n(1)+markers(n)%v(1)*markers(n)%n(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 ('VISC22')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(2)*markers(n)%n(2)+markers(n)%v(2)*markers(n)%n(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 ('VISC33')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(3)*markers(n)%n(3)+markers(n)%v(3)*markers(n)%n(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 ('VISC12')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(1)*markers(n)%n(2)+markers(n)%v(2)*markers(n)%n(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 ('VISC13')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(1)*markers(n)%n(3)+markers(n)%v(3)*markers(n)%n(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 ('VISC23')

          do n=1,this%count_
            ! Get a bump function centered
            ! on the marker
            ! -------------------------- !
            slo=markers(n)%c-this%stib/2
            shi=markers(n)%c+this%stib/2
            call markers(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=markers(n)%SA*(markers(n)%v(2)*markers(n)%n(3)+markers(n)%v(3)*markers(n)%n(2))
                  field%cell(i,j,k) = field%cell(i,j,k)+ coef*bump(i,j,k)
                end do
              end do
            end do
          end do
        end select
      end select

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

      ! Update ghostcells
      call field%AddUpGhostCells

      return
    end subroutine marker_set_Filter