multiple_anisotropy Subroutine

public subroutine multiple_anisotropy(tabin, long, larg, scale_xy, multi_fft, vec_ani)

Note

Function that returns simple_anisotropy min, max and max/min for different Gaussian filter cutoff

Arguments

Type IntentOptional Attributes Name
real(kind=R8), intent(in), dimension(1:long, 1:larg) :: tabin

surface acf array

integer(kind=I4), intent(in) :: long

surface length

integer(kind=I4), intent(in) :: larg

surface width

real(kind=R8), intent(in), dimension(1:2) :: scale_xy

lag along x and y

logical(kind=I4), intent(in) :: multi_fft

parallel ffts?

real(kind=R8), intent(out), dimension(1:9) :: vec_ani

anisotropy parameters


Calls

proc~~multiple_anisotropy~~CallsGraph proc~multiple_anisotropy multiple_anisotropy init_scal init_scal proc~multiple_anisotropy->init_scal proc~fft_filter fft_filter proc~multiple_anisotropy->proc~fft_filter proc~simple_anisotropy simple_anisotropy proc~multiple_anisotropy->proc~simple_anisotropy write_surf write_surf proc~multiple_anisotropy->write_surf calc_fftw3_real_bwd calc_fftw3_real_bwd proc~fft_filter->calc_fftw3_real_bwd calc_fftw3_real_fwd calc_fftw3_real_fwd proc~fft_filter->calc_fftw3_real_fwd extend extend proc~fft_filter->extend proc~gaussian_filter gaussian_filter proc~fft_filter->proc~gaussian_filter tab_calc_fftw3_real_bwd tab_calc_fftw3_real_bwd proc~fft_filter->tab_calc_fftw3_real_bwd tab_calc_fftw3_real_fwd tab_calc_fftw3_real_fwd proc~fft_filter->tab_calc_fftw3_real_fwd

Called by

proc~~multiple_anisotropy~~CalledByGraph proc~multiple_anisotropy multiple_anisotropy program~test_anisotropy test_anisotropy program~test_anisotropy->proc~multiple_anisotropy

Source Code

   subroutine multiple_anisotropy(tabin, long, larg, scale_xy, multi_fft, vec_ani)
   !================================================================================================
   !< @note Function that returns simple_anisotropy min, max and max/min for different Gaussian filter
   !< cutoff
   !< @endnote
   !------------------------------------------------------------------------------------------------
   implicit none
   integer(kind=I4), intent(in )                                  :: long       !! *surface length*
   integer(kind=I4), intent(in )                                  :: larg       !! *surface width*
   logical(kind=I4), intent(in )                                  :: multi_fft  !! *parallel ffts?*
   real   (kind=R8), intent(in ), dimension(1:long, 1:larg)       :: tabin      !! *surface acf array*
   real   (kind=R8), intent(in ), dimension(1:2)                  :: scale_xy   !! *lag along x and y*
   real   (kind=R8), intent(out), dimension(1:9)                  :: vec_ani    !! *anisotropy parameters*

      integer(kind=I4) :: icut, i
      real   (kind=R8) :: dx, dy, fft_cutoff,length

      real   (kind=R8) :: max_mean_pks, min_mean_pks, rat_mean_pks, max_mean_len
      real   (kind=R8) :: min_mean_len, rat_mean_len, max_mean_slp, min_mean_slp, rat_mean_slp

      character(len=2) :: str

      real   (kind=R8), dimension(0:179) :: mean_len, mean_pks, mean_slp

      real   (kind=R8), allocatable, dimension(:,:) :: bf_tab

      real   (kind=R8), allocatable, dimension(:,:) :: mat_len    !! *array containing anisotropy outputs*
      real   (kind=R8), allocatable, dimension(:,:) :: mat_pks    !! *array containing anisotropy outputs*
      real   (kind=R8), allocatable, dimension(:,:) :: mat_slp    !! *array containing anisotropy outputs*

      type(SCALE_SURF) :: scal_surf

      allocate( bf_tab(1:long, 1:larg) )

      dx = scale_xy(1) ! x lag
      dy = scale_xy(2) ! y lag

      length = dx / 1.e-6 ! base cutoff = 1 µm

      allocate( mat_len(0:179, 1:10) )
      allocate( mat_pks(0:179, 1:10) )
      allocate( mat_slp(0:179, 1:10) )

      call init_scal( scal   = scal_surf,       &  !
                      nx     = long,            &  !
                      ny     = larg,            &  !
                      lx     = long * dx,       &  !
                      ly     = larg * dy,       &  !
                      unit_z = 'nm')               !

      do icut = 1, 10

         fft_cutoff = length / icut

         call fft_filter(tab       = tabin(1:long, 1:larg),    &  ! in
                         long      = long,                     &  ! in
                         larg      = larg,                     &  ! in
                         cutoff    = fft_cutoff,               &  ! in
                         bf_tab    = bf_tab(1:long, 1:larg),   &  ! out
                         multi_fft = multi_fft)                   ! in

         write(str,'(i2.2)') icut

         if (.false.) then
            call write_surf( nom_fic = "out/bf_tab"//str//".sur",    &  !
                             tab_s   = bf_tab(1:long, 1:larg),       &  !
                             scal    = scal_surf )                      !
         endif

         call simple_anisotropy( tabin     = bf_tab(1:long, 1:larg),     &  ! in
                                 long      = long,                       &  ! in
                                 larg      = larg,                       &  ! in
                                 vec_len   = mat_len(0:179, icut),       &  ! out   path length ratios
                                 vec_pks   = mat_pks(0:179, icut),       &  ! out   peak mean width
                                 vec_slp   = mat_slp(0:179, icut),       &  ! out   RMS slopes
                                 scale_xy  = scale_xy(1:2),              &  ! in
                                 multi_fft = multi_fft)                     ! in

      enddo

      do i = 0, 179

        mean_pks(i) = sum( mat_pks(i, 06:10) ) / 5. ! mean on smoother profiles
        mean_len(i) = sum( mat_len(i, 01:05) ) / 5.
        mean_slp(i) = sum( mat_slp(i, 01:05) ) / 5.

      enddo

      deallocate( bf_tab, mat_len, mat_pks, mat_slp )

      max_mean_pks = maxval( mean_pks(0:179) )     ! robust maximum of the peak mean width
      min_mean_pks = minval( mean_pks(0:179) )     !        minimum
      rat_mean_pks = max_mean_pks / min_mean_pks   ! ratio

      max_mean_len = maxval( mean_len(0:179) )     ! robust maximum of the path length ratio
      min_mean_len = minval( mean_len(0:179) )     !        minimum
      rat_mean_len = max_mean_len / min_mean_len   ! ratio

      max_mean_slp = maxval( mean_slp(0:179) )     ! robust maximum of the RMS slopes
      min_mean_slp = minval( mean_slp(0:179) )     !        minimum
      rat_mean_slp = max_mean_slp / min_mean_slp   ! ratio

      vec_ani(1:9) = [ max_mean_pks, min_mean_pks, rat_mean_pks,  &  !
                       max_mean_len, min_mean_len, rat_mean_len,  &  !
                       max_mean_slp, min_mean_slp, rat_mean_slp ]    !

   return
   endsubroutine multiple_anisotropy