https://www.cdslab.org/paramonte/fortran/2
Current view: top level - test - test_pm_mathCompare@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 92 92 100.0 %
Date: 2024-04-08 03:18:57 Functions: 24 24 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_mathCompare](@ref pm_mathCompare).
      19             : !>
      20             : !>  \fintest
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, Tuesday 2:06 AM, September 21, 2021, Dallas, TX
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         use pm_distUnif, only: getUnifRand, setUnifRand
      28             :         use pm_except, only: setInfNeg, setInfPos, isInf
      29             :         use pm_except, only: setNAN, isNAN
      30             : #if     isClose_CK_ENABLED
      31             :         complex(CKC)                    :: x, y
      32             :         complex(CKC)    , parameter     :: ZERO = (0._CKC, 0._CKC)
      33             :         complex(CKC)    , parameter     :: low = -cmplx(sqrt(sqrt(huge(0._CKC))), sqrt(sqrt(huge(0._CKC))), CKC)
      34             :         complex(CKC)    , parameter     :: upp = +cmplx(sqrt(sqrt(huge(0._CKC))), sqrt(sqrt(huge(0._CKC))), CKC)
      35             : #elif   isClose_RK_ENABLED
      36             :         real(RKC)                       :: x, y
      37             :         real(RKC)       , parameter     :: low = -sqrt(sqrt(huge(low)))
      38             :         real(RKC)       , parameter     :: upp = +sqrt(sqrt(huge(upp)))
      39             :         real(RKC)       , parameter     :: ZERO = 0._RKC
      40             : #else
      41             : #error  "Unrecognized interface."
      42             : #endif
      43             :         integer         , parameter     :: TKC = kind(x)
      44             :         real(TKC)       , parameter     :: EPS = epsilon(1._TKC) * 100, TIN = tiny(0._TKC)
      45             :         logical(LK)                     :: close, close_def, isRelTolDef, isAbsTolDef
      46             :         real(TKC)                       :: reltol, abstol
      47          16 :         class(*)        , allocatable   :: method
      48             :         integer(IK)                     :: i, choice
      49             : 
      50           8 :         reltol = EPS
      51           8 :         abstol = TIN
      52           8 :         assertion = .true._LK
      53        8016 :         do i = 1, 1000
      54             : 
      55        8000 :             choice = getUnifRand(0_IK, 9_IK)
      56        8000 :             if (choice == 0_IK) then
      57         821 :                 call setNAN(x)
      58        7179 :             elseif (choice == 1_IK) then
      59         805 :                 call setInfNeg(x)
      60        6374 :             elseif (choice == 2_IK) then
      61         760 :                 call setInfPos(x)
      62        5614 :             elseif (choice == 3_IK) then
      63         760 :                 x = ZERO
      64             :             else
      65        4854 :                 call setUnifRand(x, low, upp)
      66             :             end if
      67             : 
      68        8000 :             isAbsTolDef = getUnifRand()
      69        8000 :             isRelTolDef = getUnifRand()
      70        8000 :             if (.not. (isInf(x) .or. isNAN(x))) reltol = merge(EPS, getUnifRand(0.1_TKC * abs(x), 10._TKC * abs(x)), isRelTolDef)
      71        8000 :             if (.not. (isInf(x) .or. isNAN(x))) abstol = merge(TIN, getUnifRand(0.2_TKC * abs(x),  5._TKC * abs(x)), isAbsTolDef)
      72             : 
      73        8000 :             choice = getUnifRand(0_IK, 9_IK)
      74        8000 :             if (choice == 0_IK) then
      75         772 :                 call setNAN(y)
      76        7228 :             elseif (choice == 1_IK) then
      77         819 :                 call setInfNeg(y)
      78        6409 :             elseif (choice == 2_IK) then
      79         783 :                 call setInfPos(y)
      80        5626 :             elseif (choice == 3_IK) then
      81         810 :                 y = ZERO
      82        4816 :             elseif (.not. (isInf(x) .or. isNAN(x))) then
      83        3356 :                 call setUnifRand(y, x - 3 * reltol, x + 3 * reltol)
      84             :             else
      85        1460 :                 call setUnifRand(y, ZERO - 3 * reltol, ZERO + 3 * reltol)
      86             :             end if
      87             : 
      88        8008 :             if (getUnifRand()) then ! method
      89        4087 :                 choice = getUnifRand(1_IK, 4_IK)
      90        4087 :                 if (choice == 1_IK) then
      91        1065 :                     method = reference_type()
      92        3022 :                 elseif (choice == 2_IK) then
      93         990 :                     method = strong_type()
      94        2032 :                 elseif (choice == 3_IK) then
      95        1055 :                     method = weak_type()
      96         977 :                 elseif (choice == 4_IK) then
      97         977 :                     method = mean_type()
      98             :                 end if
      99        4087 :                 if (getUnifRand()) then ! reltol
     100        2056 :                     if (getUnifRand()) then ! abstol
     101        1000 :                         call runTestsWith(method, reltol, abstol)
     102             :                     else ! no abstol
     103        1056 :                         call runTestsWith(method, reltol)
     104             :                     end if
     105             :                 else ! no reltol
     106        2031 :                     if (getUnifRand()) then ! abstol
     107        1003 :                         call runTestsWith(method, abstol = abstol)
     108             :                     else ! no abstol
     109        1028 :                         call runTestsWith(method)
     110             :                     end if
     111             :                 end if
     112             :             else ! no method
     113        3913 :                 if (getUnifRand()) then ! reltol
     114        1940 :                     if (getUnifRand()) then ! abstol
     115         989 :                         call runTestsWith(reltol = reltol, abstol = abstol)
     116             :                     else ! no abstol
     117         951 :                         call runTestsWith(reltol = reltol)
     118             :                     end if
     119             :                 else ! no reltol
     120        1973 :                     if (getUnifRand()) then ! abstol
     121         971 :                         call runTestsWith(abstol = abstol)
     122             :                     else ! no abstol
     123        1002 :                         call runTestsWith()
     124             :                     end if
     125             :                 end if
     126             :             endif
     127             :         end do
     128             : 
     129             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     130             : 
     131             :     contains
     132             : 
     133             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     134             : 
     135        8000 :         subroutine runTestsWith(method, reltol, abstol)
     136             :             class(*)    , intent(in), optional  :: method
     137             :             real(TKC)   , intent(in), optional  :: reltol, abstol
     138        8000 :             if (present(method)) then
     139        4087 :                 close_def = isClose_def(method, reltol, abstol)
     140             :                 select type (method)
     141             :                     type is (reference_type)
     142        1065 :                         close = isClose(x, y, method, reltol, abstol)
     143             :                     type is (strong_type)
     144         990 :                         close = isClose(x, y, method, reltol, abstol)
     145             :                     type is (weak_type)
     146        1055 :                         close = isClose(x, y, method, reltol, abstol)
     147             :                     type is (mean_type)
     148         977 :                         close = isClose(x, y, method, reltol, abstol)
     149             :                 end select
     150             :             else
     151        3913 :                 close_def = isClose_def(reltol = reltol, abstol = abstol)
     152        3913 :                 close = isClose(x, y, reltol, abstol)
     153             :             endif
     154        8000 :             call report(__LINE__, method, reltol, abstol)
     155        8000 :         end subroutine
     156             : 
     157             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     158             : 
     159        8000 :         function isClose_def(method, reltol, abstol) result(close)
     160             :             use pm_except, only: isInf, isInfNeg, isInfPos
     161             :             use pm_except, only: isNAN
     162             :             class(*)    , intent(in), optional  :: method
     163             :             real(TKC)   , intent(in), optional  :: reltol, abstol
     164             :             logical(LK) :: close
     165             : 
     166             :             real(TKC) :: reltol_def, abstol_def
     167             :             real(TKC) :: absDiff
     168             : 
     169        8000 :             if (isNAN(x) .or. isNAN(y)) then
     170             :                 close = .false._LK
     171             :                 return
     172             :             end if
     173             : 
     174        6505 :             if (.not. (isInf(x) .or. isInf(y))) then
     175        3913 :                 close = logical(x == y, LK)
     176        3913 :                 if (close) return
     177        2166 :                 if (present(reltol)) then
     178        1071 :                     reltol_def = reltol
     179             :                 else
     180         542 :                     reltol_def = epsilon(0._TKC)
     181             :                 end if
     182        2166 :                 if (present(abstol)) then
     183        1047 :                     abstol_def = abstol
     184             :                 else
     185         562 :                     abstol_def = tiny(0._TKC)
     186             :                 end if
     187        2166 :                 absDiff = abs(y - x)
     188        2166 :                 if (present(method)) then
     189             :                     select type (method)
     190             :                         type is (reference_type)
     191         291 :                             close = logical(absDiff <= abs(reltol_def * x) .or. absDiff <= abstol_def, LK)
     192             :                         type is (strong_type)
     193         292 :                             close = logical((absDiff <= abs(reltol_def * x) .and. absDiff <= abs(reltol_def * y)) .or. absDiff <= abstol_def, LK)
     194             :                         type is (weak_type)
     195         279 :                             close = logical(absDiff <= abs(reltol_def * x) .or. absDiff <= abs(reltol_def * y) .or. absDiff <= abstol_def, LK)
     196             :                         type is (mean_type)
     197         256 :                             close = logical(absDiff <= abs(reltol_def * 0.5_TKC * (x + y)) .or. absDiff <= abstol_def, LK)
     198             :                     end select
     199             :                 else ! assume `weak_type`.
     200        1048 :                     close = logical(absDiff <= abs(reltol_def * x) .or. absDiff <= abs(reltol_def * y) .or. absDiff <= abstol_def, LK)
     201             :                 end if
     202        2166 :                 return
     203             :             end if
     204        2592 :             close = (isInfNeg(x) .and. isInfNeg(y)) .or. (isInfPos(x) .and. isInfPos(y))
     205             :         end function
     206             : 
     207             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     208             : 
     209        8000 :         subroutine report(line, method, reltol, abstol)
     210             :             use pm_io, only: display_type
     211             :             use pm_val2str, only: getStr
     212             :             use pm_option, only: getOption
     213             :             integer     , intent(in) :: line
     214             :             class(*)    , intent(in), optional :: method
     215             :             real(TKC)   , intent(in), optional :: reltol, abstol
     216        8000 :             type(display_type)  :: disp
     217        8000 :             assertion = assertion .and. (close .eqv. close_def)
     218        8000 :             if (test%traceable .and. .not. assertion) then
     219             :                 ! LCOV_EXCL_START
     220             :                 call disp%skip()
     221             :                 call disp%show("[x, y]")
     222             :                 call disp%show( [x, y] )
     223             :                 call disp%show("[present(reltol), present(abstol)]")
     224             :                 call disp%show( [present(reltol), present(abstol)] )
     225             :                 call disp%show("[getOption(EPS, reltol), getOption(TIN, abstol)]")
     226             :                 call disp%show( [getOption(EPS, reltol), getOption(TIN, abstol)] )
     227             :                 call disp%show("present(method)")
     228             :                 call disp%show( present(method) )
     229             :                 if (present(method)) then
     230             :                     select type (method)
     231             :                         type is (reference_type)
     232             :                             call disp%show("method: reference_type()")
     233             :                         type is (strong_type)
     234             :                             call disp%show("method: strong_type()")
     235             :                         type is (weak_type)
     236             :                             call disp%show("method: weak_type()")
     237             :                         type is (mean_type)
     238             :                             call disp%show("method: mean_type()")
     239             :                     end select
     240             :                 end if
     241             :                 call disp%show("[close, close_def]")
     242             :                 call disp%show( [close, close_def] )
     243             :                 call disp%skip()
     244             :                 ! LCOV_EXCL_STOP
     245             :             end if
     246        8000 :             call test%assert(assertion, SK_"The proximity of the two numbers must be computed correctly.", int(line, IK))
     247        8000 :         end subroutine

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