File: test_distribution_normal.fypp

package info (click to toggle)
fortran-stdlib 0.8.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 34,008 kB
  • sloc: f90: 24,178; ansic: 1,244; cpp: 623; python: 119; makefile: 13
file content (322 lines) | stat: -rw-r--r-- 14,520 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322

#:include "common.fypp"
#:set RC_KINDS_TYPES = REAL_KINDS_TYPES + CMPLX_KINDS_TYPES
program test_distribution_normal
    use stdlib_kinds, only : sp, dp, xdp, qp
    use stdlib_error, only : check
    use stdlib_random, only : random_seed
    use stdlib_stats_distribution_uniform, only : uni => rvs_uniform
    use stdlib_stats_distribution_normal, only : nor_rvs => rvs_normal,        &
                          nor_pdf => pdf_normal, nor_cdf => cdf_normal

    implicit none
    #:for k1, t1 in REAL_KINDS_TYPES
    ${t1}$, parameter :: ${k1}$tol = 1000 * epsilon(1.0_${k1}$)
    #:endfor
    logical ::  warn = .true.
    integer :: put, get

    put = 12345678
    call random_seed(put, get)

    call test_normal_random_generator


    #:for k1, t1 in RC_KINDS_TYPES
    call test_nor_rvs_${t1[0]}$${k1}$
    #:endfor

    #:for k1, t1 in RC_KINDS_TYPES
    call test_nor_rvs_default_${t1[0]}$${k1}$
    #:endfor



    #:for k1, t1 in RC_KINDS_TYPES
    call test_nor_pdf_${t1[0]}$${k1}$
    #:endfor



    #:for k1, t1 in RC_KINDS_TYPES
    call test_nor_cdf_${t1[0]}$${k1}$
    #:endfor





