The ParaMonte Documentation Website
Current view: top level - kernel - Constants_mod.f90 (source / functions) Hit Total Coverage
Test: ParaMonte 1.5.1 :: Coarray Parallel Kernel - Code Coverage Report Lines: 7 7 100.0 %
Date: 2021-01-08 12:59:07 Functions: 2 2 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       2             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
       3             : !!!!
       4             : !!!!   MIT License
       5             : !!!!
       6             : !!!!   ParaMonte: plain powerful parallel Monte Carlo library.
       7             : !!!!
       8             : !!!!   Copyright (C) 2012-present, The Computational Data Science Lab
       9             : !!!!
      10             : !!!!   This file is part of the ParaMonte library.
      11             : !!!!
      12             : !!!!   Permission is hereby granted, free of charge, to any person obtaining a 
      13             : !!!!   copy of this software and associated documentation files (the "Software"), 
      14             : !!!!   to deal in the Software without restriction, including without limitation 
      15             : !!!!   the rights to use, copy, modify, merge, publish, distribute, sublicense, 
      16             : !!!!   and/or sell copies of the Software, and to permit persons to whom the 
      17             : !!!!   Software is furnished to do so, subject to the following conditions:
      18             : !!!!
      19             : !!!!   The above copyright notice and this permission notice shall be 
      20             : !!!!   included in all copies or substantial portions of the Software.
      21             : !!!!
      22             : !!!!   THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, 
      23             : !!!!   EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF 
      24             : !!!!   MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. 
      25             : !!!!   IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, 
      26             : !!!!   DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR 
      27             : !!!!   OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE 
      28             : !!!!   OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
      29             : !!!!
      30             : !!!!   ACKNOWLEDGMENT
      31             : !!!!
      32             : !!!!   ParaMonte is an honor-ware and its currency is acknowledgment and citations.
      33             : !!!!   As per the ParaMonte library license agreement terms, if you use any parts of 
      34             : !!!!   this library for any purposes, kindly acknowledge the use of ParaMonte in your 
      35             : !!!!   work (education/research/industry/development/...) by citing the ParaMonte 
      36             : !!!!   library as described on this page:
      37             : !!!!
      38             : !!!!       https://github.com/cdslaborg/paramonte/blob/main/ACKNOWLEDGMENT.md
      39             : !!!!
      40             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      41             : !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
      42             : 
      43             : !> \brief This module contains the mathematical and programming constants.
      44             : !> \author Amir Shahmoradi
      45             : 
      46             : module Constants_mod
      47             : 
      48             :     use, intrinsic :: iso_fortran_env, only: real32, real64, int32, int64
      49             :     use, intrinsic :: iso_c_binding, only: CIK => c_int32_t, CRK => c_double
      50             : #if defined CFI_ENABLED
      51             :     use, intrinsic :: iso_c_binding, only: IK => c_int32_t, RK => c_double
      52             : #else
      53             :     use, intrinsic :: iso_fortran_env, only: IK => int32, RK => real64
      54             : #endif
      55             : 
      56             :     implicit none
      57             : 
      58             :     character(*), parameter :: MODULE_NAME = "@Constants_mod"
      59             : 
      60             :     ! Constants for computational accuracy
      61             : 
      62             :     integer     , parameter :: SPR = real32                                             !< @public single precision real kind
      63             :     integer     , parameter :: DPR = real64                                             !< @public double precision real kind
      64             :     integer     , parameter :: SPI = int32                                              !< @public single precision integer kind
      65             :     integer     , parameter :: DPI = int64                                              !< @public double precision integer kind
      66             :     integer     , parameter :: SPC = kind((1._SPR,1._SPR))                              !< @public single-precision complex kind
      67             :     integer     , parameter :: DPC = kind((1._DPR,1._DPR))                              !< @public double-precision complex kind
      68             :     integer     , parameter :: CK = kind((1._RK,1._RK))                                 !< @public complex kind
      69             :     integer     , parameter :: RKP = precision(1._RK)                                   !< @public real kind precision
      70             :     integer(IK) , parameter :: MAX_REC_LEN = 9999                                       !< @public maximum string record length
      71             : 
      72             :     ! Mathematical constants
      73             : 
      74             :     real(RK)    , parameter :: PI = 3.141592653589793238462643383279502884197_DPR       !< @public = acos(-1._RK) : The irrational number Pi.
      75             :     real(RK)    , parameter :: INVPI = 1._DPR / PI                                      !< @public = inverse Pi.
      76             :     real(RK)    , PARAMETER :: TWOPI = 6.283185307179586476925286766559005768394_DPR    !< @public 2*PI
      77             :     real(RK)    , parameter :: LN2 = log(2._RK)                                         !< @public Natural Log of 2 (= 0.693147180559945_RK).
      78             :     real(RK)    , parameter :: INVLN2 = 1._RK / LN2                                     !< @public Inverse of the natural Log of 2 (= 0.693147180559945_RK).
      79             :     real(RK)    , parameter :: LN10 = log(1.e1_RK)                                      !< @public Natural Log of 10 (= 2.302585092994046_RK).
      80             :     real(RK)    , parameter :: LNPI = log(PI)                                           !< @public ln(2pi) (= 1.144729885849400_RK)
      81             :     real(RK)    , parameter :: LN2PI = log(2._RK*PI)                                    !< @public ln(2pi) (= 1.837877066409345_RK)
      82             :     real(RK)    , parameter :: SQRT2 = sqrt(2._RK)                                      !< @public Square root of 2.
      83             :     real(RK)    , parameter :: NAPIER = exp(1._RK)                                      !< @public Napier number e.
      84             :     real(RK)    , parameter :: SQRTPI = sqrt(PI)                                        !< @public Square root of Pi.
      85             :     real(RK)    , parameter :: SQRT2PI = sqrt(2._RK*acos(-1._RK))                       !< @public Square root of 2Pi.
      86             :     real(RK)    , parameter :: INVSQRT2 = 1._RK / SQRT2PI                               !< @public Square root of 2.
      87             :     real(RK)    , parameter :: HALFLN2PI = 0.5_RK*LN2PI                                 !< @public ln(sqrt(2pi))
      88             :     real(RK)    , parameter :: INVSQRTPI = sqrt(INVPI)                                  !< @public = inverse of the square root of Pi.
      89             :     real(RK)    , parameter :: INVSQRT2PI = 1._RK / SQRT2PI                             !< @public 1/sqrt(2*Pi) (= 0.398942280401432_RK)
      90             :     real(RK)    , parameter :: LOGINVSQRT2PI = log(INVSQRT2PI)                          !< @public Log(1/sqrt(2Pi)), used in Gaussian distribution.
      91             :     real(RK)    , parameter :: SQRT_HALF_PI = sqrt(0.5_RK*PI)                           !< @public Square root of PI/2 (= 1.2533141373155_RK)
      92             :     real(RK)    , parameter :: LOG10NAPIER = log10(NAPIER)                              !< @public Log10 of Napier constant (= 0.434294481903259_RK).
      93             :     real(RK)    , parameter :: EPS_RK = epsilon(1._RK)                                  !< @public the smallest representable real increment (highest precision) by the machine
      94             :     real(RK)    , parameter :: SQRT_EPS_RK = sqrt(EPS_RK)                               !< @public the smallest representable real increment (highest precision) by the machine
      95             :     real(RK)    , parameter :: HUGE_RK = huge(1._RK)                                    !< @public largest number of kind RK
      96             :     real(RK)    , parameter :: TINY_RK = tiny(1._RK)                                    !< @public tiniest number of kind RK
      97             :     real(RK)    , parameter :: LOGHUGE_RK = log(HUGE_RK)                                !< @public log of the largest number of kind RK
      98             :     real(RK)    , parameter :: LOGTINY_RK = log(TINY_RK)                                !< @public log of the smallest number of kind RK
      99             :     real(RK)    , parameter :: POSINF_RK =  HUGE_RK / 1.e1_RK                           !< @public positive virtual infinite real. The division is done to avoid overflow in output
     100             :     real(RK)    , parameter :: NEGINF_RK = -POSINF_RK                                   !< @public negative virtual infinite real. Defined to avoid underflow/overflow and compatibility with Dynamic languages.
     101             :     real(RK)    , parameter :: LOGINF_RK =  log(POSINF_RK)                              !< @public represents the logarithm of the largest representable number
     102             :     real(RK)    , parameter :: NEGLOGINF_RK = -LOGINF_RK                                !< @public represents the logarithm of the smallest representable number
     103             :     integer(IK) , parameter :: HUGE_IK = huge(1_IK)                                     !< @public largest number of kind RK
     104             :     integer(IK) , parameter :: POSINF_IK =  HUGE_IK / 2_IK                              !< @public positive virtually-infinite integer. the division is done to avoid overflow in output
     105             :     integer(IK) , parameter :: NEGINF_IK = -POSINF_IK                                   !< @public negative virtually-infinite integer.
     106             :     real(RK)    , parameter :: LOGINF_IK =  log(real(POSINF_IK,RK))                     !< @public the natural logarithm of the positive virtually-infinite integer.
     107             : 
     108             :     real(RK)    , parameter :: NULL_RK = -HUGE_RK                                       !< @public the value used to represent unassigned variables.
     109             :     integer(IK) , parameter :: NULL_IK = -HUGE_IK
     110             :     character(1), parameter :: NULL_SK = achar(30)                                      ! This must remain a single character as it is assumed to be so in multiple routines: Record separator
     111             : 
     112             :     character(1), parameter :: NLC = achar(10)                                          ! the New Line Character
     113             :     character(1), parameter :: TAB = achar(9)                                           ! the TAB Character
     114             :     character(*), parameter :: UNDEFINED = "UNDEFINED"
     115             : 
     116             :     ! null values
     117             : 
     118             :     type, private  :: NullType
     119             :         real(RK)     :: RK = NULL_RK
     120             :         integer(IK)  :: IK = NULL_IK
     121             :         character(1) :: SK = NULL_SK
     122             :     end type NullType
     123             :   
     124             :     type(NullType), protected :: NullVal
     125             : 
     126             :     ! Physical constants
     127             : 
     128             :     real(RK), parameter :: ERG2KEV = 6.241509125883258e8_RK                      !< @public 1 (erg) = ERG2KEV (keV)
     129             :     real(RK), parameter :: KEV2ERG = 1.60217662080000e-9_RK                      !< @public 1 (keV) = KEV2ERG (erg)
     130             :     real(RK), parameter :: LOG_ERG2KEV = log(ERG2KEV)
     131             :     real(RK), parameter :: LOG_KEV2ERG = log(KEV2ERG)
     132             : 
     133             :     ! Cosmological constants
     134             : 
     135             :     !real(RK), parameter :: LIGHT_SPEED = 3.e5_RK                                ! LIGHT_SPEED is the speed of light (Km/s).
     136             :     !real(RK), parameter :: HUBBLE_TIME_GYRS = 13.8_RK                               ! hubble time (liddle 2003, page 57) in units of gyrs
     137             :     !real(RK), parameter :: HUBBLE_CONST = 7.1e1_RK                              ! HUBBLE_CONST is the Hubble constant in units of km/s/MPc.
     138             :     !real(RK), parameter :: LS2HC = LIGHT_SPEED / HUBBLE_CONST                   ! the speed of light in units of km/s divided by the Hubble constant.
     139             :     !real(RK), parameter :: MPC2CM = 3.09e24_RK                                  ! 1 Mega Parsec = MPC2CM centimeters.
     140             :     !real(RK), parameter :: LOG10MPC2CMSQ4PI = log10(4._RK*PI) + 2*log10(MPC2CM) ! Log10(MPC2CM centimeters.
     141             :     !real(RK), parameter :: OMEGA_DE = 0.7_RK                                    ! Dark Energy density.
     142             :     !real(RK), parameter :: OMEGA_DM = 0.3_RK                                    ! Dark Matter density.
     143             : 
     144             :     character(len=1), parameter :: CARRIAGE_RETURN = achar(13)
     145             :     character(len=1), parameter :: ESCAPE = achar(27)
     146             :     character(len=1), parameter :: ESC = achar(27)
     147             :     character(len=1), parameter :: CLOCK_TICK(4) = [ "|" , "/" , "-" , "\" ]
     148             : 
     149             :     interface getPosInf
     150             :         module procedure :: getPosInf_RK
     151             :     end interface getPosInf
     152             : 
     153             :     interface getNegInf
     154             :         module procedure :: getNegInf_RK
     155             :     end interface getNegInf
     156             : 
     157             :     ! file extentions
     158             : 
     159             :     type, private       :: FileType_type
     160             :         character(6)    :: binary  = "binary"
     161             :         character(6)    :: matlab  = "MATLAB"
     162             :         character(6)    :: python  = "Python"
     163             :         character(5)    :: julia   = "Julia"
     164             :         character(5)    :: ascii   = "ASCII"
     165             :         character(1)    :: rlang   = "R"
     166             :     end type FileType_type
     167             : 
     168             :     type, private       :: FileExt_type
     169             :         character(4)    :: binary  = ".bin"
     170             :         character(2)    :: matlab  = ".m"
     171             :         character(3)    :: python  = ".py"
     172             :         character(3)    :: julia   = ".jl"
     173             :         character(4)    :: ascii   = ".txt"
     174             :         character(2)    :: r       = ".R"
     175             :     end type FileExt_type   
     176             : 
     177             :     type(FileExt_type)  , parameter :: FILE_EXT = FileExt_type()
     178             :     type(FileType_type) , parameter :: FILE_TYPE = FileType_type()
     179             : 
     180             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     181             : ! ParaMonte methods
     182             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     183             : 
     184             :     type, private       :: ParaMonteSamplingMethod_type
     185             :         integer(IK)     :: count = 4_IK
     186             :         character(8)    :: ParaDRAM = "ParaDRAM"
     187             :         character(8)    :: ParaDISE = "ParaDISE"
     188             :         character(8)    :: ParaHDMC = "ParaHDMC"
     189             :         character(8)    :: ParaTemp = "ParaTemp"
     190             :         character(8)    :: ParaNest = "ParaNest"
     191             :     end type ParaMonteSamplingMethod_type
     192             : 
     193             :     type(ParaMonteSamplingMethod_type), parameter :: PMSM = ParaMonteSamplingMethod_type()
     194             : 
     195             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     196             : 
     197             : contains
     198             : 
     199             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     200             : 
     201           6 :     pure function getPosInf_RK() result(posInf)
     202             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     203             :         !DEC$ ATTRIBUTES DLLEXPORT :: getPosInf_RK
     204             : #endif
     205             :         !> \brief Return IEEE positive infinity.
     206             :         !> @param[out] posInf The positive infinity.
     207             :         use, intrinsic :: ieee_arithmetic, only: ieee_value, ieee_positive_inf
     208             :         implicit none
     209             :         real(RK) :: posInf
     210           3 :         posInf = ieee_value(0._RK, ieee_positive_inf)
     211           3 :     end function getPosInf_RK
     212             : 
     213             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     214             : 
     215           3 :     pure function getNegInf_RK() result(negInf)
     216             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     217             :         !DEC$ ATTRIBUTES DLLEXPORT :: getNegInf_RK
     218             : #endif
     219             :         !> \brief Return IEEE negative infinity.
     220             :         !> @param[out] posInf The negative infinity.
     221           3 :         use, intrinsic :: ieee_arithmetic, only: ieee_value, ieee_negative_inf
     222             :         implicit none
     223             :         real(RK) :: negInf
     224           3 :         negInf = ieee_value(0._RK, ieee_negative_inf)
     225           3 :     end function getNegInf_RK
     226             : 
     227             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     228             : 
     229             : end module Constants_mod ! LCOV_EXCL_LINE

ParaMonte: Plain Powerful Parallel Monte Carlo Library 
The Computational Data Science Lab
© Copyright 2012 - 2021