https://www.cdslab.org/paramonte/fortran/2
Current view: top level - test - test_pm_distExp@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 95 95 100.0 %
Date: 2024-04-08 03:18:57 Functions: 32 32 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       2             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       3             : !!!!                                                                                                                            !!!!
       4             : !!!!    ParaMonte: Parallel Monte Carlo and Machine Learning Library.                                                           !!!!
       5             : !!!!                                                                                                                            !!!!
       6             : !!!!    Copyright (C) 2012-present, The Computational Data Science Lab                                                          !!!!
       7             : !!!!                                                                                                                            !!!!
       8             : !!!!    This file is part of the ParaMonte library.                                                                             !!!!
       9             : !!!!                                                                                                                            !!!!
      10             : !!!!    LICENSE                                                                                                                 !!!!
      11             : !!!!                                                                                                                            !!!!
      12             : !!!!       https://github.com/cdslaborg/paramonte/blob/main/LICENSE.md                                                          !!!!
      13             : !!!!                                                                                                                            !!!!
      14             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      15             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      16             : 
      17             : !>  \brief
      18             : !>  This include file contains procedure implementations of the tests of [pm_distExp](@ref pm_distExp).
      19             : !>
      20             : !>  \fintest
      21             : !>
      22             : !>  \author
      23             : !>  \FatemehBagheri, 12:27 AM Tuesday, February 22, 2022, Dallas, TX
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      28             : #if     getExpLogPDF_ENABLED || setExpLogPDF_ENABLED
      29             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      30             : 
      31             : #if     RK_ENABLED
      32             :         real(RKC)   , parameter :: TOL = epsilon(0._RKC) * 10
      33             :         real(RKC)   , parameter :: ZERO = 0._RKC, ONE = 1._RKC
      34             :         real(RKC)               :: x, mu, invSigma, logInvSigma
      35             :         real(RKC)               :: logPDF, logPDF_ref
      36             : #else
      37             : #error  "Unrecognized interface."
      38             : #endif
      39             :         integer(IK) :: itry
      40           8 :         assertion = .true._LK
      41         808 :         do itry = 1, 100
      42        2400 :             mu = getChoice([ZERO, getUnifRand(-10._RKC, 10._RKC)])
      43        2400 :             invSigma = getChoice([ONE, getUnifRand(0.1_RKC, 10._RKC)])
      44         800 :             logInvSigma = log(invSigma)
      45         800 :             x = getUnifRand(mu, mu + 10._RKC)
      46         800 :             logPDF_ref = getLogPDF_ref(x, mu, invSigma)
      47             : #if         getExpLogPDF_ENABLED
      48         800 :             logPDF = getExpLogPDF(x, mu, invSigma)
      49         400 :             call report(__LINE__, mu, invSigma)
      50         400 :             if (mu == ZERO) then
      51         200 :                 logPDF = getExpLogPDF(x, invSigma = invSigma)
      52         200 :                 call report(__LINE__, invSigma = invSigma)
      53             :             end if
      54         400 :             if (invSigma == ONE) then
      55         201 :                 logPDF = getExpLogPDF(x, mu)
      56         201 :                 call report(__LINE__, mu)
      57             :             end if
      58         404 :             if (mu == ZERO .and. invSigma == ONE) then
      59         107 :                 logPDF = getExpLogPDF(x)
      60         107 :                 call report(__LINE__)
      61             :             end if
      62             : #elif       setExpLogPDF_ENABLED
      63         400 :             call setExpLogPDF(logPDF, x, mu, invSigma, logInvSigma)
      64         400 :             call report(__LINE__, mu, invSigma)
      65         404 :             if (mu == ZERO) then
      66         174 :                 call setExpLogPDF(logPDF, x, invSigma, logInvSigma)
      67         174 :                 call report(__LINE__, invSigma = invSigma)
      68         174 :                 if (invSigma == ONE) then
      69          86 :                     call setExpLogPDF(logPDF, x)
      70          86 :                     call report(__LINE__, mu)
      71             :                 end if
      72             :             end if
      73             : #else
      74             : #error      "Unrecognized interface."
      75             : #endif
      76             :         end do
      77             : 
      78             :     contains
      79             : 
      80             :         pure elemental function getLogPDF_ref(x, mu, invSigma) result(logPDF)
      81             :             real(RKC), intent(in) :: x, mu, invSigma
      82             :             real(RKC) :: logPDF
      83         800 :             logPDF = logInvSigma - (x - mu) * invSigma
      84             :         end function
      85             : 
      86        1568 :         subroutine report(line, mu, invSigma)
      87             :             integer(IK) , intent(in)            :: line
      88             :             real(RKC)   , intent(in), optional  :: mu, invSigma
      89        1568 :             assertion = assertion .and. isClose(logPDF, logPDF_ref, abstol = TOL)
      90        1568 :             call test%assert(assertion, SK_"The PDF must be computed correctly.", int(line, IK))
      91        1568 :             if (test%traceable .and. .not. assertion) then
      92             :                 ! LCOV_EXCL_START
      93             :                 call test%disp%skip()
      94             :                 call test%disp%show("x")
      95             :                 call test%disp%show( x )
      96             :                 call test%disp%show("present(mu)")
      97             :                 call test%disp%show( present(mu) )
      98             :                 if (present(mu)) then
      99             :                 call test%disp%show("mu")
     100             :                 call test%disp%show( mu )
     101             :                 end if
     102             :                 call test%disp%show("present(invSigma)")
     103             :                 call test%disp%show( present(invSigma) )
     104             :                 if (present(invSigma)) then
     105             :                 call test%disp%show("logInvSigma")
     106             :                 call test%disp%show( logInvSigma )
     107             :                 call test%disp%show("invSigma")
     108             :                 call test%disp%show( invSigma )
     109             :                 end if
     110             :                 call test%disp%show("logPDF_ref")
     111             :                 call test%disp%show( logPDF_ref )
     112             :                 call test%disp%show("logPDF")
     113             :                 call test%disp%show( logPDF )
     114             :                 call test%disp%skip()
     115             :                 ! LCOV_EXCL_STOP
     116             :             end if
     117        1568 :         end subroutine
     118             : 
     119             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     120             : #elif   getExpCDF_ENABLED || setExpCDF_ENABLED
     121             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     122             : 
     123             : #if     RK_ENABLED
     124             :         real(RKC)   , parameter :: TOL = epsilon(0._RKC) * 10
     125             :         real(RKC)   , parameter :: ZERO = 0._RKC, ONE = 1._RKC
     126             :         real(RKC)               :: invSigma, mu, x
     127             :         real(RKC)               :: cdf, cdf_ref
     128             : #else
     129             : #error  "Unrecognized interface."
     130             : #endif
     131             :         integer(IK) :: itry
     132           8 :         assertion = .true._LK
     133         808 :         do itry = 1, 100
     134        2400 :             mu = getChoice([ZERO, getUnifRand(-10._RKC, 10._RKC)])
     135        2400 :             invSigma = getChoice([ONE, getUnifRand(1._RKC, 10._RKC)])
     136         800 :             x = getUnifRand(mu, mu + 10._RKC)
     137         800 :             cdf_ref = getCDF_ref(x, mu, invSigma)
     138             : #if         getExpCDF_ENABLED
     139         800 :             cdf = getExpCDF(x, mu, invSigma)
     140         400 :             call report(__LINE__, mu, invSigma)
     141         400 :             if (mu == ZERO) then
     142         201 :                 cdf = getExpCDF(x, invSigma = invSigma)
     143         201 :                 call report(__LINE__, invSigma = invSigma)
     144             :             end if
     145         400 :             if (invSigma == ONE) then
     146         211 :                 cdf = getExpCDF(x, mu)
     147         211 :                 call report(__LINE__, mu)
     148             :             end if
     149         404 :             if (mu == ZERO .and. invSigma == ONE) then
     150         115 :                 cdf = getExpCDF(x)
     151         115 :                 call report(__LINE__)
     152             :             end if
     153             : #elif       setExpCDF_ENABLED
     154         400 :             call setExpCDF(cdf, x, mu, invSigma)
     155         400 :             call report(__LINE__, mu, invSigma)
     156         404 :             if (mu == ZERO) then
     157         195 :                 call setExpCDF(cdf, x, invSigma)
     158         195 :                 call report(__LINE__, invSigma = invSigma)
     159         195 :                 if (invSigma == ONE) then
     160         107 :                     call setExpCDF(cdf, x)
     161         107 :                     call report(__LINE__, mu)
     162             :                 end if
     163             :             end if
     164             : #else
     165             : #error      "Unrecognized interface."
     166             : #endif
     167             :         end do
     168             : 
     169             :     contains
     170             : 
     171             :         pure elemental function getCDF_ref(x, mu, invSigma) result(cdf)
     172             :             real(RKC), intent(in) :: x, mu, invSigma
     173             :             real(RKC) :: cdf
     174         800 :             cdf = 1._RKC - exp(-(x - mu) * invSigma)
     175             :         end function
     176             : 
     177        1629 :         subroutine report(line, mu, invSigma)
     178             :             integer(IK) , intent(in)            :: line
     179             :             real(RKC)   , intent(in), optional  :: mu, invSigma
     180        1629 :             assertion = assertion .and. isClose(cdf, cdf_ref, abstol = TOL)
     181        1629 :             call test%assert(assertion, SK_"The CDF must be computed correctly.", int(line, IK))
     182        1629 :             if (test%traceable .and. .not. assertion) then
     183             :                 ! LCOV_EXCL_START
     184             :                 call test%disp%skip()
     185             :                 call test%disp%show("x")
     186             :                 call test%disp%show( x )
     187             :                 call test%disp%show("present(mu)")
     188             :                 call test%disp%show( present(mu) )
     189             :                 if (present(mu)) then
     190             :                 call test%disp%show("mu")
     191             :                 call test%disp%show( mu )
     192             :                 end if
     193             :                 call test%disp%show("present(invSigma)")
     194             :                 call test%disp%show( present(invSigma) )
     195             :                 if (present(invSigma)) then
     196             :                 call test%disp%show("invSigma")
     197             :                 call test%disp%show( invSigma )
     198             :                 end if
     199             :                 call test%disp%show("cdf_ref")
     200             :                 call test%disp%show( cdf_ref )
     201             :                 call test%disp%show("cdf")
     202             :                 call test%disp%show( cdf )
     203             :                 call test%disp%skip()
     204             :                 ! LCOV_EXCL_STOP
     205             :             end if
     206        1629 :         end subroutine
     207             : 
     208             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     209             : #elif   getExpRand_ENABLED || setExpRand_ENABLED
     210             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     211             : 
     212             :         use pm_kind, only: IK
     213             :         use pm_val2str, only: getStr
     214             :         use pm_option, only: getOption
     215             :         use pm_sampleMean, only: getMean
     216             :         use pm_distUnif, only: getUnifRand
     217             :         use pm_distUnif, only: setUnifRand
     218             : 
     219             :         integer(IK) :: i
     220             :         integer(IK) , parameter :: NP = 10000_IK
     221             :         real(RKC)   , parameter :: TOL = 0.1_RKC
     222             :         real(RKC)   , parameter :: SIGMA_DEF = 1._RKC, MU_DEF = 0._RKC
     223             :         real(RKC)               :: rand(NP), diff, reltol
     224             : 
     225           8 :         diff = huge(0._RKC)
     226           8 :         reltol = -huge(0._RKC)
     227           4 :         assertion = .true._LK
     228             : 
     229             : #if     setExpRand_ENABLED
     230           4 :         call runTestsWith()
     231             : #endif
     232          48 :         do i = 1, 5
     233          40 :             call runTestsWith(sigma = getUnifRand(.5_RKC, 5._RKC))
     234          48 :             call runTestsWith(sigma = getUnifRand(.5_RKC, 5._RKC), mu = getUnifRand(-.5_RKC, 5._RKC))
     235             :         end do
     236             : 
     237             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     238             : 
     239             :     contains
     240             : 
     241             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     242             : 
     243          84 :         subroutine runTestsWith(sigma, mu)
     244             :             real(RKC), intent(in), optional :: sigma, mu
     245             :             real(RKC) :: sigmac
     246             : 
     247             :             !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     248             : 
     249             : #if         getExpRand_ENABLED
     250             :             integer :: i
     251          40 :             if (present(mu)) then
     252      200020 :                 do i = 1, size(rand)
     253      200020 :                     rand(i) = getExpRand(sigma, mu)
     254             :                 end do
     255             :             else
     256      200020 :                 do i = 1, size(rand)
     257      200020 :                     rand(i) = getExpRand(sigma)
     258             :                 end do
     259             :             end if
     260             : #elif       setExpRand_ENABLED
     261      440044 :             call setUnifRand(rand)
     262          44 :             if (present(sigma) .and. present(mu)) then
     263      200020 :                 call setExpRand(rand, sigma, mu)
     264          24 :             elseif (present(sigma)) then
     265      200020 :                 call setExpRand(rand, sigma)
     266             :             else
     267       40004 :                 call setExpRand(rand)
     268             :             end if
     269             : #else
     270             : #error      "Unrecognized interface."
     271             : #endif
     272      840084 :             assertion = assertion .and. all(getOption(MU_DEF, mu) <= rand)
     273          84 :             call report(sigma, mu)
     274         252 :             call test%assert(assertion, SK_"The condition `all(mu_ref <= rand)` must hold. Repeating the test may resolve the failure. present(sigma), present(mu) = "//getStr([present(sigma), present(mu)]), int(__LINE__, IK))
     275             : 
     276          84 :             sigmac = getMean(rand) - getOption(MU_DEF, mu)
     277          84 :             diff = abs(sigmac - getOption(SIGMA_DEF, sigma))
     278          84 :             reltol = getOption(SIGMA_DEF, sigma) * TOL
     279          84 :             assertion = assertion .and. diff < reltol
     280          84 :             call report(sigma, mu, sigmac)
     281         252 :             call test%assert(assertion, SK_"The condition `abs(sigmac - getOption(SIGMA_DEF, sigma)) < getOption(SIGMA_DEF, sigma) * TOL` must hold. Repeating the test may resolve the failure. present(sigma), present(mu) = "//getStr([present(sigma), present(mu)]), int(__LINE__, IK))
     282             : 
     283             :             !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     284             : 
     285          84 :         end subroutine
     286             : 
     287             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     288             : 
     289         168 :         subroutine report(sigma, mu, sigmac)
     290             :             integer :: i
     291             :             real(RKC), intent(in), optional :: sigma, mu, sigmac
     292         168 :             if (test%traceable .and. .not. assertion) then
     293             :                 ! LCOV_EXCL_START
     294             :                 call test%disp%skip()
     295             :                 call test%disp%show("present(mu)")
     296             :                 call test%disp%show( present(mu) )
     297             :                 if (present(mu)) then
     298             :                     call test%disp%show("mu")
     299             :                     call test%disp%show( mu )
     300             :                 end if
     301             :                 call test%disp%show("present(sigma)")
     302             :                 call test%disp%show( present(sigma) )
     303             :                 if (present(sigma)) then
     304             :                     call test%disp%show("sigma")
     305             :                     call test%disp%show( sigma )
     306             :                 end if
     307             :                 call test%disp%show("present(sigmac)")
     308             :                 call test%disp%show( present(sigmac) )
     309             :                 if (present(sigmac)) then
     310             :                     call test%disp%show("sigmac")
     311             :                     call test%disp%show( sigmac )
     312             :                 end if
     313             :                 call test%disp%show("reltol")
     314             :                 call test%disp%show( reltol )
     315             :                 call test%disp%show("diff")
     316             :                 call test%disp%show( diff )
     317             :                 do i = 1, NP
     318             :                     if (rand(i) < getOption(MU_DEF, mu)) then
     319             :                         call test%disp%show("[real(i, RKC), getOption(MU_DEF, mu), rand(i)]")
     320             :                         call test%disp%show( [real(i, RKC), getOption(MU_DEF, mu), rand(i)] )
     321             :                         exit
     322             :                     end if
     323             :                 end do
     324             :                 call test%disp%skip()
     325             :                 ! LCOV_EXCL_STOP
     326             :             end if
     327         168 :         end subroutine
     328             : 
     329             : #else
     330             :             !%%%%%%%%%%%%%%%%%%%%%%%%
     331             : #error      "Unrecognized interface."
     332             :             !%%%%%%%%%%%%%%%%%%%%%%%%
     333             : #endif

ParaMonte: Parallel Monte Carlo and Machine Learning Library 
The Computational Data Science Lab
© Copyright 2012 - 2024