https://www.cdslab.org/paramonte/fortran/2
Current view: top level - main - pm_sampling.F90 (source / functions) Hit Total Coverage
Test: ParaMonte 2.0.0 :: Serial Fortran - Code Coverage Report Lines: 1 3 33.3 %
Date: 2024-04-08 03:18:57 Functions: 1 3 33.3 %
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 module contains procedures and generic interfaces for the ParaMonte library sampler routines.
      19             : !>
      20             : !>  \test
      21             : !>  [test_pm_sampling](@ref test_pm_sampling)<br>
      22             : !>
      23             : !>  \finmain
      24             : !>
      25             : !>  \author
      26             : !>  \AmirShahmoradi, Monday 00:01 AM, January 1, 2018, Institute for Computational Engineering and Sciences, University of Texas Austin
      27             : 
      28             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      29             : 
      30             : module pm_sampling
      31             : 
      32             :     use pm_err, only: err_type
      33             :     use pm_str, only: getCharSeq
      34             :     use, intrinsic :: iso_c_binding, only: c_char, c_funptr, c_null_char, c_f_procpointer
      35             :     use pm_kind, only: SK, IK, LK, RKL, RKD, RKH, RKALL, SKC => SK ! LCOV_EXCL_LINE
      36             :     use pm_sysPath, only: PATHLEN => MAX_LEN_FILE_PATH
      37             : 
      38             :     implicit none
      39             : 
      40             :     character(*, SK), parameter :: MODULE_NAME = "@pm_sampling"
      41             : 
      42             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      43             : 
      44             :     !>  \brief
      45             :     !>  This is a derived type for constructing objects containing
      46             :     !>  the **optional** simulation properties of the ParaMonte library explorers and samplers.<br>
      47             :     !>
      48             :     !>  \details
      49             :     !>  Objects of this derived type are not meant to be used directly by the end users.<br>
      50             :     !>  Instead use one of the appropriate derived types:<br>
      51             :     !>  <ol>
      52             :     !>      <li>    [paradram_type](@ref pm_sampling::paradram_type)
      53             :     !>      <li>    [paranest_type](@ref pm_sampling::paranest_type)
      54             :     !>      <li>    [paradise_type](@ref pm_sampling::paradise_type)
      55             :     !>  </ol>
      56             :     !>
      57             :     !>  \note
      58             :     !>  For more information on any of the components of this derived type,
      59             :     !>  see the output `_report.txt` files from the corresponding sampler simulations.<br>
      60             :     !>  Alternatively, see the ParaMonte library cross-language sampler documentation
      61             :     !>  for detailed descriptions of the input simulation specifications at:<br>
      62             :     !>  <ol>
      63             :     !>      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
      64             :     !>  </ol>
      65             :     !>
      66             :     !>  \see
      67             :     !>  [paradram_type](@ref pm_sampling::paradram_type)<br>
      68             :     !>  [paranest_type](@ref pm_sampling::paranest_type)<br>
      69             :     !>  [paradise_type](@ref pm_sampling::paradise_type)<br>
      70             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
      71             :     !>
      72             :     !>  \finmain{sampler_type}
      73             :     !>
      74             :     !>  \author
      75             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
      76             :     type                                        :: sampler_type
      77             :         character(:,SKC)        , allocatable   :: description                      !<  \specdram{description}
      78             :         character(:,SKC)        , allocatable   :: domain                           !<  \specdram{domain}
      79             :         character(:,SKC)        , allocatable   :: domainAxisName(:)                !<  \specdram{domainaxisname}
      80             :         real(RKH)               , allocatable   :: domainBallCenter(:)              !<  \specdram{domainballcenter}
      81             :         real(RKH)               , allocatable   :: domainBallCorMat(:,:)            !<  \specdram{domainballcormat}
      82             :         real(RKH)               , allocatable   :: domainBallCovMat(:,:)            !<  \specdram{domainballcovmat}
      83             :         real(RKH)               , allocatable   :: domainBallStdVec(:)              !<  \specdram{domainballstdvec}
      84             :         real(RKH)               , allocatable   :: domainCubeLimitLower(:)          !<  \specdram{domaincubelimitlower}
      85             :         real(RKH)               , allocatable   :: domainCubeLimitUpper(:)          !<  \specdram{domaincubelimitupper}
      86             :         integer(IK)             , allocatable   :: domainErrCount                   !<  \specdram{domainerrcount}
      87             :         integer(IK)             , allocatable   :: domainErrCountMax                !<  \specdram{domainerrcountmax}
      88             :         character(:,SKC)        , allocatable   :: inputFile                        !<  \public The input scalar `character` of default kind \SK which be either,
      89             :                                                                                     !!  <ol>
      90             :                                                                                     !!      <li>    The path to an external input `namelist` file containing the simulation settings.<br>
      91             :                                                                                     !!      <li>    A string containing the simulation settings in `namelist` format.<br>
      92             :                                                                                     !!              The following rules apply to variable assignments in namelist-style input:<br>
      93             :                                                                                     !!              <ol>
      94             :                                                                                     !!                  <li>    Comments must begin with an exclamation mark `!`.
      95             :                                                                                     !!                  <li>    Comments can appear anywhere on an empty line or, after a variable assignment
      96             :                                                                                     !!                          (but not in the middle of a variable assignment whether in single or multiple lines).
      97             :                                                                                     !!                  <li>    All variable assignments are optional and can be commented out. In such cases, appropriate default values will be used.
      98             :                                                                                     !!                  <li>    All variable assignments must appear within a `namelist` group which starts with either
      99             :                                                                                     !!                          <ol>
     100             :                                                                                     !!                              <li>    the (case-INsensitive) keyword `&paradise` to be used for specifying [paradram_type](@ref pm_sampling::paradram_type) simulation specifications.<br>
     101             :                                                                                     !!                              <li>    the (case-INsensitive) keyword `&paradram` to be used for specifying [paradram_type](@ref pm_sampling::paradram_type) simulation specifications.<br>
     102             :                                                                                     !!                              <li>    the (case-INsensitive) keyword `&paranest` to be used for specifying [paradram_type](@ref pm_sampling::paradram_type) simulation specifications.<br>
     103             :                                                                                     !!                          </ol>
     104             :                                                                                     !!                          and ends with a forward slash `/` on the last line of the namelist group.
     105             :                                                                                     !!                  <li>    Every time a sampler opens an input file, it will first search for the namelist group with the corresponding sampler name.<br>
     106             :                                                                                     !!                  <li>    The order of the input variables in the namelist groups is irrelevant and unimportant.
     107             :                                                                                     !!                  <li>    Variables can be defined multiple times, but **only the last definition** will be considered as input.
     108             :                                                                                     !!                  <li>    All variable names are case-insensitive.
     109             :                                                                                     !!                          However, for clarity, this software follows the camelCase code-writing practice for simulation specification names.
     110             :                                                                                     !!                  <li>    String values must be enclosed with either single or double quotation marks.
     111             :                                                                                     !!                  <li>    Logical values are case-insensitive and can be either<br>
     112             :                                                                                     !!                          <ol>
     113             :                                                                                     !!                              <li>    `.true.`, `true`, or `t` (all case-insensitive) representing the truth, or,
     114             :                                                                                     !!                              <li>    `.false.`, `false`, or `f` (all case-insensitive) representing the falsehood.
     115             :                                                                                     !!                          </ol>
     116             :                                                                                     !!                  <li>    All vectors and arrays specified within the input file start with a lower-bound of `1` as opposed to `0` in the C family of languages.<br>
     117             :                                                                                     !!                          This follows the convention of the majority of science-oriented programming languages: Fortran, Julia, Mathematica, MATLAB, and R.
     118             :                                                                                     !!                  <li>    Specifying multiple values is as simple as listing the elements as comma or space -separated values.
     119             :                                                                                     !!                  <li>    Individual vector or element values can be specified by an explicit index or index range:<br>
     120             :                                                                                     !!                          \code{.F90}
     121             :                                                                                     !!                              &paranest
     122             :                                                                                     !!                                  domainAxisName(1 : 2) = "sampleVariable1", "sampleVariable2"
     123             :                                                                                     !!                                  domainAxisName(3) = "sampleVariable3"
     124             :                                                                                     !!                              /
     125             :                                                                                     !!                          \endcode
     126             :                                                                                     !!                  <li>    If a sequence of elements of a vector or array are equal,
     127             :                                                                                     !!                          the assignment can be succinctly expressed via the following convention:<br>
     128             :                                                                                     !!                          \code{.F90}
     129             :                                                                                     !!                              &paradram
     130             :                                                                                     !!                                  proposalStart = 4 * 0. ! assign value `0.` to the first four elements of proposalStart.
     131             :                                                                                     !!                              /
     132             :                                                                                     !!                          \endcode
     133             :                                                                                     !!                  <li>    Matrix assignments follow the [column-major](https://en.wikipedia.org/wiki/Row-_and_column-major_order) storage scheme in memory.<br>
     134             :                                                                                     !!                          \code{.F90}
     135             :                                                                                     !!                              &paradram
     136             :                                                                                     !!                                  proposalCorMat = 1, 0, 0, 0
     137             :                                                                                     !!                                                 , 0, 1, 0, 0
     138             :                                                                                     !!                                                 , 0, 0, 1, 0
     139             :                                                                                     !!                                                 , 0, 0, 0, 1
     140             :                                                                                     !!                              /
     141             :                                                                                     !!                          \endcode
     142             :                                                                                     !!              </ol>
     143             :                                                                                     !!  </ol>
     144             :                                                                                     !!  For detailed descriptions of the input simulation settings see,
     145             :                                                                                     !!  <ol>
     146             :                                                                                     !!      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
     147             :                                                                                     !!      <li>    [ParaDISE simulation options](\pmdoc_usage_sampling/paradise/specifications/)
     148             :                                                                                     !!      <li>    [ParaNest simulation options](\pmdoc_usage_sampling/paranest/specifications/)
     149             :                                                                                     !!  </ol>
     150             :                                                                                     !!  (**optional**. If missing, the default simulation settings will be used.)
     151             :        !logical(LK)             , allocatable   :: inputFileHasPriority             !<
     152             :         character(:,SKC)        , allocatable   :: outputChainFileFormat            !<  \specdram{outputchainfileformat}
     153             :         integer(IK)             , allocatable   :: outputColumnWidth                !<  \specdram{outputcolumnwidth}
     154             :         character(:,SKC)        , allocatable   :: outputFileName                   !<  \specdram{outputfilename}
     155             :         character(:,SKC)        , allocatable   :: outputStatus                     !<  \specdram{outputstatus}
     156             :         integer(IK)             , allocatable   :: outputPrecision                  !<  \specdram{outputprecision}
     157             :         integer(IK)             , allocatable   :: outputReportPeriod               !<  \specdram{outputreportperiod}
     158             :         character(:,SKC)        , allocatable   :: outputRestartFileFormat          !<  \specdram{outputrestartfileformat}
     159             :         integer(IK)             , allocatable   :: outputSampleSize                 !<  \specdram{outputsamplesize}
     160             :         character(:,SKC)        , allocatable   :: outputSeparator                  !<  \specdram{outputseparator}
     161             :         character(:,SKC)        , allocatable   :: outputSplashMode                 !<  \specdram{outputsplashmode}
     162             :         character(:,SKC)        , allocatable   :: parallelism                      !<  \specdram{parallelism}
     163             :         logical(LK)             , allocatable   :: parallelismMpiFinalizeEnabled    !<  \specdram{parallelismmpifinalizeenabled}
     164             :         integer(IK)             , allocatable   :: parallelismNumThread             !<  \specdram{parallelismnumthread}
     165             :         integer(IK)             , allocatable   :: randomSeed                       !<  \specdram{randomseed}
     166             :         character(:,SKC)        , allocatable   :: sysInfoFilePath                  !<  \specdram{sysinfofilepath}
     167             :         real(RKH)               , allocatable   :: targetAcceptanceRate(:)          !<  \specdram{targetacceptancerate}
     168             :     end type
     169             : 
     170             :     !>  \brief
     171             :     !>  This is a derived type for constructing objects containing
     172             :     !>  the **optional** simulation properties of the ParaMonte library MCMC explorers and samplers.<br>
     173             :     !>
     174             :     !>  \details
     175             :     !>  Objects of this derived type are not meant to be used directly by the end users.<br>
     176             :     !>  Instead use one of the appropriate derived types:<br>
     177             :     !>  <ol>
     178             :     !>      <li>    [paradram_type](@ref pm_sampling::paradram_type)
     179             :     !>      <li>    [paranest_type](@ref pm_sampling::paranest_type)
     180             :     !>      <li>    [paradise_type](@ref pm_sampling::paradise_type)
     181             :     !>  </ol>
     182             :     !>
     183             :     !>  \note
     184             :     !>  For more information on any of the components of this derived type,
     185             :     !>  see the output `_report.txt` files from the corresponding sampler simulations.<br>
     186             :     !>  Alternatively, see the ParaMonte library cross-language sampler documentation
     187             :     !>  for detailed descriptions of the input simulation specifications at:<br>
     188             :     !>  <ol>
     189             :     !>      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
     190             :     !>      <li>    [ParaDISE simulation options](\pmdoc_usage_sampling/paradise/specifications/)
     191             :     !>      <li>    [ParaNest simulation options](\pmdoc_usage_sampling/paranest/specifications/)
     192             :     !>  </ol>
     193             :     !>
     194             :     !>  \see
     195             :     !>  [paradram_type](@ref pm_sampling::paradram_type)<br>
     196             :     !>  [paranest_type](@ref pm_sampling::paranest_type)<br>
     197             :     !>  [paradise_type](@ref pm_sampling::paradise_type)<br>
     198             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
     199             :     !>
     200             :     !>  \finmain{paramcmc_type}
     201             :     !>
     202             :     !>  \author
     203             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     204             :     type, extends(sampler_type)                 :: paramcmc_type
     205             :         integer(IK)             , allocatable   :: outputChainSize                      !<  \specdram{outputchainsize}
     206             :         integer(IK)             , allocatable   :: outputSampleRefinementCount          !<  \specdram{outputsamplerefinementcount}
     207             :         character(:,SKC)        , allocatable   :: outputSampleRefinementMethod         !<  \specdram{outputsamplerefinementmethod}
     208             :         character(:,SKC)        , allocatable   :: proposal                             !<  \specdram{proposal}
     209             :         real(RKH)               , allocatable   :: proposalCorMat(:,:)                  !<  \specdram{proposalcormat}
     210             :         real(RKH)               , allocatable   :: proposalCovMat(:,:)                  !<  \specdram{proposalcovmat}
     211             :         character(:,SKC)        , allocatable   :: proposalScaleFactor                  !<  \specdram{proposalscalefactor}
     212             :         real(RKH)               , allocatable   :: proposalStart(:)                     !<  \specdram{proposalstart}
     213             :         real(RKH)               , allocatable   :: proposalStartDomainCubeLimitLower(:) !<  \specdram{proposalstartdomaincubelimitlower}
     214             :         real(RKH)               , allocatable   :: proposalStartDomainCubeLimitUpper(:) !<  \specdram{proposalstartdomaincubelimitupper}
     215             :         logical(LK)             , allocatable   :: proposalStartRandomized              !<  \specdram{proposalstartrandomized}
     216             :         real(RKH)               , allocatable   :: proposalStdVec(:)                    !<  \specdram{proposalstdvec}
     217             :     end type
     218             : 
     219             :     !>  \brief
     220             :     !>  This is a derived type for constructing objects containing the **optional** simulation
     221             :     !>  properties of the ParaMonte library Delayed-Rejection Adaptive Metropolis (DRAM) MCMC explorers and samplers.<br>
     222             :     !>
     223             :     !>  \note
     224             :     !>  For more information on any of the components of this derived type,
     225             :     !>  see the output `_report.txt` files from the corresponding sampler simulations.<br>
     226             :     !>  Alternatively, see the ParaMonte library cross-language sampler documentation
     227             :     !>  for detailed descriptions of the input simulation specifications at:<br>
     228             :     !>  <ol>
     229             :     !>      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
     230             :     !>      <li>    [ParaDISE simulation options](\pmdoc_usage_sampling/paradise/specifications/) (work in progress)
     231             :     !>      <li>    [ParaNest simulation options](\pmdoc_usage_sampling/paranest/specifications/) (work in progress)
     232             :     !>  </ol>
     233             :     !>
     234             :     !>  \see
     235             :     !>  [paradram_type](@ref pm_sampling::paradram_type)<br>
     236             :     !>  [paranest_type](@ref pm_sampling::paranest_type)<br>
     237             :     !>  [paradise_type](@ref pm_sampling::paradise_type)<br>
     238             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
     239             :     !>
     240             :     !>  \finmain{paradram_type}
     241             :     !>
     242             :     !>  \author
     243             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     244             :     type, extends(paramcmc_type)                :: paradram_type
     245             :         real(RKH)               , allocatable   :: burninAdaptationMeasure                  !<  \specdram{burninadaptationmeasure}
     246             :         integer(IK)             , allocatable   :: proposalAdaptationCount                  !<  \specdram{proposaladaptationcount}
     247             :         integer(IK)             , allocatable   :: proposalAdaptationCountGreedy            !<  \specdram{proposaladaptationcountgreedy}
     248             :         integer(IK)             , allocatable   :: proposalAdaptationPeriod                 !<  \specdram{proposaladaptationperiod}
     249             :         integer(IK)             , allocatable   :: proposalDelayedRejectionCount            !<  \specdram{proposaldelayedrejectioncount}
     250             :         real(RKH)               , allocatable   :: proposalDelayedRejectionScaleFactor(:)   !<  \specdram{proposaldelayedrejectionscalefactor}
     251             :     end type
     252             : 
     253             :     !>  \brief
     254             :     !>  This is a derived type for constructing objects containing the **optional** simulation
     255             :     !>  properties of the ParaMonte library DRAM-enhanced MCMC explorers and samplers.<br>
     256             :     !>
     257             :     !>  \note
     258             :     !>  For more information on any of the components of this derived type,
     259             :     !>  see the output `_report.txt` files from the corresponding sampler simulations.<br>
     260             :     !>  Alternatively, see the ParaMonte library cross-language sampler documentation
     261             :     !>  for detailed descriptions of the input simulation specifications at:<br>
     262             :     !>  <ol>
     263             :     !>      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
     264             :     !>      <li>    [ParaDISE simulation options](\pmdoc_usage_sampling/paradise/specifications/) (work in progress)
     265             :     !>      <li>    [ParaNest simulation options](\pmdoc_usage_sampling/paranest/specifications/) (work in progress)
     266             :     !>  </ol>
     267             :     !>
     268             :     !>  \see
     269             :     !>  [paradram_type](@ref pm_sampling::paradram_type)<br>
     270             :     !>  [paranest_type](@ref pm_sampling::paranest_type)<br>
     271             :     !>  [paradise_type](@ref pm_sampling::paradise_type)<br>
     272             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
     273             :     !>
     274             :     !>  \finmain{paradise_type}
     275             :     !>
     276             :     !>  \author
     277             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     278             :     type, extends(paradram_type)                :: paradise_type
     279             :     end type
     280             : 
     281             :     !>  \brief
     282             :     !>  This is a derived type for constructing objects containing the **optional** simulation
     283             :     !>  properties of the ParaMonte library Nested explorers and samplers.<br>
     284             :     !>
     285             :     !>  \note
     286             :     !>  For more information on any of the components of this derived type,
     287             :     !>  see the output `_report.txt` files from the corresponding sampler simulations.<br>
     288             :     !>  Alternatively, see the ParaMonte library cross-language sampler documentation
     289             :     !>  for detailed descriptions of the input simulation specifications at:<br>
     290             :     !>  <ol>
     291             :     !>      <li>    [ParaDRAM simulation options](\pmdoc_usage_sampling/paradram/specifications/)
     292             :     !>      <li>    [ParaDISE simulation options](\pmdoc_usage_sampling/paradise/specifications/) (work in progress)
     293             :     !>      <li>    [ParaNest simulation options](\pmdoc_usage_sampling/paranest/specifications/) (work in progress)
     294             :     !>  </ol>
     295             :     !>
     296             :     !>  \see
     297             :     !>  [paradram_type](@ref pm_sampling::paradram_type)<br>
     298             :     !>  [paranest_type](@ref pm_sampling::paranest_type)<br>
     299             :     !>  [paradise_type](@ref pm_sampling::paradise_type)<br>
     300             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
     301             :     !>  [getErrSampling](@ref pm_sampling::getErrSampling)<br>
     302             :     !>
     303             :     !>  \finmain{paradise_type}
     304             :     !>
     305             :     !>  \author
     306             :     !>  \AmirShahmoradi, September 1, 2017, 12:00 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin
     307             :     type, extends(sampler_type)                 :: paranest_type
     308             :         integer(IK)             , allocatable   :: domainPartitionAdaptationCount
     309             :         integer(IK)             , allocatable   :: domainPartitionAdaptationPeriod
     310             :         logical(LK)             , allocatable   :: domainPartitionBiasCorrectionEnabled
     311             :         integer(IK)             , allocatable   :: domainPartitionCountMax
     312             :         real(RKH)               , allocatable   :: domainPartitionFactorExpansion
     313             :         real(RKH)               , allocatable   :: domainPartitionFactorShrinkage
     314             :         integer(IK)             , allocatable   :: domainPartitionKmeansClusterCountMax
     315             :         integer(IK)             , allocatable   :: domainPartitionKmeansClusterSizeMin
     316             :         logical(LK)             , allocatable   :: domainPartitionKmeansNormalizationEnabled
     317             :         integer(IK)             , allocatable   :: domainPartitionKmeansNumFailMax
     318             :         integer(IK)             , allocatable   :: domainPartitionKmeansNumRecursionMax
     319             :         integer(IK)             , allocatable   :: domainPartitionKmeansNumTry
     320             :         integer(IK)             , allocatable   :: domainPartitionKvolumeNumRecursionMax
     321             :         real(RKH)               , allocatable   :: domainPartitionKvolumeWeightExponent
     322             :         character(:,SKC)        , allocatable   :: domainPartitionMethod
     323             :         character(:,SKC)        , allocatable   :: domainPartitionObject
     324             :         logical(LK)             , allocatable   :: domainPartitionOptimizationScaleEnabled
     325             :         logical(LK)             , allocatable   :: domainPartitionOptimizationShapeEnabled
     326             :         logical(LK)             , allocatable   :: domainPartitionOptimizationShapeScaleEnabled
     327             :         real(RKH)               , allocatable   :: domainPartitionScaleFactor
     328             :         character(:,SKC)        , allocatable   :: domainSampler
     329             :         integer(IK)             , allocatable   :: liveSampleSize
     330             :         real(RKH)               , allocatable   :: tolerance
     331             :     end type
     332             : 
     333             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     334             : 
     335             :     !   This abstract interface defines the set of objective function interfaces
     336             :     !   that can be passed to the ParaMonte samplers (e.g., ParaDRAM, ParaNest, ...).
     337             :     !   These interfaces only have internal library specific use cases.
     338             :     !>  \cond excluded
     339             :     abstract interface
     340             : 
     341             : #if OMP_ENABLED && (MATLAB_ENABLED || PYTHON_ENABLED || R_ENABLED)
     342             : 
     343             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     344             : 
     345             : #if RK5_ENABLED
     346             :     function getLogFunc_proc_RK5(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     347             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     348             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK5
     349             : #endif
     350             :         use pm_kind, only: RKD, RKC => RK5
     351             :         real(RKC), intent(inout), contiguous :: logFuncState(0:,:)
     352             :         real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     353             :         real(RKC) :: mold
     354             :     end function
     355             : #endif
     356             : 
     357             : #if RK4_ENABLED
     358             :     function getLogFunc_proc_RK4(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     359             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     360             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK4
     361             : #endif
     362             :         use pm_kind, only: RKD, RKC => RK4
     363             :         real(RKC), intent(inout), contiguous :: logFuncState(0:,:)
     364             :         real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     365             :         real(RKC) :: mold
     366             :     end function
     367             : #endif
     368             : 
     369             : #if RK3_ENABLED
     370             :     function getLogFunc_proc_RK3(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     371             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     372             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK3
     373             : #endif
     374             :         use pm_kind, only: RKD, RKC => RK3
     375             :         real(RKC), intent(inout), contiguous :: logFuncState(0:,:)
     376             :         real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     377             :         real(RKC) :: mold
     378             :     end function
     379             : #endif
     380             : 
     381             : #if RK2_ENABLED
     382             :     function getLogFunc_proc_RK2(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     383             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     384             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK2
     385             : #endif
     386             :         use pm_kind, only: RKD, RKC => RK2
     387             :         real(RKC), intent(inout), contiguous :: logFuncState(0:,:)
     388             :         real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     389             :         real(RKC) :: mold
     390             :     end function
     391             : #endif
     392             : 
     393             : #if RK1_ENABLED
     394             :     function getLogFunc_proc_RK1(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     395             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     396             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK1
     397             : #endif
     398             :         use pm_kind, only: RKD, RKC => RK1
     399             :         real(RKC), intent(inout), contiguous :: logFuncState(0:,:)
     400             :         real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     401             :         real(RKC) :: mold
     402             :     end function
     403             : #endif
     404             : 
     405             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     406             : 
     407             : #else
     408             : 
     409             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     410             : 
     411             : #if RK5_ENABLED
     412             :     function getLogFunc_proc_RK5(state) result(logFunc)
     413             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     414             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK5
     415             : #endif
     416             :         use pm_kind, only: RKC => RK5
     417             :         real(RKC), intent(in), contiguous :: state(:)
     418             :         real(RKC) :: logFunc
     419             :     end function
     420             : #endif
     421             : 
     422             : #if RK4_ENABLED
     423             :     function getLogFunc_proc_RK4(state) result(logFunc)
     424             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     425             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK4
     426             : #endif
     427             :         use pm_kind, only: RKC => RK4
     428             :         real(RKC), intent(in), contiguous :: state(:)
     429             :         real(RKC) :: logFunc
     430             :     end function
     431             : #endif
     432             : 
     433             : #if RK3_ENABLED
     434             :     function getLogFunc_proc_RK3(state) result(logFunc)
     435             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     436             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK3
     437             : #endif
     438             :         use pm_kind, only: RKC => RK3
     439             :         real(RKC), intent(in), contiguous :: state(:)
     440             :         real(RKC) :: logFunc
     441             :     end function
     442             : #endif
     443             : 
     444             : #if RK2_ENABLED
     445             :     function getLogFunc_proc_RK2(state) result(logFunc)
     446             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     447             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK2
     448             : #endif
     449             :         use pm_kind, only: RKC => RK2
     450             :         real(RKC), intent(in), contiguous :: state(:)
     451             :         real(RKC) :: logFunc
     452             :     end function
     453             : #endif
     454             : 
     455             : #if RK1_ENABLED
     456             :     function getLogFunc_proc_RK1(state) result(logFunc)
     457             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     458             :         !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK1
     459             : #endif
     460             :         use pm_kind, only: RKC => RK1
     461             :         real(RKC), intent(in), contiguous :: state(:)
     462             :         real(RKC) :: logFunc
     463             :     end function
     464             : #endif
     465             : 
     466             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     467             : 
     468             : #endif
     469             : 
     470             :     end interface
     471             :     !>  \endcond excluded
     472             : 
     473             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     474             : 
     475             :     !>  \brief
     476             :     !>  Generate and return `.true.` if the procedure fails to fully accomplish the task of generating a
     477             :     !>  Monte Carlo sample of the specified input mathematical objective function, otherwise, return `.false.`.
     478             :     !>
     479             :     !>  \details
     480             :     !>  This generic interface is the entry point to all ParaMonte Monte Carlo samplers of mathematical density functions.<br>
     481             :     !>  Although the procedures of this generic interface return a single scalar object of type [err_type](@ref pm_err::err_type),
     482             :     !>  the procedures generate massive amounts of information about each simulation which are stored in appropriate external hard drive files.<br>
     483             :     !>
     484             :     !>  \param[in]  sampler     :   The input scalar object that can be of,<br>
     485             :     !>                              <ol>
     486             :     !>                                  <li>    type [paradram_type](@ref pm_sampling::paradram_type) indicating the use of the ParaDRAM sampler,<br>
     487             :     !>                                  <li>    type [paradise_type](@ref pm_sampling::paradise_type) indicating the use of the ParaDISE sampler,<br>
     488             :     !>                                  <li>    type [paranest_type](@ref pm_sampling::paranest_type) indicating the use of the ParaNest sampler,<br>
     489             :     !>                              </ol>
     490             :     !>  \param[in]  getLogFunc  :   The input user-specified `function` that takes an input `contiguous` vector of shape `state(1:ndim)` of type `real` of kind \RKALL.<br>
     491             :     !>                              On input, `state` represents a state (point) from within the domain of the user-specified target density function whose function value must be returned.<br>
     492             :     !>                              On output, the user-specified procedure `getLogFunc()` must return the function value corresponding to the input `state(1:ndim)`.<br>
     493             :     !>                              The following illustrate the generic interface of `getLogFunc(state)`,
     494             :     !>                              \code{.F90}
     495             :     !>                                  function getLogFunc(state) result(logFunc)
     496             :     !>                                      real(RKC), intent(in), contiguous :: state(:) ! (1 : ndim)
     497             :     !>                                      real(RKC) :: logFunc
     498             :     !>                                  end function
     499             :     !>                              \endcode
     500             :     !>                              where the condition `all(shape(state) == [ndim])` holds where `ndim` is the number of dimensions of the target density function.<br>
     501             :     !>                              **If you are an end user, this is where you stop reading about this input argument**.<br>
     502             :     !>                              If you are a developer, continue reading about the alternative internal interfaces of this argument for dynamic (MATLAB/Python/R) languages.<br>
     503             :     !>                              Setting the preprocessing flag `CFI_ENABLED` to a non-zero value at the time of library compilation
     504             :     !>                              activates the C-Fortran Interoperability (CFI) interface of this generic interface.<br>
     505             :     !>                              This requires `getLogFunc()` to have the following interfaces:<br>
     506             :     !>                              <ol>
     507             :     !>                                  <li>    When the library is built for OpenMP parallelism (i.e., `OMP_ENABLED=1` holds) for MATLAB/Python/R programming languages,
     508             :     !>                                          the input `getLogFunc()` must have the following interface:<br>
     509             :     !>                                          \code{.F90}
     510             :     !>                                              function getLogFunc(logFuncState, ndim, njob, avgTimePerFunCall, avgCommPerFunCall) result(mold)
     511             :     !>                                                  integer(IK), value :: ndim, njob
     512             :     !>                                                  real(RKC), intent(inout) :: logFuncState(ndim + 1, njob)
     513             :     !>                                                  real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
     514             :     !>                                                  real(RKC) :: mold
     515             :     !>                                              end function
     516             :     !>                                          \endcode
     517             :     !>                                          where,
     518             :     !>                                          <ol>
     519             :     !>                                              <li>    the variable `ndim` is the number of dimensions of the target density function and
     520             :     !>                                                      the variable `njob` is the number of states for which the function must be evaluated concurrently.
     521             :     !>                                              <li>    the constant \RKD refers to the `real` kind `double precision` supported by the compiler.
     522             :     !>                                              <li>    the output `mold` merely dictates the `real` kind of the interface and its value is ignored.
     523             :     !>                                              <li>    on input, the subset `logFuncState(2:ndim+1, 1:njob)` represents the set of `njob` states (points)
     524             :     !>                                                      from within the domain of the user-specified target density function whose function values must be returned.<br>
     525             :     !>                                                      On output, the user-specified procedure `getLogFunc()` must return `njob` function values in elements `logFuncState(1, 1:njob)`
     526             :     !>                                                      corresponding to the input states `logFuncState(2:ndim+1, 1:njob)`.<br>
     527             :     !>                                              <li>    the two extra `intent(inout)` arguments `avgTimePerFunCall, avgCommPerFunCall` represent:<br>
     528             :     !>                                                      <ol>
     529             :     !>                                                          <li>    the average time (in seconds) it takes to compute the user-specified function for a **single** input state.
     530             :     !>                                                          <li>    the average time (in seconds) it takes to create images/processes/threads for parallel computation of
     531             :     !>                                                                  the user-specified function and communicate information among them, collectively called **parallelization overhead**.
     532             :     !>                                                      </ol>
     533             :     !>                                                      These two extra arguments are essential for accurate post-processing the performance of the parallel computations in the simulation.<br>
     534             :     !>                                                      If there is no parallelization within he user-specified `getLogFunc`, then the parallelization overhead can be set to zero.<br>
     535             :     !>                                          </ol>
     536             :     !>                                          This procedure interface is entirely different from the density function interfaces in the previous versions of the ParaMonte library.<br>
     537             :     !>                                          The new interface has been primarily developed to facilitate the **concurrent or parallel** evaluation
     538             :     !>                                          of target density function in higher-level languages for a set of `njob` input states.<br>
     539             :     !>                                          This interface is particularly vital in higher-level programming languages such as MATLAB, Python, and R where ParaMonte thread-level
     540             :     !>                                          parallelism is impossible due to the [global interpreter lock (GIL)](https://en.wikipedia.org/wiki/Global_interpreter_lock).<br>
     541             :     !>                                          This interface is currently activated **if and only if** the preprocessors `OMP_ENABLED=1` and either `MATLAB_ENABLED=1` or `PYTHON_ENABLED=1` or `R_ENABLED=1` are set.<br>
     542             :     !>                                  <li>    When the library is built for any build configuration other than the above with preprocessor `CDF_ENABLED=1`, particularly, for the C and C++ programming languages,
     543             :     !>                                          the input `getLogFunc()` must have the following interface:<br>
     544             :     !>                                          \code{.F90}
     545             :     !>                                              function getLogFunc(state, ndim) result(logFunc)
     546             :     !>                                                  integer(IK), value :: ndim
     547             :     !>                                                  real(RKC), intent(in) :: state(ndim)
     548             :     !>                                                  real(RKC) :: logFunc
     549             :     !>                                              end function
     550             :     !>                                          \endcode
     551             :     !>                                          where,
     552             :     !>                                          <ol>
     553             :     !>                                              <li>    the variable `ndim` is the number of dimensions of the target density function.
     554             :     !>                                              <li>    the variable `state` is an input point from within the domain of the density function.
     555             :     !>                                              <li>    the constant `RKC` refers to the `real` kind used for computations.
     556             :     !>                                              <li>    the output `logFunc` is the computed value of the density function at the specified input `state`.
     557             :     !>                                          </ol>
     558             :     !>                              </ol>
     559             :     !>  \param[in]  ndim        :   The input scalar `integer` of default kind \IK representing the number of
     560             :     !>                              dimensions of the domain of the objective function that is to be sampled.<br>
     561             :     !>
     562             :     !>  \return
     563             :     !>  `err`                   :   The output scalar of type [err_type](@ref pm_err::err_type) whose component `occur` is set to the `.true.` **if and only if** the
     564             :     !>                              procedure fails to accomplish the sampling task, in which case, the `msg` component of `err` is set to a description of the error origin.<br>
     565             :     !>                              In such a case, the description of the possible cause(s) of the error will be written to the output *report* file that is automatically generated for all simulations.<br>
     566             :     !>                              The primary reason for simulation failures is a wrong syntactic or semantic input value(s) for the simulation settings specified within the components of the input
     567             :     !>                              argument `sampler` or within the external output file specified by the `inputFile` component of the input `sampler`.<br>
     568             :     !>                              On output, the `msg` component of `err` is an empty string if the sampling completes successfully.<br>
     569             :     !>
     570             :     !>  \interface{getErrSampling}
     571             :     !>  \code{.F90}
     572             :     !>
     573             :     !>      use pm_sampling, only: getErrSampling
     574             :     !>
     575             :     !>      err = getErrSampling(sampler, getLogFunc, ndim)
     576             :     !>
     577             :     !>  \endcode
     578             :     !>
     579             :     !>  \warning
     580             :     !>  The condition `0 < ndim` must hold for the corresponding input arguments.<br>
     581             :     !>  \vericons
     582             :     !>
     583             :     !>  \impure
     584             :     !>
     585             :     !>  \example{getErrSampling}
     586             :     !>  \include{lineno} example/pm_sampling/getErrSampling/main.F90
     587             :     !>  \compilef{getErrSampling}
     588             :     !>  \output{getErrSampling}
     589             :     !>  \include{lineno} example/pm_sampling/getErrSampling/main.out.F90
     590             :     !>
     591             :     !>  \example{mvn}
     592             :     !>  \include{lineno} example/pm_sampling/mvn/main.F90
     593             :     !>  \inputfile{mvn, input.nml, example specifications namelist input file}
     594             :     !>  \include{lineno} example/pm_sampling/mvn/input.nml
     595             :     !>  \compilef{mvn}
     596             :     !>  \postproc{mvn}
     597             :     !>  \include{lineno} example/pm_sampling/mvn/main.py
     598             :     !>  \vis{mvn}
     599             :     !>  \image html example/pm_sampling/mvn/mvn.traceplot.png width=700
     600             :     !>  \image html example/pm_sampling/mvn/mvn.scatterplot.png width=700
     601             :     !>  \image html example/pm_sampling/mvn/mvn.adaptationMeasure.png width=700
     602             :     !>
     603             :     !>  \example{himmelblau}
     604             :     !>  \include{lineno} example/pm_sampling/himmelblau/main.F90
     605             :     !>  \inputfile{himmelblau, input.nml, example specifications namelist input file}
     606             :     !>  \include{lineno} example/pm_sampling/himmelblau/input.nml
     607             :     !>  \compilef{himmelblau}
     608             :     !>  \postproc{himmelblau}
     609             :     !>  \include{lineno} example/pm_sampling/himmelblau/main.py
     610             :     !>  \vis{himmelblau}
     611             :     !>  \image html example/pm_sampling/himmelblau/himmelblau.traceplot.png width=700
     612             :     !>  \image html example/pm_sampling/himmelblau/himmelblau.scatterplot.png width=700
     613             :     !>  \image html example/pm_sampling/himmelblau/himmelblau.adaptationMeasure.png width=700
     614             :     !>
     615             :     !>  \example{mixture}
     616             :     !>  \include{lineno} example/pm_sampling/mixture/main.F90
     617             :     !>  \compilef{mixture}
     618             :     !>  \postproc{mixture}
     619             :     !>  \include{lineno} example/pm_sampling/mixture/main.py
     620             :     !>  \vis{mixture}
     621             :     !>  \image html example/pm_sampling/mixture/mixLogNormLogNorm.hist.png width=700
     622             :     !>  \image html example/pm_sampling/mixture/mixLogNormLogNorm.traceplot.png width=700
     623             :     !>  \image html example/pm_sampling/mixture/mixLogNormLogNorm.adaptationMeasure.png width=700
     624             :     !>  <br>
     625             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoFlatPowetoTapered.hist.png width=700
     626             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoFlatPowetoTapered.traceplot.png width=700
     627             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoFlatPowetoTapered.adaptationMeasure.png width=700
     628             :     !>  <br>
     629             :     !>  \image html example/pm_sampling/mixture/mixLogNormFlatPowetoTapered.hist.png width=700
     630             :     !>  \image html example/pm_sampling/mixture/mixLogNormFlatPowetoTapered.traceplot.png width=700
     631             :     !>  \image html example/pm_sampling/mixture/mixLogNormFlatPowetoTapered.adaptationMeasure.png width=700
     632             :     !>  <br>
     633             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoLogNorm.hist.png width=700
     634             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoLogNorm.traceplot.png width=700
     635             :     !>  \image html example/pm_sampling/mixture/mixFlatPowetoLogNorm.adaptationMeasure.png width=700
     636             :     !>
     637             :     !>  \test
     638             :     !>  [test_pm_sampling](@ref test_pm_sampling)
     639             :     !>
     640             :     !>  \finmain{getErrSampling}
     641             :     !>
     642             :     !>  \author
     643             :     !>  \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
     644             :     interface getErrSampling
     645             : 
     646             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     647             : 
     648             : #if RK5_ENABLED
     649             :     impure module function getErrParaDRAM_RK5(sampler, getLogFunc, ndim) result(err)
     650             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     651             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK5
     652             : #endif
     653             :         use pm_kind, only: RKC => RK5
     654             :         procedure(getLogFunc_proc_RK5) :: getLogFunc
     655             :         type(paradram_type), intent(in) :: sampler
     656             :         integer(IK), intent(in) :: ndim
     657             :         type(err_type) :: err
     658             :     end function
     659             : #endif
     660             : 
     661             : #if RK4_ENABLED
     662             :     impure module function getErrParaDRAM_RK4(sampler, getLogFunc, ndim) result(err)
     663             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     664             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK4
     665             : #endif
     666             :         use pm_kind, only: RKC => RK4
     667             :         procedure(getLogFunc_proc_RK4) :: getLogFunc
     668             :         type(paradram_type), intent(in) :: sampler
     669             :         integer(IK), intent(in) :: ndim
     670             :         type(err_type) :: err
     671             :     end function
     672             : #endif
     673             : 
     674             : #if RK3_ENABLED
     675             :     impure module function getErrParaDRAM_RK3(sampler, getLogFunc, ndim) result(err)
     676             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     677             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK3
     678             : #endif
     679             :         use pm_kind, only: RKC => RK3
     680             :         procedure(getLogFunc_proc_RK3) :: getLogFunc
     681             :         type(paradram_type), intent(in) :: sampler
     682             :         integer(IK), intent(in) :: ndim
     683             :         type(err_type) :: err
     684             :     end function
     685             : #endif
     686             : 
     687             : #if RK2_ENABLED
     688             :     impure module function getErrParaDRAM_RK2(sampler, getLogFunc, ndim) result(err)
     689             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     690             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK2
     691             : #endif
     692             :         use pm_kind, only: RKC => RK2
     693             :         procedure(getLogFunc_proc_RK2) :: getLogFunc
     694             :         type(paradram_type), intent(in) :: sampler
     695             :         integer(IK), intent(in) :: ndim
     696             :         type(err_type) :: err
     697             :     end function
     698             : #endif
     699             : 
     700             : #if RK1_ENABLED
     701             :     impure module function getErrParaDRAM_RK1(sampler, getLogFunc, ndim) result(err)
     702             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     703             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK1
     704             : #endif
     705             :         use pm_kind, only: RKC => RK1
     706             :         procedure(getLogFunc_proc_RK1) :: getLogFunc
     707             :         type(paradram_type), intent(in) :: sampler
     708             :         integer(IK), intent(in) :: ndim
     709             :         type(err_type) :: err
     710             :     end function
     711             : #endif
     712             : 
     713             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     714             : 
     715             : #if RK5_ENABLED
     716             :     impure module function getErrParaDISE_RK5(sampler, getLogFunc, ndim) result(err)
     717             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     718             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK5
     719             : #endif
     720             :         use pm_kind, only: RKC => RK5
     721             :         procedure(getLogFunc_proc_RK5) :: getLogFunc
     722             :         type(paradise_type), intent(in) :: sampler
     723             :         integer(IK), intent(in) :: ndim
     724             :         type(err_type) :: err
     725             :     end function
     726             : #endif
     727             : 
     728             : #if RK4_ENABLED
     729             :     impure module function getErrParaDISE_RK4(sampler, getLogFunc, ndim) result(err)
     730             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     731             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK4
     732             : #endif
     733             :         use pm_kind, only: RKC => RK4
     734             :         procedure(getLogFunc_proc_RK4) :: getLogFunc
     735             :         type(paradise_type), intent(in) :: sampler
     736             :         integer(IK), intent(in) :: ndim
     737             :         type(err_type) :: err
     738             :     end function
     739             : #endif
     740             : 
     741             : #if RK3_ENABLED
     742             :     impure module function getErrParaDISE_RK3(sampler, getLogFunc, ndim) result(err)
     743             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     744             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK3
     745             : #endif
     746             :         use pm_kind, only: RKC => RK3
     747             :         procedure(getLogFunc_proc_RK3) :: getLogFunc
     748             :         type(paradise_type), intent(in) :: sampler
     749             :         integer(IK), intent(in) :: ndim
     750             :         type(err_type) :: err
     751             :     end function
     752             : #endif
     753             : 
     754             : #if RK2_ENABLED
     755             :     impure module function getErrParaDISE_RK2(sampler, getLogFunc, ndim) result(err)
     756             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     757             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK2
     758             : #endif
     759             :         use pm_kind, only: RKC => RK2
     760             :         procedure(getLogFunc_proc_RK2) :: getLogFunc
     761             :         type(paradise_type), intent(in) :: sampler
     762             :         integer(IK), intent(in) :: ndim
     763             :         type(err_type) :: err
     764             :     end function
     765             : #endif
     766             : 
     767             : #if RK1_ENABLED
     768             :     impure module function getErrParaDISE_RK1(sampler, getLogFunc, ndim) result(err)
     769             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     770             :         !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK1
     771             : #endif
     772             :         use pm_kind, only: RKC => RK1
     773             :         procedure(getLogFunc_proc_RK1) :: getLogFunc
     774             :         type(paradise_type), intent(in) :: sampler
     775             :         integer(IK), intent(in) :: ndim
     776             :         type(err_type) :: err
     777             :     end function
     778             : #endif
     779             : 
     780             :     !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     781             : 
     782             :     end interface
     783             : 
     784             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     785             : 
     786             : contains
     787             : 
     788             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     789             : 
     790             :     !>  \name runParaDRAM
     791             :     !>  \{
     792             :     !>
     793             :     !>  \cfi
     794             :     !>
     795             :     !>  \brief
     796             :     !>  Generate and return a non-zero value (`1`) if the procedure fails to fully accomplish the task of generating
     797             :     !>  a Monte Carlo sample of the specified input mathematical objective function, otherwise, return `0`.
     798             :     !>
     799             :     !>  \details
     800             :     !>  This interface group is the entry point to all **C-style interfaces** to the ParaDRAM samplers of mathematical density functions.<br>
     801             :     !>  Although the procedures of this generic interface return a single scalar of type `int32_t`, the procedures generate
     802             :     !>  massive amounts of information about each simulation which are stored in appropriate external hard drive files.<br>
     803             :     !>
     804             :     !>  \param[in]  getLogFunc  :   The input user-specified procedure pointer to the natural logarithm of the target density function.<br>
     805             :     !>                              On input, the function must take an input vector `state` of size `ndim` of floating-point type of kind \C_RKALL
     806             :     !>                              representing a state (point) from within the domain of the user-specified target density function whose function value must be returned.<br>
     807             :     !>                              On output, the user-specified procedure `getLogFunc()` must return the function value corresponding to the input `state[ndim]`.<br>
     808             :     !>                              The following illustrate the generic interface of input function pointer `getLogFunc(state, ndim)`,
     809             :     !>                              \code{.F90}
     810             :     !>                                  REAL getLogFunc(REAL[] state, int32_t ndim);
     811             :     !>                              \endcode
     812             :     !>                              where `REAL` can be a floating-point type of kind \C_RKALL for the corresponding varying-precision sampler interfaces:<br>
     813             :     !>                              <ol>
     814             :     !>                                  <li>    [runParaDRAMF](@ref pm_sampling::runParaDRAMF),
     815             :     !>                                  <li>    [runParaDRAMD](@ref pm_sampling::runParaDRAMD),
     816             :     !>                                  <li>    [runParaDRAML](@ref pm_sampling::runParaDRAML).
     817             :     !>                              </ol>
     818             :     !>  \param[in]  ndim        :   The input scalar constant of type `int32_t` representing the number of dimensions of the domain of the objective function.
     819             :     !>  \param[in]  input       :   The input scalar pointer of type `char` representing the null-terminated C string
     820             :     !>                              that contains either,
     821             :     !>                              <ol>
     822             :     !>                                  <li>    the path to an external input file containing the namelist group of ParaDRAM sampler specifications
     823             :     !>                                          as outlined in the corresponding page of [ParaMonte library generic documentation website](\pmdoc_usage_sampling/paradram/specifications/).<br>
     824             :     !>                                  <li>    the namelist group of ParaDRAM sampler specifications as the can appear in an external input specification file.<br>
     825             :     !>                              </ol>
     826             :     !>                              While all input simulation specifications are optional, it is highly recommended to pay
     827             :     !>                              attention to the default settings of the domain boundaries and sampler starting point.<br>
     828             :     !>                              (**optional**. It is considered as missing if set to `null()`.)
     829             :     !>
     830             :     !>  \return
     831             :     !>  `stat`                  :   The output scalar of type `int32_t` that is `0` **if and only if**
     832             :     !>                              the sampler succeeds in sampling the specified density function.<br>
     833             :     !>
     834             :     !>  \interface{runParaDRAM}
     835             :     !>  \code{.F90}
     836             :     !>
     837             :     !>      stat = runParaDRAMF(getLogFunc, ndim, input)
     838             :     !>      stat = runParaDRAMD(getLogFunc, ndim, input)
     839             :     !>      stat = runParaDRAML(getLogFunc, ndim, input)
     840             :     !>
     841             :     !>  \endcode
     842             :     !>
     843             :     !>  \warning
     844             :     !>  Beware that the definition of extended precision real type `long double` is compiler and platform dependent
     845             :     !>  making the use of `long double` with precompiled ParaMonte libraries problematic and non-functional.<br>
     846             :     !>
     847             :     !>  \warning
     848             :     !>  The condition `0 < ndim` must hold for the corresponding input arguments.<br>
     849             :     !>  \vericon
     850             :     !>
     851             :     !>  \example{runParaDRAM}
     852             :     !>  \include{lineno} example/pm_sampling/runParaDRAM/main.F90
     853             :     !>  \compilef{runParaDRAM}
     854             :     !>  \postproc{runParaDRAM}
     855             :     !>  \include{lineno} example/pm_sampling/runParaDRAM/main.py
     856             :     !>  \vis{runParaDRAM}
     857             :     !>  \image html example/pm_sampling/runParaDRAM/runParaDRAM.traceplot.png width=700
     858             :     !>  \image html example/pm_sampling/runParaDRAM/runParaDRAM.scatterplot.png width=700
     859             :     !>  \image html example/pm_sampling/runParaDRAM/runParaDRAM.adaptationMeasure.png width=700
     860             :     !>
     861             :     !>  \test
     862             :     !>  [test_pm_sampling](@ref test_pm_sampling)
     863             :     !>
     864             :     !>  \bug
     865             :     !>  \status \unresolved
     866             :     !>  \source \ifx{2024.0.0 20231017}
     867             :     !>  \desc
     868             :     !>  \ifx This interface yields a segmentation fault error for all `real` types supported when linked with the ParaMonte library built with \ifx.
     869             :     !>  \remedy
     870             :     !>  For now, only \ifort will be used.<br>
     871             :     !>
     872             :     !>  \bug
     873             :     !>  \status \unresolved
     874             :     !>  \source \ifort{2021.11.0 20231010}
     875             :     !>  \desc
     876             :     !>  The `runParaDRAML` interface for `long double` yields a segmentation fault error.
     877             :     !>  \remedy
     878             :     !>  None as of today.<br>
     879             :     !>
     880             :     !>  \todo
     881             :     !>  \pvhigh
     882             :     !>  The current tests for `long double` interface fail with \ifx and \icx compilers.<br>
     883             :     !>  Apparently, there is a mismatch between `long double` and `c_long_double` storage.<br>
     884             :     !>  This issue does not exist with GNU compilers, although the GNU definition of `long double`
     885             :     !>  appears to yield incorrect values in some calculations, e.g., in `isFailedGeomCyclicFit()` of the ParaMonte library.<br>
     886             :     !>
     887             :     !>  \finminc{runParaDRAM}
     888             :     !>
     889             :     !>  \author
     890             :     !>  \AmirShahmoradi, May 16 2016, 9:03 AM, Oden Institute for Computational Engineering and Sciences (ICES), UT Austin
     891           0 :     function runParaDRAML(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAML")
     892             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     893             :         !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAML
     894             : #endif
     895             :         use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_long_double
     896             :         integer, parameter :: RKC = merge(c_long_double, RKH, any(c_long_double == RKALL))
     897             :         character(1, SK), intent(in), optional :: input(*)
     898             :         type(c_funptr), intent(in), value :: getLogFunc
     899             :         integer(IK), intent(in), value :: ndim
     900             :         integer(IK) :: stat
     901             : #define runParaDRAM_ENABLED 1
     902             : #include "pm_sampling@routines.inc.F90"
     903             : #undef  runParaDRAM_ENABLED
     904             :     end function
     905             : 
     906           1 :     function runParaDRAMD(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMD")
     907             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     908             :         !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMD
     909             : #endif
     910             :         use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_double
     911             :         integer, parameter :: RKC = merge(c_double, RKD, any(c_double == RKALL))
     912             :         character(1, SK), intent(in), optional :: input(*)
     913             :         type(c_funptr), intent(in), value :: getLogFunc
     914             :         integer(IK), intent(in), value :: ndim
     915             :         integer(IK) :: stat
     916             : #define runParaDRAM_ENABLED 1
     917             : #include "pm_sampling@routines.inc.F90"
     918             : #undef  runParaDRAM_ENABLED
     919             :     end function
     920             : 
     921           0 :     function runParaDRAMF(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMF")
     922             : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
     923             :         !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMF
     924             : #endif
     925             :         use pm_kind, only: RKALL
     926             :         use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_float
     927             :         integer, parameter :: RKC = merge(c_float, RKL, any(c_float == RKALL))
     928             :         character(1, SK), intent(in), optional :: input(*)
     929             :         type(c_funptr), intent(in), value :: getLogFunc
     930             :         integer(IK), intent(in), value :: ndim
     931             :         integer(IK) :: stat
     932             : #define runParaDRAM_ENABLED 1
     933             : #include "pm_sampling@routines.inc.F90"
     934             : #undef  runParaDRAM_ENABLED
     935             :     end function
     936             :     !>  \}
     937             : 
     938             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     939             : 
     940             : end module pm_sampling ! LCOV_EXCL_LINE

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