https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_distLogUnif@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 31 32 96.9 %
Date: 2024-04-08 03:18:57 Functions: 0 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 the implementation of procedures in [pm_distLogUnif](@ref pm_distLogUnif).
      19             : !>
      20             : !>  \finmain
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, Oct 16, 2009, 11:14 AM, Michigan
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         !%%%%%%%%%%%%%%%%%%%%%%
      28             : #if     getLogUnifPDFNF_ENABLED
      29             :         !%%%%%%%%%%%%%%%%%%%%%%
      30             : 
      31       65661 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getLogUnifPDFNF(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      32       21887 :         pdfnf = 1._RKC / (logMaxX - logMinX)
      33             : 
      34             :         !%%%%%%%%%%%%%%%%%%%%
      35             : #elif   getLogUnifPDF_ENABLED
      36             :         !%%%%%%%%%%%%%%%%%%%%
      37             : 
      38        1013 :         CHECK_ASSERTION(__LINE__, 0._RKC < minx, SK_"@getLogUnifPDF(): The condition `0._RKC < minx` must hold. minx = "//getStr(minx))
      39        3039 :         CHECK_ASSERTION(__LINE__, minx <= x, SK_"@getLogUnifPDF(): The condition `minx <= x` must hold. minx, x = "//getStr([minx, x]))
      40        3039 :         CHECK_ASSERTION(__LINE__, x <= maxx, SK_"@getLogUnifPDF(): The condition `x <= maxx` must hold. x, exp(logMaxX) = "//getStr([x, maxx]))
      41        1013 :         call setLogUnifPDF(pdf, x, getLogUnifPDFNF(log(minx), log(maxx)))
      42             : 
      43             :         !%%%%%%%%%%%%%%%%%%%%
      44             : #elif   setLogUnifPDF_ENABLED
      45             :         !%%%%%%%%%%%%%%%%%%%%
      46             : 
      47        2026 :         CHECK_ASSERTION(__LINE__, 0._RKC <= x, SK_"@setLogUnifPDF(): The condition `0._RKC <= x` must hold. x = "//getStr(x))
      48        2026 :         CHECK_ASSERTION(__LINE__, 0._RKC < pdfnf, SK_"@setLogUnifPDF(): The condition `0._RKC <= pdfnf` must hold. pdfnf = "//getStr(pdfnf))
      49        2026 :         pdf = pdfnf / x
      50             : 
      51             :         !%%%%%%%%%%%%%%%%%%%%
      52             : #elif   getLogUnifCDF_ENABLED
      53             :         !%%%%%%%%%%%%%%%%%%%%
      54             : 
      55        5754 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@getLogUnifCDF(): The condition `logx <= logMinX` must hold. logMinX, logx = "//getStr([logMinX, logx]))
      56        5754 :         CHECK_ASSERTION(__LINE__, logx <= logMaxX, SK_"@getLogUnifCDF(): The condition `logx <= logMaxX` must hold. logx, logMaxX = "//getStr([logx, logMaxX]))
      57        1918 :         call setLogUnifCDF(cdf, logx, logMinX, getLogUnifPDFNF(logMinX, logMaxX))
      58             : 
      59             :         !%%%%%%%%%%%%%%%%%%%%
      60             : #elif   setLogUnifCDF_ENABLED
      61             :         !%%%%%%%%%%%%%%%%%%%%
      62             : 
      63       11508 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@setLogUnifCDF(): The condition `logMinX <= logx` must hold. logMinX, logx = "//getStr([logMinX, logx]))
      64        3836 :         if (logMinX < logx) then
      65        3836 :             cdf = (logx - logMinX) * pdfnf
      66       15344 :             CHECK_ASSERTION(__LINE__, cdf < 1._RKC + sqrt(epsilon(0._RKC)), SK_"@setLogUnifCDF(): The condition `cdf <= 1._RKC` must hold. The input arguments are inconsistent or `logx` is out of support. cdf, logx, logMinX = "//getStr([cdf, logx, logMinX]))
      67             :         else
      68           0 :             cdf = 0._RKC
      69             :         end if
      70             : 
      71             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      72             : #elif   getLogUnifLogQuan_ENABLED
      73             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      74             : 
      75       12012 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getLogUnifLogQuan(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      76        4004 :         call setLogUnifLogQuan(logx, cdf, logMinX, getLogUnifPDFNF(logMinX, logMaxX))
      77             : 
      78             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      79             : #elif   setLogUnifLogQuan_ENABLED
      80             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      81             : 
      82        8008 :         CHECK_ASSERTION(__LINE__, 0._RKC <= cdf, SK_"@setLogUnifLogQuan(): The condition `0._RKC <= cdf` must hold. cdf = "//getStr(cdf))
      83        8008 :         CHECK_ASSERTION(__LINE__, 1._RKC >= cdf, SK_"@setLogUnifLogQuan(): The condition `1._RKC >= cdf` must hold. cdf = "//getStr(cdf))
      84        8008 :         logx = logMinX + cdf / pdfnf
      85       24024 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@setLogUnifLogQuan(): The condition `logMinX <= logx` must hold. The input parameters are inconsistent. logMinX, logx = "//getStr([logMinX, logx]))
      86             : 
      87             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      88             : #elif   getLogUnifRand_ENABLED && MM_ENABLED && IK_ENABLED
      89             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      90             : 
      91             :         integer, parameter :: RKC = selected_real_kind(r = range(rand), radix = radix(rand))
      92         520 :         rand = int(getLogUnifRand(real(minx, RKC), real(maxx, RKC)), kind = IKC)
      93             :         !use pm_kind, only: RKS, RKD, RKQ, RKHR
      94             :         !if (real(huge(rand), RKH) < real(huge(0._RKS), RKH)) then
      95             :         !    rand = int(getLogUnifRand(real(minx, RKS), real(maxx, RKS)), kind = IKC)
      96             :         !    return
      97             :         !elseif (real(huge(rand), RKH) < real(huge(0._RKD), RKH)) then
      98             :         !    rand = int(getLogUnifRand(real(minx, RKD), real(maxx, RKD)), kind = IKC)
      99             :         !    return
     100             :         !else
     101             :         !    if (0 < RKH) then ! LCOV_EXCL_LINE
     102             :         !        if (real(huge(rand), RKH) < huge(0._RKH)) then ! LCOV_EXCL_LINE
     103             :         !            rand = int(getLogUnifRand(real(minx, RKH), real(maxx, RKH)), kind = IKC) ! LCOV_EXCL_LINE
     104             :         !            return ! LCOV_EXCL_LINE
     105             :         !        end if
     106             :         !    end if
     107             :         !end if
     108             :         !error stop "The kind of the input integer requires real numbers that can handle much larger values than what is supported by the available single, double, or quad precision real kinds. You must be living in the year 3000!" ! LCOV_EXCL_LINE
     109             : 
     110             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     111             : #elif   getLogUnifRand_ENABLED && MM_ENABLED && RK_ENABLED
     112             :         !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     113             : 
     114             :         real(RKC) :: logMinX, logRand
     115       13572 :         CHECK_ASSERTION(__LINE__, minx < maxx, SK_"@getLogUnifRand(): The condition `minx < maxx` must hold. minx, maxx = "//getStr([minx, maxx]))
     116        4524 :         logMinX = log(minx)
     117        4524 :         call setLogUnifLogRand(logRand, getUnifRand(0._RKC, 1._RKC), logMinX, getLogUnifPDFNF(logMinX, log(maxx)))
     118        4524 :         rand = exp(logRand)
     119             : 
     120             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     121             : #elif   setLogUnifLogRand_ENABLED
     122             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     123             : 
     124        8528 :         CHECK_ASSERTION(__LINE__, 0._RKC <= urand, SK_"@setLogUnifLogRand(): The condition `0._RKC <= urand` must hold. urand = "//getStr(urand))
     125        8528 :         CHECK_ASSERTION(__LINE__, 1._RKC >= urand, SK_"@setLogUnifLogRand(): The condition `1._RKC >= urand` must hold. urand = "//getStr(urand))
     126        8528 :         logRand = logMinX + urand / pdfnf
     127       25584 :         CHECK_ASSERTION(__LINE__, logMinX <= logRand, SK_"@setLogUnifLogRand(): The condition `logMinX <= logRand` must hold. The input parameters are inconsistent. logMinX, logRand = "//getStr([logMinX, logRand]))
     128             : 
     129             : #else
     130             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     131             : #error  "Unrecognized interface."
     132             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     133             : #endif

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