morpho_filter Subroutine

public subroutine morpho_filter(tabin, tabou, long, larg, scale_xyz, ray, omp, nb_div, mtype)

Note

Morphological filter: uses combinations of roll_smooth to provide all kind of transformation :

  • closing
  • opening
  • dilation
  • erosion

Arguments

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

2D array in

real(kind=R8), intent(out), dimension(1:long, 1:larg) :: tabou

2D array out

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

2D array length

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

2D array width

real(kind=R8), intent(in), dimension(1:3) :: scale_xyz

lag along x, y and scale z

real(kind=R8), intent(in) :: ray

roll radius

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

if multithreading

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

number of macro elements along an axis

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

closing, opening, dilation or erosion


Calls

proc~~morpho_filter~~CallsGraph proc~morpho_filter morpho_filter proc~roll_smooth roll_smooth proc~morpho_filter->proc~roll_smooth omp_get_num_procs omp_get_num_procs proc~roll_smooth->omp_get_num_procs

Called by

proc~~morpho_filter~~CalledByGraph proc~morpho_filter morpho_filter program~test_smooth test_smooth program~test_smooth->proc~morpho_filter

Source Code

   subroutine morpho_filter(tabin, tabou, long, larg, scale_xyz, ray, omp, nb_div, mtype)
   !================================================================================================
   !< @note
   !<
   !< Morphological filter: uses combinations of [[roll_smooth]] to provide all kind of transformation :
   !<
   !< + closing
   !< + opening
   !< + dilation
   !< + erosion
   !<
   !< @endnote
   !------------------------------------------------------------------------------------------------
   implicit none
   integer(kind=I4), intent(in )                            :: long        !! *2D array length*
   integer(kind=I4), intent(in )                            :: larg        !! *2D array width*
   integer(kind=I4), intent(in )                            :: nb_div      !! *number of macro elements along an axis*
   logical(kind=I4), intent(in )                            :: omp         !! *if multithreading*
   real   (kind=R8), intent(in )                            :: ray         !! *roll radius*
   character(len=*), intent(in )                            :: mtype       !! *closing, opening, dilation or erosion*
   real   (kind=R8), intent(in ), dimension(1:3)            :: scale_xyz   !! *lag along x, y and scale z*
   real   (kind=R8), intent(in ), dimension(1:long, 1:larg) :: tabin       !! *2D array in*
   real   (kind=R8), intent(out), dimension(1:long, 1:larg) :: tabou       !! *2D array out*

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

      logical(kind=I4) :: op1, op2
      integer(kind=I4) :: k

      select case ( mtype(1:7) )

         case("closing")

            op1 = .true. ; op2 = .true. ; k = +1

         case("opening")

            op1 = .true. ; op2 = .true. ; k = -1

         case("dilatio")

            op1 = .true. ; op2 = .false. ; k = +1

         case("erosion")

            op1 = .true. ; op2 = .false. ; k = -1

         case default

            stop "bad choice, morpho_filter"

      endselect

      if (op1) then

         call roll_smooth( tabin       = tabin(1:long, 1:larg),      &  ! IN
                           tabou       = tabou(1:long, 1:larg),      &  ! OUT
                           long        = long,                       &  ! IN
                           larg        = larg,                       &  ! IN
                           scale_xyz   = scale_xyz,                  &  ! IN
                           sgn         = +k,                         &  ! IN
                           ray         = ray,                        &  ! IN
                           omp         = omp,                        &  ! IN
                           nb_div      = nb_div )                       ! IN

      endif

      if (op2) then

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

         call roll_smooth( tabin       = tabou(1:long, 1:larg),      &  ! IN
                           tabou       = tab_tmp(1:long, 1:larg),    &  ! OUT
                           long        = long,                       &  ! IN
                           larg        = larg,                       &  ! IN
                           scale_xyz   = scale_xyz,                  &  ! IN
                           sgn         = -k,                         &  ! IN
                           ray         = ray,                        &  ! IN
                           omp         = omp,                        &  ! IN
                           nb_div      = nb_div )                       ! IN

         tabou(1:long, 1:larg) = tab_tmp(1:long, 1:larg)

         deallocate( tab_tmp )

      endif

   return
   endsubroutine morpho_filter