contains

    subroutine test_normal_random_generator

        integer, parameter :: num = 10000000, array_size = 1000
        integer :: i, j, freq(0:array_size)
        real(dp) :: chisq, expct

        print *, ""
        print *, "Test normal random generator with chi-squared"
        freq = 0
        do i = 1, num
            j = array_size * (1 + erf(nor_rvs(0.0, 1.0) / sqrt(2.0))) / 2.0
            freq(j) = freq(j) + 1
        end do
        chisq = 0.0_dp
        expct = num / array_size
        do i = 0, array_size - 1
           chisq = chisq + (freq(i) - expct) ** 2 / expct
        end do
        write(*,*) "The critical values for chi-squared with 1000 d. of f. is" &
            //" 1143.92"
        write(*,*) "Chi-squared for normal random generator is : ", chisq
        call check((chisq < 1143.9),                                           &
            msg = "normal randomness failed chi-squared test", warn = warn)
    end subroutine test_normal_random_generator




    #:for k1, t1 in RC_KINDS_TYPES
    subroutine test_nor_rvs_${t1[0]}$${k1}$
        ${t1}$ :: res(10), loc, scale
        integer, parameter :: k = 5
        integer :: i
        integer :: seed, get
        #:if t1[0] == "r"
            #! for real type
        ${t1}$, parameter :: ans(10) =                                         &
                            [2.66708039318040679432897377409972250_${k1}$,     &
                             2.36030794936128329730706809641560540_${k1}$,     &
                             1.27712218793084242296487218482070602_${k1}$,     &
                            -2.39132544130814794769435138732660562_${k1}$,     &
                             1.72566595106028652928387145948363468_${k1}$,     &
                            -1.50621775537767632613395107910037041_${k1}$,     &
                             2.13518835158352082714827702147886157_${k1}$,     &
                           -0.636788253742142318358787633769679815_${k1}$,     &
                             2.48600787778845799813609573902795091_${k1}$,     &
                            -3.03711473837981227319460231228731573_${k1}$]
        #:else
            #! for complex type
        ${t1}$, parameter :: ans(10) =                                         &
                            [(2.12531029488530509574673033057479188_${k1}$,    &
                              1.46507698734032082432676702410390135_${k1}$),   &
                             (1.08284164094813181722365413861552952_${k1}$,    &
                             0.277168639672963013076412153168348595_${k1}$),   &
                             (1.41924946329521489696290359461272601_${k1}$,    &
                             0.498445561155580918466512230224907398_${k1}$),   &
                             (1.72639126368764062036120776610914618_${k1}$,    &
                             0.715802936564464420410303091557580046_${k1}$),   &
                             (1.98950590834134349860207180427096318_${k1}$,    &
                             0.115721315405046931701349421928171068_${k1}$),   &
                            (-1.16929014824793620075382705181255005_${k1}$,    &
                             0.250744737486995217246033007540972903_${k1}$),   &
                             (1.57160542831869509683428987045772374_${k1}$,    &
                             0.638282596371312238581197107123443857_${k1}$),   &
                            (-1.36106107654239116833139178197598085_${k1}$,    &
                             0.166259201494369124318950525776017457_${k1}$),   &
                             (1.13403096805387920698038328737311531_${k1}$,    &
                              1.04232618148691447146347854868508875_${k1}$),   &
                            (-1.68220535920475811053620418533682823_${k1}$,    &
                              1.63361446685040256898702182297711261_${k1}$)]
        #:endif

        print *, "Test normal_distribution_rvs_${t1[0]}$${k1}$"
        seed = 25836914
        call random_seed(seed, get)
        #:if t1[0] == "r"
            #! for real type
        loc = 0.5_${k1}$; scale = 2.0_${k1}$
        #:else
            #! for complex type
        loc = (0.5_${k1}$, 1.0_${k1}$); scale = (1.5_${k1}$, 0.5_${k1}$)
        #:endif

        do i = 1, k
           res(i) = nor_rvs(loc, scale)     ! 2 dummies
        end do
        res(6:10) = nor_rvs(loc, scale, k)  ! 3 dummies
        call check(all(abs(res - ans) < ${k1}$tol),                            &
            msg="normal_distribution_rvs_${t1[0]}$${k1}$ failed", warn=warn)
    end subroutine test_nor_rvs_${t1[0]}$${k1}$

    #:endfor


    #:for k1, t1 in RC_KINDS_TYPES
    subroutine test_nor_rvs_default_${t1[0]}$${k1}$
        ${t1}$ :: a1(10), a2(10), mold
        integer :: i
        integer :: seed, get

        print *, "Test normal_distribution_rvs_default_${t1[0]}$${k1}$"
        seed = 25836914
        call random_seed(seed, get)

        ! explicit form with loc=0, scale=1
        #:if t1[0] == "r"
            a1 = nor_rvs(0.0_${k1}$, 1.0_${k1}$, 10)
        #:else
            a1 = nor_rvs((0.0_${k1}$, 0.0_${k1}$), (1.0_${k1}$, 1.0_${k1}$), 10)
        #:endif

        ! reset seed to reproduce same random sequence
        seed = 25836914
        call random_seed(seed, get)

        ! default mold form: mold used only to disambiguate kind
        ! For real(dp), mold is optional; for other types (including complex), it's required
        #:if t1[0] == "r" and k1 == "dp"
            a2 = nor_rvs(10)  ! mold optional for rdp only, defaults to real(dp)
        #:else
            #! mold required for all other types including complex and non-dp kinds
            #:if t1[0] == "r"
                mold = 0.0_${k1}$
            #:else
                mold = (0.0_${k1}$, 0.0_${k1}$)
            #:endif
            a2 = nor_rvs(10, mold)
        #:endif

        call check(all(a1 == a2), msg="normal_distribution_rvs_default_${t1[0]}$${k1}$ failed", warn=warn)
    end subroutine test_nor_rvs_default_${t1[0]}$${k1}$

    #:endfor




    #:for k1, t1 in RC_KINDS_TYPES
    subroutine test_nor_pdf_${t1[0]}$${k1}$

        ${t1}$ :: x1, x2(3,4), loc, scale
        integer, parameter :: k = 5
        integer :: seed, get
        real(${k1}$) :: res(3,5)
        #:if t1[0] == "r"
            #! for real type
        real(${k1}$), parameter :: ans(15) =                                   &
                         [0.215050766989949083210785218076278553_${k1}$,       &
                          0.215050766989949083210785218076278553_${k1}$,       &
                          0.215050766989949083210785218076278553_${k1}$,       &
                          0.200537626692880596839299818454439236_${k1}$,       &
                          5.66161527403268434368575024104022496E-0002_${k1}$,  &
                          0.238986957612021514867582359518579138_${k1}$,       &
                          0.265935969411942911029638292783132425_${k1}$,       &
                          0.262147558654079961109031890374902943_${k1}$,       &
                          0.249866408914952245533320687656894701_${k1}$,       &
                          3.98721117498705317877792757313696510E-0002_${k1}$,  &
                          0.265902369803533466897906694845094995_${k1}$,       &
                          0.161311603170650092038944290133124635_${k1}$,       &
                          0.249177740354276111998717092437037695_${k1}$,       &
                          0.237427217242213206474603807278971527_${k1}$,       &
                          0.155696086384122017518186260628090478_${k1}$]
        #:else
            #! for complex type
        real(${k1}$), parameter :: ans(15) =                                   &
                         [0.129377311291944176372137325120411497_${k1}$,       &
                          0.129377311291944176372137325120411497_${k1}$,       &
                          0.129377311291944176372137325120411497_${k1}$,       &
                          4.05915662853246811934977653001971736E-0002_${k1}$,  &
                          0.209143395418940756076861773161637924_${k1}$,       &
                          2.98881041363874672676853084975547667E-0002_${k1}$,  &
                          0.128679412679649127469385460133445161_${k1}$,       &
                          0.177484732473055532384223611956177231_${k1}$,       &
                          3.82205306942578982084957100753849738E-0002_${k1}$,  &
                          7.09915714309796034515515428785324918E-0002_${k1}$,  &
                          4.56126582912124629544443072483362352E-0002_${k1}$,  &
                          6.57454133967021123696499056531595921E-0002_${k1}$,  &
                          0.165161039915667041643464172210282279_${k1}$,       &
                          3.86104822953520989775015755966798359E-0002_${k1}$,  &
                          0.196922947431391188040943672442575686_${k1}$]
        #:endif

        print *, "Test normal_distribution_pdf_${t1[0]}$${k1}$"
        seed = 741852963
        call random_seed(seed, get)
        #:if t1[0] == "r"
            #! for real type
        loc = -0.5_${k1}$; scale = 1.5_${k1}$
        #:else
            #! for complex type
        loc = (-0.5_${k1}$, 0.5_${k1}$); scale = (0.5_${k1}$, 1.5_${k1}$)
        #:endif

        x1 = nor_rvs(loc, scale)
        x2 = reshape(nor_rvs(loc, scale, 12), [3,4])
        res(:,1) = nor_pdf(x1, loc, scale)
        res(:, 2:5) = nor_pdf(x2, loc, scale)
        call check(all(abs(res - reshape(ans, [3,5])) < ${k1}$tol),            &
            msg="normal_distribution_pdf_${t1[0]}$${k1}$ failed", warn=warn)
    end subroutine test_nor_pdf_${t1[0]}$${k1}$

    #:endfor




    #:for k1, t1 in RC_KINDS_TYPES
    subroutine test_nor_cdf_${t1[0]}$${k1}$

        ${t1}$ :: x1, x2(3,4), loc, scale
        integer :: seed, get
        real(${k1}$) :: res(3,5)
        #:if t1[0] == "r"
            #!for real type
        real(${k1}$), parameter :: ans(15) =                                   &
                          [7.50826305038441048487991102776953948E-0002_${k1}$, &
                           7.50826305038441048487991102776953948E-0002_${k1}$, &
                           7.50826305038441048487991102776953948E-0002_${k1}$, &
                           0.143119834108717983250834016885129863_${k1}$,      &
                           0.241425421525703182028420560765471735_${k1}$,      &
                           0.284345878626039240974266199229875972_${k1}$,      &
                           0.233239836366015928845367994433532757_${k1}$,      &
                           0.341059506137219171082517155967522896_${k1}$,      &
                           0.353156850199835111081038166086606192_${k1}$,      &
                           6.81066766396638231790017005897813244E-0002_${k1}$, &
                           4.38792331441682923984716366123285346E-0002_${k1}$, &
                           0.763679637882860826030745070304416929_${k1}$,      &
                           0.363722187587355040667876190724308059_${k1}$,      &
                           0.868187114884980488672309198087692444_${k1}$,      &
                           0.626506799809652872401992867475200722_${k1}$]
        #:else
            #! for complex type
        real(${k1}$), parameter :: ans(15) =                                   &
                          [1.07458136221563368133842063954746170E-0002_${k1}$, &
                           1.07458136221563368133842063954746170E-0002_${k1}$, &
                           1.07458136221563368133842063954746170E-0002_${k1}$, &
                           6.86483236063879585051085536740820057E-0002_${k1}$, &
                           7.95486634025192048896990048539218724E-0002_${k1}$, &
                           2.40523393996423661445007940057223384E-0002_${k1}$, &
                           3.35096768781160662250307446207445131E-0002_${k1}$, &
                           0.315778916661119434962814841317323376_${k1}$,      &
                           0.446311293878359175362094845206410428_${k1}$,      &
                           0.102010220821382542292905161748120877_${k1}$,      &
                           7.66919007012121545175655727052974512E-0002_${k1}$, &
                           0.564690968410069125818268877247699603_${k1}$,      &
                           0.708769523556518785240723539383512333_${k1}$,      &
                           6.40553790808161720088070925562830659E-0002_${k1}$, &
                           5.39999153072107729358158443133850711E-0002_${k1}$]
        #:endif

        print *, "Test normal_distribution_cdf_${t1[0]}$${k1}$"
        seed = 369147582
        call random_seed(seed, get)
        #:if t1[0] == "r"
            #! for real type
        loc = -1.0_${k1}$; scale = 2.0_${k1}$
        #:else
            #! for complex type
        loc = (-1.0_${k1}$, 1.0_${k1}$); scale = (1.7_${k1}$, 2.4_${k1}$)
        #:endif

        x1 = nor_rvs(loc, scale)
        x2 = reshape(nor_rvs(loc, scale, 12), [3,4])
        res(:,1) = nor_cdf(x1, loc, scale)
        res(:, 2:5) = nor_cdf(x2, loc, scale)
        call check(all(abs(res - reshape(ans,[3,5])) < ${k1}$tol),             &
            msg="normal_distribution_cdf_${t1[0]}$${k1}$ failed", warn=warn)
    end subroutine test_nor_cdf_${t1[0]}$${k1}$

    #:endfor

end program test_distribution_normal