https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_distPower@routines.inc.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 41 41 100.0 %
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_distPower](@ref pm_distPower).
      19             : !>
      20             : !>  \finmain
      21             : !>
      22             : !>  \author
      23             : !>  \AmirShahmoradi, Oct 16, 2009, 11:14 AM, Michigan
      24             : 
      25             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      26             : 
      27             :         !%%%%%%%%%%%%%%%%%%%%%%%
      28             : #if     getPowerLogPDFNF_ENABLED
      29             :         !%%%%%%%%%%%%%%%%%%%%%%%
      30             : 
      31      652884 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@getPowerLogPDFNF(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
      32             : #if     ALD_ENABLED
      33      582590 :         logPDFNF = log(alpha) - alpha * logMaxX
      34             : #elif   ALL_ENABLED
      35      210882 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getPowerLogPDFNF(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      36       70294 :         logPDFNF = log(alpha) - getLogSubExp(smaller = alpha * logMinX, larger = alpha * logMaxX)
      37             : #else
      38             : #error  "Unrecognized interface."
      39             : #endif
      40             : 
      41             :         !%%%%%%%%%%%%%%%%%%%%%
      42             : #elif   getPowerLogPDF_ENABLED
      43             :         !%%%%%%%%%%%%%%%%%%%%%
      44             : 
      45     1934784 :         CHECK_ASSERTION(__LINE__, logx <= logMaxX, SK_"@getPowerLogPDF(): The condition `logx <= logMaxX` must hold. logx, logMaxX = "//getStr([logx, logMaxX]))
      46             : #if     ALD_ENABLED
      47      578358 :         call setPowerLogPDF(logPDF, logx, alpha, getPowerLogPDFNF(alpha, logMaxX))
      48             : #elif   ALL_ENABLED
      49      199710 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@getPowerLogPDF(): The condition `logMinX <= logx` must hold. logMinX, logx = "//getStr([logMinX, logx]))
      50       66570 :         call setPowerLogPDF(logPDF, logx, alpha, getPowerLogPDFNF(alpha, logMinX, logMaxX))
      51             : #else
      52             : #error  "Unrecognized interface."
      53             : #endif
      54             : 
      55             :         !%%%%%%%%%%%%%%%%%%%%%
      56             : #elif   setPowerLogPDF_ENABLED
      57             :         !%%%%%%%%%%%%%%%%%%%%%
      58             : 
      59      758648 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@setPowerLogPDF(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
      60      758648 :         logPDF = logPDFNF + (alpha - 1._RKC) * logx
      61             : 
      62             :         !%%%%%%%%%%%%%%%%%%%%%%%
      63             : #elif   getPowerLogCDFNF_ENABLED
      64             :         !%%%%%%%%%%%%%%%%%%%%%%%
      65             : 
      66             : #if     ALD_ENABLED
      67       41472 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@getPowerLogCDFNF(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
      68       41472 :         logCDFNF = - alpha * logMaxX
      69             : #elif   ALL_ENABLED
      70             :         real(RKC) :: alphaLogMinX
      71       48968 :         alphaLogMinX = alpha * logMinX
      72       48968 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@getPowerLogCDFNF(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
      73      146904 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getPowerLogCDFNF(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
      74       48968 :         logCDFNF = alphaLogMinX - getLogSubExp(smaller = alphaLogMinX, larger = alpha * logMaxX)
      75             : #else
      76             : #error  "Unrecognized interface."
      77             : #endif
      78             : 
      79             :         !%%%%%%%%%%%%%%%%%%%%%
      80             : #elif   getPowerLogCDF_ENABLED
      81             :         !%%%%%%%%%%%%%%%%%%%%%
      82             : 
      83       20208 :         CHECK_ASSERTION(__LINE__, logx <= logMaxX, SK_"@getPowerLogCDF(): The condition `logx <= logMaxX` must hold. logx, logMaxX = "//getStr([logx, logMaxX]))
      84             : #if     ALD_ENABLED
      85        5106 :         call setPowerLogCDF(logCDF, logx, alpha, getPowerLogCDFNF(alpha, logMaxX))
      86             : #elif   ALL_ENABLED
      87        1630 :         call setPowerLogCDF(logCDF, logx, alpha, logMinX, getPowerLogCDFNF(alpha, logMinX, logMaxX))
      88             : #else
      89             : #error  "Unrecognized interface."
      90             : #endif
      91             : 
      92             :         !%%%%%%%%%%%%%%%%%%%%%
      93             : #elif   setPowerLogCDF_ENABLED
      94             :         !%%%%%%%%%%%%%%%%%%%%%
      95             : 
      96       62376 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@setPowerLogCDF(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
      97             : #if     ALD_ENABLED
      98       33812 :         logCDF = logCDFNF + alpha * logx
      99      169060 :         CHECK_ASSERTION(__LINE__, exp(logCDF) <= 1._RKC + sqrt(epsilon(0._RKC)), SK_"@setPowerLogCDF(): The condition `logCDF <= 0._RKC` must hold. The input arguments are inconsistent or `logx` is out of support. logCDF, logx, alpha, logCDFNF = "//getStr([logCDF, logx, alpha, logCDFNF]))
     100             : #elif   ALL_ENABLED
     101       85692 :         CHECK_ASSERTION(__LINE__, logMinX <= logx, SK_"@setPowerLogCDF(): The condition `logMinX <= logx` must hold. logMinX, logx = "//getStr([logMinX, logx]))
     102       28564 :         if (logx /= logMinX) then
     103       27363 :             logCDF = logCDFNF + getLogSubExp(smaller = 0._RKC, larger = alpha * (logx - logMinX))
     104      191541 :             CHECK_ASSERTION(__LINE__, exp(logCDF) <= 1._RKC + sqrt(epsilon(0._RKC)), SK_"@setPowerLogCDF(): The condition `logCDF <= 0._RKC` must hold. The input arguments are inconsistent or `logx` is out of support. epsilon(0._RKC), logCDF, logx, alpha, logMinX = "//getStr([epsilon(0._RKC), logCDF, logx, alpha, logMinX, logCDFNF]))
     105             :         else
     106        1201 :             logCDF = -log(huge(0._RKC))
     107             :         end if
     108             : #else
     109             : #error  "Unrecognized interface."
     110             : #endif
     111             : 
     112             :         !%%%%%%%%%%%%%%%%%%%%%%
     113             : #elif   getPowerLogQuan_ENABLED
     114             :         !%%%%%%%%%%%%%%%%%%%%%%
     115             : 
     116             : #if     ALD_ENABLED
     117        6004 :         call setPowerLogQuan(logx, logCDF, alpha, getPowerLogCDFNF(alpha, logMaxX))
     118             : #elif   ALL_ENABLED
     119       18012 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getPowerLogQuan(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
     120        6004 :         call setPowerLogQuan(logx, logCDF, alpha, logMinX, getPowerLogCDFNF(alpha, logMinX, logMaxX))
     121             : #else
     122             : #error  "Unrecognized interface."
     123             : #endif
     124             : 
     125             :         !%%%%%%%%%%%%%%%%%%%%%%
     126             : #elif   setPowerLogQuan_ENABLED
     127             :         !%%%%%%%%%%%%%%%%%%%%%%
     128             : 
     129       84048 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@setPowerLogQuan(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
     130      252144 :         CHECK_ASSERTION(__LINE__, logCDF <= epsilon(0._RKC), SK_"@setPowerLogQuan(): The condition `logCDF <= epsilon(0._RKC)` must hold. logCDF, epsilon(0._RKC) = "//getStr([logCDF, epsilon(0._RKC)]))
     131             : #if         LLALD_ENABLED
     132       36024 :             logx = (logCDF - logCDFNF) / alpha
     133             : #elif       LLALL_ENABLED
     134       48024 :             logx = logMinX + getLogAddExp(getMinMax(logCDF - logCDFNF, 0._RKC)) / alpha
     135             : #else
     136             : #error      "Unrecognized interface."
     137             : #endif
     138             : 
     139             :         !%%%%%%%%%%%%%%%%%%%%%%
     140             : #elif   getPowerLogRand_ENABLED
     141             :         !%%%%%%%%%%%%%%%%%%%%%%
     142             : 
     143             : #if     ALD_ENABLED
     144       12004 :         call setPowerLogRand(logRand, getNegExpRand(1._RKC), alpha, getPowerLogCDFNF(alpha, logMaxX))
     145             : #elif   ALL_ENABLED
     146       30012 :         CHECK_ASSERTION(__LINE__, logMinX < logMaxX, SK_"@getPowerLogRand(): The condition `logMinX < logMaxX` must hold. logMinX, logMaxX = "//getStr([logMinX, logMaxX]))
     147       10004 :         call setPowerLogRand(logRand, getNegExpRand(1._RKC), alpha, logMinX, getPowerLogCDFNF(alpha, logMinX, logMaxX))
     148             : #else
     149             : #error  "Unrecognized interface."
     150             : #endif
     151             : 
     152             :         !%%%%%%%%%%%%%%%%%%%%%%
     153             : #elif   setPowerLogRand_ENABLED
     154             :         !%%%%%%%%%%%%%%%%%%%%%%
     155             : 
     156       44016 :         CHECK_ASSERTION(__LINE__, alpha > 0._RKC, SK_"@setPowerLogRand(): The condition `alpha > 0._RKC` must hold. alpha = "//getStr(alpha))
     157       44016 :         CHECK_ASSERTION(__LINE__, negExpRand <= 0._RKC, SK_"@setPowerLogRand(): The condition `negExpRand <= 0._RKC` must hold. alpha = "//getStr(negExpRand))
     158             : #if     LNALD_ENABLED
     159       24008 :         call setPowerLogQuan(logRand, negExpRand, alpha, logCDFNF)
     160             : #elif   LNALL_ENABLED
     161       20008 :         call setPowerLogQuan(logRand, negExpRand, alpha, logMinX, logCDFNF)
     162             : #else
     163             : #error  "Unrecognized interface."
     164             : #endif
     165             : 
     166             : #else
     167             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     168             : #error  "Unrecognized interface."
     169             :         !%%%%%%%%%%%%%%%%%%%%%%%%
     170             : #endif

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