stdlib_stats_moment.fypp Source File


Source Code

#:include "common.fypp"
#:set RANKS = range(1, MAXRANK + 1)
#:set REDRANKS = range(2, MAXRANK + 1)
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
submodule (stdlib_stats) stdlib_stats_moment

  use, intrinsic:: ieee_arithmetic, only: ieee_value, ieee_quiet_nan
  use stdlib_error, only: error_stop
  use stdlib_optval, only: optval
  implicit none

contains

  #:for k1, t1 in RC_KINDS_TYPES
    #:for rank in RANKS
      #:set RName = rname("moment",rank, t1, k1)
      module function ${RName}$(x, order, dim, center, mask) result(res)
        ${t1}$, intent(in) :: x${ranksuffix(rank)}$
        integer, intent(in) :: order
        integer, intent(in) :: dim
        ${t1}$, intent(in), optional :: center${reduced_shape('x', rank, 'dim')}$
        logical, intent(in), optional :: mask
        ${t1}$ :: res${reduced_shape('x', rank, 'dim')}$

        integer :: i
        real(${k1}$) :: n
        ${t1}$, allocatable :: mean_${ranksuffix(rank-1)}$

        if (.not.optval(mask, .true.)) then
          res = ieee_value(1._${k1}$, ieee_quiet_nan)
          return
        end if

        n = real(size(x, dim), ${k1}$)

        res = 0
        select case(dim)
          #:for fi in range(1, rank+1)
          case(${fi}$)
            if (present(center)) then
              do i = 1, size(x, ${fi}$)
                res = res + (x${select_subarray(rank, [(fi, 'i')])}$ - center)**order
              end do
            else
              allocate(mean_, source = mean(x, ${fi}$))
              do i = 1, size(x, ${fi}$)
                res = res + (x${select_subarray(rank, [(fi, 'i')])}$ - mean_)**order
              end do
              deallocate(mean_)
            end if
          #:endfor
          case default
            call error_stop("ERROR (moment): wrong dimension")
        end select
        res = res / n

      end function ${RName}$
    #:endfor
  #:endfor


  #:for k1, t1 in INT_KINDS_TYPES
    #:for rank in RANKS
      #:set RName = rname("moment",rank, t1, k1, 'dp')
      module function ${RName}$(x, order, dim, center, mask) result(res)
        ${t1}$, intent(in) :: x${ranksuffix(rank)}$
        integer, intent(in) :: order
        integer, intent(in) :: dim
        real(dp),intent(in), optional :: center${reduced_shape('x', rank, 'dim')}$
        logical, intent(in), optional :: mask
        real(dp) :: res${reduced_shape('x', rank, 'dim')}$

        integer :: i
        real(dp) :: n
        real(dp), allocatable :: mean_${ranksuffix(rank-1)}$

        if (.not.optval(mask, .true.)) then
          res = ieee_value(1._dp, ieee_quiet_nan)
          return
        end if

        n = real(size(x, dim), dp)

        res = 0
        select case(dim)
          #:for fi in range(1, rank+1)
          case(${fi}$)
            if (present(center)) then
              do i = 1, size(x, ${fi}$)
                res = res + (real(x${select_subarray(rank, [(fi, 'i')])}$, dp) -&
                  center)**order
              end do
            else
              allocate(mean_, source = mean(x, ${fi}$))
              do i = 1, size(x, ${fi}$)
                res = res + (real(x${select_subarray(rank, [(fi, 'i')])}$, dp) - mean_)**order
              end do
              deallocate(mean_)
            end if
          #:endfor
          case default
            call error_stop("ERROR (moment): wrong dimension")
        end select
        res = res / n

      end function ${RName}$
    #:endfor
  #:endfor

end submodule