check_if_in_profiles_domain_p Subroutine

private subroutine check_if_in_profiles_domain_p(Y_R, Y_PHI, Y_Z, flag)

Arguments

Type IntentOptional AttributesName
real(kind=rp), intent(in), DIMENSION(8):: Y_R
real(kind=rp), intent(in), DIMENSION(8):: Y_PHI
real(kind=rp), intent(in), DIMENSION(8):: Y_Z
integer(kind=is), intent(inout), DIMENSION(8):: flag

Flag that determines whether particles are followed in the simulation (flag=1), or not (flag=0).


Contents


Source Code

  subroutine check_if_in_profiles_domain_p(Y_R,Y_PHI,Y_Z,flag)
    REAL(rp), DIMENSION(8),  INTENT(IN)      :: Y_R,Y_PHI,Y_Z

    INTEGER(is), DIMENSION(8), INTENT(INOUT)  :: flag
    !! Flag that determines whether particles are followed
    !! in the simulation (flag=1), or not (flag=0).
    INTEGER                                                :: IR
    !! @param IR Variable used to localize the grid cell in
    !! the \((R,\phi,Z)\) or \((R,Z)\) grid containing the fields data that
    !! corresponds to the radial position of the particles.
    INTEGER                                                :: IPHI
    !! @param IPHI Variable used to localize the grid cell in
    !! the \((R,\phi,Z)\) or \((R,Z)\) grid containing the fields data that
    !! corresponds to the azimuthal position of the particles.
    INTEGER                                                :: IZ
    !! @param IZ Variable used to localize the grid cell in the
    !! \((R,\phi,Z)\) or \((R,Z)\) grid containing the fields data that
    !! corresponds to the vertical position of the particles.
    INTEGER(ip)                                            :: pp
    !! @param pp Particle iterator.
    INTEGER(ip)                                            :: ss
    !! @param ss Species iterator.
    

    if (ALLOCATED(profiles_domain%FLAG3D)) then
       !$OMP SIMD
!       !$OMP& aligned(IR,IPHI,IZ)
       do pp=1_idef,8_idef
          IR = INT(FLOOR((Y_R(pp)  - profiles_domain%Ro + &
               0.5_rp*profiles_domain%DR)/profiles_domain%DR) + 1.0_rp,idef)
          IPHI = INT(FLOOR((Y_PHI(pp)  + 0.5_rp*profiles_domain%DPHI)/ &
               profiles_domain%DPHI) + 1.0_rp,idef)
          IZ = INT(FLOOR((Y_Z(pp)  + ABS(profiles_domain%Zo) + &
               0.5_rp*profiles_domain%DZ)/profiles_domain%DZ) + 1.0_rp,idef)

          if ((profiles_domain%FLAG3D(IR,IPHI,IZ).NE.1_is).OR. &
               ((IR.GT.profiles_3d%NR).OR.(IZ.GT.profiles_3d%NZ))) then
             flag(pp) = 0_is
          end if
       end do
       !$OMP END SIMD
    else
       !$OMP SIMD
!       !$OMP& aligned(IR,IZ)
       do pp=1_idef,8_idef
          IR = INT(FLOOR((Y_R(pp)  - profiles_domain%Ro + &
               0.5_rp*profiles_domain%DR)/profiles_domain%DR) + 1.0_rp,idef)
          IZ = INT(FLOOR((Y_Z(pp)  + ABS(profiles_domain%Zo) + &
               0.5_rp*profiles_domain%DZ)/profiles_domain%DZ) + 1.0_rp,idef)

          if ((profiles_domain%FLAG2D(IR,IZ).NE.1_is).OR. &
               ((IR.GT.profiles_2d%NR).OR.(IZ.GT.profiles_2d%NZ))) then
             flag(pp) = 0_is

!             write(6,'("Shit''s fucked.")')
          end if
       end do      
       !$OMP END SIMD
!       write(6,'("Shit''s not fucked.")')
    end if
  end subroutine check_if_in_profiles_domain_p