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_mathCumPropExp](@ref pm_mathCumPropExp).
19 : !>
20 : !> \fintest
21 : !>
22 : !> \author
23 : !> \AmirShahmoradi, Tuesday 2:06 AM, September 21, 2021, Dallas, TX
24 :
25 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
26 :
27 : integer(IK) :: itry
28 : #if getCumPropExp_ENABLED
29 : logical(LK), parameter :: isnew = .true.
30 : #elif setCumPropExp_ENABLED
31 : logical(LK) :: isnew
32 : #else
33 : #error "Unrecognized interface."
34 : #endif
35 : #if Sel_ENABLED
36 : type(selection_type), parameter :: control = selection_type()
37 : #elif Seq_ENABLED
38 : type(sequence_type), parameter :: control = sequence_type()
39 : #else
40 : #error "Unrecognized interface."
41 : #endif
42 : real(RKC), parameter :: TOL = epsilon(1._RKC) * 10._RKC
43 : !real(RKC), parameter :: LB = 1._RKC, UB = 2._RKC
44 : real(RKC), parameter :: LB = minexponent(0._RKC), UB = maxexponent(0._RKC)
45 : real(RKC), allocatable :: cumPropExp_ref(:), cumPropExp(:), array(:), diff(:)
46 : logical(LK) :: isbackward, isreverse
47 :
48 16 : assertion = .true._LK
49 :
50 4816 : do itry = 1, 300
51 :
52 : #if setCumPropExp_ENABLED
53 2400 : isnew = getUnifRand()
54 : #endif
55 4800 : isreverse = getUnifRand()
56 4800 : isbackward = getUnifRand()
57 35774 : array = getUnifRand(LB, UB, getUnifRand(1_IK, 10_IK))
58 4800 : call setResized(cumPropExp, size(array, 1, IK))
59 :
60 4800 : if (isbackward .and. isreverse) then
61 24369 : cumPropExp_ref = getCumPropExp_ref(array, maxval(array), backward, reverse)
62 : #if setCumPropExp_ENABLED
63 658 : if (isnew) then
64 2451 : call setCumPropExp(cumPropExp, array, maxval(array), control, backward, reverse)
65 : else
66 2420 : cumPropExp = array
67 2420 : call setCumPropExp(cumPropExp, array, maxval(array), control, backward, reverse)
68 : end if
69 : #elif getCumPropExp_ENABLED
70 11730 : cumPropExp = getCumPropExp(array, maxval(array), control, backward, reverse)
71 : #else
72 : #error "Unrecognized interface."
73 : #endif
74 3535 : elseif (isbackward) then
75 22725 : cumPropExp_ref = getCumPropExp_ref(array, maxval(array), backward, nothing)
76 : #if setCumPropExp_ENABLED
77 580 : if (isnew) then
78 2308 : call setCumPropExp(cumPropExp, array, maxval(array), control, backward, nothing)
79 : else
80 1986 : cumPropExp = array
81 1986 : call setCumPropExp(cumPropExp, maxval(array), control, backward, nothing)
82 : end if
83 : #elif getCumPropExp_ENABLED
84 11583 : cumPropExp = getCumPropExp(array, maxval(array), control, backward, nothing)
85 : #endif
86 2355 : elseif (isreverse) then
87 22737 : cumPropExp_ref = getCumPropExp_ref(array, maxval(array), forward, reverse)
88 : #if setCumPropExp_ENABLED
89 568 : if (isnew) then
90 2017 : call setCumPropExp(cumPropExp, array, maxval(array), control, forward, reverse)
91 : else
92 2172 : cumPropExp = array
93 2172 : call setCumPropExp(cumPropExp, maxval(array), control, forward, reverse)
94 : end if
95 : #elif getCumPropExp_ENABLED
96 11874 : cumPropExp = getCumPropExp(array, maxval(array), control, forward, reverse)
97 : #endif
98 : else
99 23091 : cumPropExp_ref = getCumPropExp_ref(array, maxval(array), forward, nothing)
100 : #if setCumPropExp_ENABLED
101 594 : if (isnew) then
102 2170 : call setCumPropExp(cumPropExp, array, maxval(array), control, forward, nothing)
103 : else
104 2254 : cumPropExp = array
105 2254 : call setCumPropExp(cumPropExp, maxval(array), control, forward, nothing)
106 : end if
107 : #elif getCumPropExp_ENABLED
108 11601 : cumPropExp = getCumPropExp(array, maxval(array), control, forward, nothing)
109 : #endif
110 : end if
111 4808 : call report(__LINE__)
112 :
113 : #if getCumPropExp_ENABLED
114 2400 : if (isbackward .and. .not. isreverse) call runTestsWith(direction = backward)
115 2400 : if (isreverse .and. .not. isbackward) call runTestsWith(action = reverse)
116 2408 : if (.not. (isbackward .or. isreverse)) call runTestsWith()
117 : #endif
118 :
119 : end do
120 :
121 : contains
122 :
123 : #if getCumPropExp_ENABLED
124 1793 : subroutine runTestsWith(direction, action)
125 : class(action_type), intent(in), optional :: action
126 : class(direction_type), intent(in), optional :: direction
127 35058 : cumPropExp = getCumPropExp(array, maxval(array), control, direction, action)
128 1793 : call report(__LINE__)
129 25165 : cumPropExp = getCumPropExp(array, maxval(array), direction, action)
130 1793 : call report(__LINE__)
131 1793 : end subroutine
132 : #endif
133 :
134 4800 : function getCumPropExp_ref(array, maxArray, direction, action) result(cumPropExp)
135 : class(direction_type), intent(in), optional :: direction
136 : class(action_type), intent(in), optional :: action
137 14400 : class(direction_type), allocatable :: direction_def
138 9600 : class(action_type), allocatable :: action_def
139 : real(RKC), intent(in) :: array(:), maxArray
140 : real(RKC) :: cumPropExp(size(array, 1, IK))
141 4800 : action_def = nothing
142 4800 : direction_def = forward
143 4800 : if (present(action)) action_def = action
144 4800 : if (present(direction)) direction_def = direction
145 4800 : if (same_type_as(direction_def, backward) .and. same_type_as(action_def, reverse)) then
146 8123 : cumPropExp = getCumSum(exp(array - maxArray), backward, reverse)
147 3535 : elseif (same_type_as(direction_def, backward) .and. same_type_as(action_def, nothing)) then
148 7575 : cumPropExp = getCumSum(exp(array - maxArray), backward, nothing)
149 2355 : elseif (same_type_as(direction_def, forward) .and. same_type_as(action_def, reverse)) then
150 7579 : cumPropExp = getCumSum(exp(array - maxArray), forward, reverse)
151 1182 : elseif (same_type_as(direction_def, forward) .and. same_type_as(action_def, nothing)) then
152 7697 : cumPropExp = getCumSum(exp(array - maxArray), forward, nothing)
153 : end if
154 66748 : cumPropExp = cumPropExp / maxval(cumPropExp)
155 4800 : end function
156 :
157 8386 : subroutine report(line)
158 : integer :: line
159 62732 : diff = abs(cumPropExp - cumPropExp_ref)
160 54346 : assertion = all(diff < TOL)
161 8386 : if (test%traceable .and. .not. assertion) then
162 : ! LCOV_EXCL_START
163 : call test%disp%skip
164 : call test%disp%show("isnew")
165 : call test%disp%show( isnew )
166 : call test%disp%show("isreverse")
167 : call test%disp%show( isreverse )
168 : call test%disp%show("isbackward")
169 : call test%disp%show( isbackward )
170 : call test%disp%show("cumPropExp_ref")
171 : call test%disp%show( cumPropExp_ref )
172 : call test%disp%show("cumPropExp")
173 : call test%disp%show( cumPropExp )
174 : call test%disp%show("array")
175 : call test%disp%show( array )
176 : call test%disp%show("diff")
177 : call test%disp%show( diff )
178 : call test%disp%show("TOL")
179 : call test%disp%show( TOL )
180 : call test%disp%skip
181 : ! LCOV_EXCL_STOP
182 : end if
183 8386 : call test%assert(assertion, SK_"The output `cumPropExp` must be correctly computed.", line)
184 8386 : end subroutine
|