ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_mathFisher.F90
Go to the documentation of this file.
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
80
81!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
82
84
85 use pm_kind, only: IK, RK, SK
86
87 implicit none
88
89 character(*, SK), parameter :: MODULE_NAME = "@pm_mathFisher"
90
91!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
92
153 interface getFisher
154
155 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
156
157#if RK5_ENABLED
158 PURE elemental module function getFisherFDD_RK5(val) result(fisherz)
159#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
160 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFDD_RK5
161#endif
162 use pm_kind, only: RKG => RK5
163 real(RKG) , intent(in) :: val
164 real(RKG) :: fisherz
165 end function
166#endif
167
168#if RK4_ENABLED
169 PURE elemental module function getFisherFDD_RK4(val) result(fisherz)
170#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
171 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFDD_RK4
172#endif
173 use pm_kind, only: RKG => RK4
174 real(RKG) , intent(in) :: val
175 real(RKG) :: fisherz
176 end function
177#endif
178
179#if RK3_ENABLED
180 PURE elemental module function getFisherFDD_RK3(val) result(fisherz)
181#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
182 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFDD_RK3
183#endif
184 use pm_kind, only: RKG => RK3
185 real(RKG) , intent(in) :: val
186 real(RKG) :: fisherz
187 end function
188#endif
189
190#if RK2_ENABLED
191 PURE elemental module function getFisherFDD_RK2(val) result(fisherz)
192#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
193 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFDD_RK2
194#endif
195 use pm_kind, only: RKG => RK2
196 real(RKG) , intent(in) :: val
197 real(RKG) :: fisherz
198 end function
199#endif
200
201#if RK1_ENABLED
202 PURE elemental module function getFisherFDD_RK1(val) result(fisherz)
203#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
204 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFDD_RK1
205#endif
206 use pm_kind, only: RKG => RK1
207 real(RKG) , intent(in) :: val
208 real(RKG) :: fisherz
209 end function
210#endif
211
212 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
213
214#if RK5_ENABLED
215 PURE elemental module function getFisherFLU_RK5(val, lb, ub) result(fisherz)
216#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
217 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFLU_RK5
218#endif
219 use pm_kind, only: RKG => RK5
220 real(RKG) , intent(in) :: val, lb, ub
221 real(RKG) :: fisherz
222 end function
223#endif
224
225#if RK4_ENABLED
226 PURE elemental module function getFisherFLU_RK4(val, lb, ub) result(fisherz)
227#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
228 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFLU_RK4
229#endif
230 use pm_kind, only: RKG => RK4
231 real(RKG) , intent(in) :: val, lb, ub
232 real(RKG) :: fisherz
233 end function
234#endif
235
236#if RK3_ENABLED
237 PURE elemental module function getFisherFLU_RK3(val, lb, ub) result(fisherz)
238#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
239 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFLU_RK3
240#endif
241 use pm_kind, only: RKG => RK3
242 real(RKG) , intent(in) :: val, lb, ub
243 real(RKG) :: fisherz
244 end function
245#endif
246
247#if RK2_ENABLED
248 PURE elemental module function getFisherFLU_RK2(val, lb, ub) result(fisherz)
249#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
250 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFLU_RK2
251#endif
252 use pm_kind, only: RKG => RK2
253 real(RKG) , intent(in) :: val, lb, ub
254 real(RKG) :: fisherz
255 end function
256#endif
257
258#if RK1_ENABLED
259 PURE elemental module function getFisherFLU_RK1(val, lb, ub) result(fisherz)
260#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
261 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherFLU_RK1
262#endif
263 use pm_kind, only: RKG => RK1
264 real(RKG) , intent(in) :: val, lb, ub
265 real(RKG) :: fisherz
266 end function
267#endif
268
269 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
270
271 end interface
272
273!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
274
335 interface getFisherInv
336
337 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
338
339#if RK5_ENABLED
340 pure elemental module function getFisherInvFDD_RK5(fisherz) result(val)
341#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
342 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFDD_RK5
343#endif
344 use pm_kind, only: RKG => RK5
345 real(RKG) , intent(in) :: fisherz
346 real(RKG) :: val
347 end function
348#endif
349
350#if RK4_ENABLED
351 pure elemental module function getFisherInvFDD_RK4(fisherz) result(val)
352#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
353 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFDD_RK4
354#endif
355 use pm_kind, only: RKG => RK4
356 real(RKG) , intent(in) :: fisherz
357 real(RKG) :: val
358 end function
359#endif
360
361#if RK3_ENABLED
362 pure elemental module function getFisherInvFDD_RK3(fisherz) result(val)
363#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
364 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFDD_RK3
365#endif
366 use pm_kind, only: RKG => RK3
367 real(RKG) , intent(in) :: fisherz
368 real(RKG) :: val
369 end function
370#endif
371
372#if RK2_ENABLED
373 pure elemental module function getFisherInvFDD_RK2(fisherz) result(val)
374#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
375 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFDD_RK2
376#endif
377 use pm_kind, only: RKG => RK2
378 real(RKG) , intent(in) :: fisherz
379 real(RKG) :: val
380 end function
381#endif
382
383#if RK1_ENABLED
384 pure elemental module function getFisherInvFDD_RK1(fisherz) result(val)
385#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
386 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFDD_RK1
387#endif
388 use pm_kind, only: RKG => RK1
389 real(RKG) , intent(in) :: fisherz
390 real(RKG) :: val
391 end function
392#endif
393
394 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
395
396#if RK5_ENABLED
397 PURE elemental module function getFisherInvFLU_RK5(fisherz, lb, ub) result(val)
398#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
399 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFLU_RK5
400#endif
401 use pm_kind, only: RKG => RK5
402 real(RKG) , intent(in) :: fisherz, lb, ub
403 real(RKG) :: val
404 end function
405#endif
406
407#if RK4_ENABLED
408 PURE elemental module function getFisherInvFLU_RK4(fisherz, lb, ub) result(val)
409#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
410 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFLU_RK4
411#endif
412 use pm_kind, only: RKG => RK4
413 real(RKG) , intent(in) :: fisherz, lb, ub
414 real(RKG) :: val
415 end function
416#endif
417
418#if RK3_ENABLED
419 PURE elemental module function getFisherInvFLU_RK3(fisherz, lb, ub) result(val)
420#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
421 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFLU_RK3
422#endif
423 use pm_kind, only: RKG => RK3
424 real(RKG) , intent(in) :: fisherz, lb, ub
425 real(RKG) :: val
426 end function
427#endif
428
429#if RK2_ENABLED
430 PURE elemental module function getFisherInvFLU_RK2(fisherz, lb, ub) result(val)
431#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
432 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFLU_RK2
433#endif
434 use pm_kind, only: RKG => RK2
435 real(RKG) , intent(in) :: fisherz, lb, ub
436 real(RKG) :: val
437 end function
438#endif
439
440#if RK1_ENABLED
441 PURE elemental module function getFisherInvFLU_RK1(fisherz, lb, ub) result(val)
442#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
443 !DEC$ ATTRIBUTES DLLEXPORT :: getFisherInvFLU_RK1
444#endif
445 use pm_kind, only: RKG => RK1
446 real(RKG) , intent(in) :: fisherz, lb, ub
447 real(RKG) :: val
448 end function
449#endif
450
451 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
452
453 end interface
454
455!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
456
457end module pm_mathFisher ! LCOV_EXCL_LINE
Generate and return the inverse Fisher transformation of the input Fisher z value.
Generate and return the Fisher transformation of the input Fisher z value.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK5
Definition: pm_kind.F90:478
integer, parameter RK4
Definition: pm_kind.F90:489
integer, parameter RK2
Definition: pm_kind.F90:511
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
integer, parameter RK3
Definition: pm_kind.F90:500
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RK1
Definition: pm_kind.F90:522
This module contains procedures and generic interfaces for evaluating the Fisher transformation and i...
character(*, SK), parameter MODULE_NAME