https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_distPareto@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 43 44 97.7 %
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_distPareto](@ref pm_distPareto).
      19             : !>
      20             : !>  \finmain
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, Oct 16, 2009, 11:14 AM, Michigan
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             : 
      28             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      29             : #if     getParetoLogPDFNF_ENABLED
      30             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      31             : 
      32      391109 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@getParetoLogPDFNF(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
      33             : #if     ALD_ENABLED
      34      244855 :         logPDFNF = log(-alpha) - alpha * logMinX
      35             : #elif   ALL_ENABLED
      36      438762 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getParetoLogPDFNF(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      37      146254 :         logPDFNF = log(-alpha) - getLogSubExp(smaller = alpha * logMaxX, larger = alpha * logMinX)
      38             : #else
      39             : #error  "Unrecognized interface."
      40             : #endif
      41             : 
      42             :         !%%%%%%%%%%%%%%%%%%%%%%
      43             : #elif   getParetoLogPDF_ENABLED
      44             :         !%%%%%%%%%%%%%%%%%%%%%%
      45             : 
      46     1159179 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@getParetoLogPDF(): The condition `logMinX <= logx` must hold. logx, logMaxX = "//getStr([logMinX, logx]))
      47             : #if     ALD_ENABLED
      48      243623 :         call setParetoLogPDF(logPDF, logx, alpha, getParetoLogPDFNF(alpha, logMinX))
      49             : #elif   ALL_ENABLED
      50      428310 :         CHECK_ASSERTION(__LINE__, logx <= logMaxX, SK_"@getParetoLogPDF(): The condition `logx <= logMaxX` must hold. logMinX, logx = "//getStr([logx, logMaxX]))
      51      142770 :         call setParetoLogPDF(logPDF, logx, alpha, getParetoLogPDFNF(alpha, logMinX, logMaxX))
      52             : #else
      53             : #error  "Unrecognized interface."
      54             : #endif
      55             : 
      56             :         !%%%%%%%%%%%%%%%%%%%%%%
      57             : #elif   setParetoLogPDF_ENABLED
      58             :         !%%%%%%%%%%%%%%%%%%%%%%
      59             : 
      60      464057 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@setParetoLogPDF(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
      61      464057 :         logPDF = logPDFNF + (alpha - 1._RKC) * logx
      62             : 
      63             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      64             : #elif   getParetoLogCDFNF_ENABLED
      65             :         !%%%%%%%%%%%%%%%%%%%%%%%%
      66             : 
      67             : #if     ALD_ENABLED
      68          12 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@getParetoLogCDFNF(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
      69           0 :         logCDFNF = 0._RKC
      70             : #elif   ALL_ENABLED
      71             :         real(RKC) :: alphaLogMinX
      72       52976 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@getParetoLogCDFNF(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
      73      158928 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getParetoLogCDFNF(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      74       52976 :         alphaLogMinX = alpha * logMinX
      75       52976 :         logCDFNF = alphaLogMinX - getLogSubExp(smaller = alpha * logMaxX, larger = alphaLogMinX)
      76             : #else
      77             : #error  "Unrecognized interface."
      78             : #endif
      79             : 
      80             :         !%%%%%%%%%%%%%%%%%%%%%%
      81             : #elif   getParetoLogCDF_ENABLED
      82             :         !%%%%%%%%%%%%%%%%%%%%%%
      83             : 
      84       10083 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@getParetoLogCDF(): The condition `logx <= logMinX` must hold. logMinX, logx = "//getStr([logMinX, logx]))
      85             : #if     ALD_ENABLED
      86        1731 :         call setParetoLogCDF(logCDF, logx, alpha, logMinX)
      87             : #elif   ALL_ENABLED
      88        4890 :         CHECK_ASSERTION(__LINE__, logx <= logMaxX, SK_"@getParetoLogCDF(): The condition `logx <= logMaxX` must hold. logx, logMaxX = "//getStr([logx, logMaxX]))
      89        1630 :         call setParetoLogCDF(logCDF, logx, alpha, logMinX, getParetoLogCDFNF(alpha, logMinX, logMaxX))
      90             : #else
      91             : #error  "Unrecognized interface."
      92             : #endif
      93             : 
      94             :         !%%%%%%%%%%%%%%%%%%%%%%
      95             : #elif   setParetoLogCDF_ENABLED
      96             :         !%%%%%%%%%%%%%%%%%%%%%%
      97             : 
      98       55639 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@setParetoLogCDF(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
      99      166917 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@setParetoLogCDF(): The condition `logMinX <= logx` must hold. logMinX, logx = "//getStr([logMinX, logx]))
     100       55639 :         if (logx /= logMinX) then
     101       53238 :             logCDF = getLogSubExp(smaller = alpha * (logx - logMinX), larger = 0._RKC)
     102             : #if         ALD_ENABLED
     103      129315 :             CHECK_ASSERTION(__LINE__, exp(logCDF) <= 1._RKC + sqrt(epsilon(0._RKC)), \
     104             :             SK_"@setParetoLogCDF(): The condition `logCDF <= 0._RKC` must hold. The input arguments are inconsistent or `logx` is out of support. logCDF, logx, alpha, logMinX = "//getStr([logCDF, logx, alpha, logMinX]))
     105             : #elif       ALL_ENABLED
     106       27375 :             logCDF = logCDF + logCDFNF
     107      164250 :             CHECK_ASSERTION(__LINE__, exp(logCDF) <= 1._RKC + sqrt(epsilon(0._RKC)), \
     108             :             SK_"@setParetoLogCDF(): The condition `logCDF <= 0._RKC` must hold. The input arguments are inconsistent or `logx` is out of support. logCDF, logx, alpha, logMinX, logCDFNF = "//getStr([logCDF, logx, alpha, logMinX, logCDFNF]))
     109             : #else
     110             : #error      "Unrecognized interface."
     111             : #endif
     112             :         else
     113        2401 :             logCDF = -log(huge(0._RKC))
     114             :         end if
     115             : 
     116             :         !%%%%%%%%%%%%%%%%%%%%%%%
     117             : #elif   getParetoLogQuan_ENABLED
     118             :         !%%%%%%%%%%%%%%%%%%%%%%%
     119             : 
     120             : #if     ALD_ENABLED
     121        6004 :         call setParetoLogQuan(logx, logCDF, alpha, logMinX)
     122             : #elif   ALL_ENABLED
     123       18012 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getParetoLogQuan(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
     124        6004 :         call setParetoLogQuan(logx, logCDF, alpha, logMinX, getParetoLogCDFNF(alpha, logMinX, logMaxX))
     125             : #else
     126             : #error  "Unrecognized interface."
     127             : #endif
     128             : 
     129             :         !%%%%%%%%%%%%%%%%%%%%%%%
     130             : #elif   setParetoLogQuan_ENABLED
     131             :         !%%%%%%%%%%%%%%%%%%%%%%%
     132             : 
     133       84064 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@setParetoLogQuan(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
     134             : #if     LLALD_ENABLED
     135       30024 :         CHECK_ASSERTION(__LINE__, logCDF < 0._RKC, SK_"@setParetoLogQuan(): The condition `logCDF < 0._RKC` must hold. logCDF = "//getStr(logCDF))
     136       30024 :         logx = logMinX + getLogSubExp(smaller = logCDF, larger = 0._RKC) / alpha
     137             : #elif   LLALL_ENABLED
     138       54040 :         CHECK_ASSERTION(__LINE__, logCDF <= 0._RKC, SK_"@setParetoLogQuan(): The condition `logCDF < 0._RKC` must hold. logCDF = "//getStr(logCDF))
     139       54040 :         logx = logMinX + getLogSubExp(smaller = logCDF - logCDFNF, larger = 0._RKC) / alpha
     140             : #else
     141             : #error  "Unrecognized interface."
     142             : #endif
     143             : 
     144             :         !%%%%%%%%%%%%%%%%%%%%%%%
     145             : #elif   getParetoLogRand_ENABLED
     146             :         !%%%%%%%%%%%%%%%%%%%%%%%
     147             : 
     148             : #if     ALD_ENABLED
     149             :         !call setParetoLogRand(logRand, getNegExpRand(sigma = 1._RKC), alpha, logMinX)
     150       10004 :         call setParetoLogRand(logRand, log(1._RKC - getUnifRand(0._RKC, 1._RKC)), alpha, logMinX)
     151             : #elif   ALL_ENABLED
     152       36012 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getParetoLogRand(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
     153             :         !call setParetoLogRand(logRand, getNegExpRand(sigma = 1._RKC), alpha, logMinX, getParetoLogCDFNF(alpha, logMinX, logMaxX))
     154       12004 :         call setParetoLogRand(logRand, log(1._RKC - getUnifRand(0._RKC, 1._RKC)), alpha, logMinX, getParetoLogCDFNF(alpha, logMinX, logMaxX))
     155             : #else
     156             : #error  "Unrecognized interface."
     157             : #endif
     158             : 
     159             :         !%%%%%%%%%%%%%%%%%%%%%%%
     160             : #elif   setParetoLogRand_ENABLED
     161             :         !%%%%%%%%%%%%%%%%%%%%%%%
     162             : 
     163       44016 :         CHECK_ASSERTION(__LINE__, alpha < 0._RKC, SK_"@setParetoLogRand(): The condition `alpha < 0._RKC` must hold. alpha = "//getStr(alpha))
     164             : #if     LNALD_ENABLED
     165       18008 :         CHECK_ASSERTION(__LINE__, negExpRand < 0._RKC, SK_"@setParetoLogRand(): The condition `negExpRand < 0._RKC` must hold. alpha = "//getStr(negExpRand))
     166       18008 :         call setParetoLogQuan(logRand, negExpRand, alpha, logMinX)
     167             : #elif   LNALL_ENABLED
     168       26008 :         CHECK_ASSERTION(__LINE__, negExpRand <= 0._RKC, SK_"@setParetoLogRand(): The condition `negExpRand <= 0._RKC` must hold. alpha = "//getStr(negExpRand))
     169       26008 :         call setParetoLogQuan(logRand, negExpRand, alpha, logMinX, logCDFNF)
     170             : #else
     171             : #error  "Unrecognized interface."
     172             : #endif
     173             : 
     174             : #else
     175             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     176             : #error  "Unrecognized interface."
     177             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     178             : #endif

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