ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_sampling Module Reference

This module contains procedures and generic interfaces for the ParaMonte library sampler routines. More...

Data Types

interface  getErrSampling
 Generate and return .true. if the procedure fails to fully accomplish the task of Monte Carlo sampling of the specified input mathematical objective function, otherwise, return .false..
More...
 
type  paradise_type
 This is a derived type for constructing objects containing the optional simulation properties of the ParaMonte library DRAM-enhanced MCMC explorers and samplers.
More...
 
type  paradram_type
 This is a derived type for constructing objects containing the optional simulation properties of the ParaMonte library Delayed-Rejection Adaptive Metropolis (DRAM) MCMC explorers and samplers.
More...
 
type  paramcmc_type
 This is a derived type for constructing objects containing the optional simulation properties of the ParaMonte library MCMC explorers and samplers.
More...
 
type  paranest_type
 This is a derived type for constructing objects containing the optional simulation properties of the ParaMonte library Nested explorers and samplers.
More...
 
type  sampler_type
 This is a derived type for constructing objects containing the optional simulation properties of the ParaMonte library explorers and samplers.
More...
 

Functions/Subroutines

runParaDRAM

C-Fortran Interoperation

Generate and return a non-zero value (1) if the procedure fails to fully accomplish the task of Monte Carlo sampling of the specified input mathematical objective function, otherwise, return 0.

See the documentation of getErrSampling for example usage in Fortran programming language.
This interface group is the entry point to all C-style interfaces to the ParaDRAM samplers of mathematical density functions.
Although the procedures of this generic interface return a single scalar of type int32_t, the procedures generate massive amounts of information about each simulation which are stored in appropriate external hard drive files.

See,

  1. this generic documentation page for more information on the generated output files for samplings performed using the ParaDRAM sampler.
Parameters
[in]getLogFunc: The input user-specified procedure pointer to the natural logarithm of the target density function.
On input, the function must take an input vector state of size ndim of floating-point type of kind float, double, or long double representing a state (point) from within the domain of the user-specified target density function whose function value must be returned.
On output, the user-specified procedure getLogFunc() must return the function value corresponding to the input state[ndim].
The following illustrate the generic interface of input function pointer getLogFunc(state, ndim),
REAL getLogFunc(REAL[] state, int32_t ndim);
where REAL can be a floating-point type of kind float, double, or long double for the corresponding varying-precision sampler interfaces:
  1. runParaDRAMF,
  2. runParaDRAMD,
  3. runParaDRAML.
[in]ndim: The input scalar constant of type int32_t representing the number of dimensions of the domain of the objective function.
[in]input: The input scalar pointer of type char representing the null-terminated C string that contains either,
  1. the path to an external input file containing the namelist group of ParaDRAM sampler specifications as outlined in the corresponding page of ParaMonte library generic documentation website.
  2. the namelist group of ParaDRAM sampler specifications as the can appear in an external input specification file.
While all input simulation specifications are optional, it is highly recommended to pay attention to the default settings of the domain boundaries and sampler starting point.
(optional. It is considered as missing if set to null().)
Returns
stat : The output scalar of type int32_t that is 0 if and only if the sampler succeeds in sampling the specified density function.


Possible calling interfaces

stat = runParaDRAMF(getLogFunc, ndim, input)
stat = runParaDRAMD(getLogFunc, ndim, input)
stat = runParaDRAML(getLogFunc, ndim, input)
Warning
Beware that the definition of extended precision real type long double is compiler and platform dependent.
This makes the use of long double with precompiled ParaMonte libraries problematic and non-functional.
The condition 0 < ndim must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.


Example usage

1module logfunc
2 use iso_c_binding, only: SK => c_char, IK => c_int32_t, RKG => c_double ! All real kinds are supported.
3 implicit none
4 integer(IK) , parameter :: NDIM = 4_IK
5 real(RKG) , parameter :: MEAN(NDIM) = [-6, -2, 2, 6] ! The mean vector of the MVN
6 real(RKG) , parameter :: LOG_INVERSE_SQRT_TWO_PI = log(1./sqrt(2.*acos(-1.))) ! MVN distribution coefficient: log(1/sqrt(2*Pi)^ndim).
7 real(RKG) , parameter :: COVMAT(NDIM, NDIM) = reshape([ +1.0, +0.5, +0.5, +0.5 &
8 , +0.5, +1.0, +0.5, +0.5 &
9 , +0.5, +0.5, +1.0, +0.5 &
10 , +0.5, +0.5, +0.5, +1.0 ], shape(COVMAT)) ! The covariance matrix of the MVN
11 real(RKG) , parameter :: INVCOV(NDIM, NDIM) = reshape([ +1.6, -0.4, -0.4, -0.4 &
12 , -0.4, +1.6, -0.4, -0.4 &
13 , -0.4, -0.4, +1.6, -0.4 &
14 , -0.4, -0.4, -0.4, +1.6 ], shape(INVCOV)) ! The inverse covariance matrix of the MVN
15 real(RKG) , parameter :: MVN_COEF = NDIM * LOG_INVERSE_SQRT_TWO_PI + 0.581575404902840 ! log(sqrt(det(INVCOV)))
16contains
17 function getLogFunc(state, ndim) result(logFunc) bind(C)
18 ! Return the negative natural logarithm of MVN distribution
19 ! evaluated at the input vector `state` of length `ndim`.
20 integer(IK), value :: ndim
21 real(RKG), intent(in) :: state(ndim)
22 real(RKG) :: stateNormed(size(state))
23 real(RKG) :: logFunc
24 stateNormed = state - MEAN
25 logFunc = MVN_COEF - 0.5_RKG * (dot_product(stateNormed, matmul(INVCOV, stateNormed)))
26 end function
27end module logfunc
28
29program example
30 use pm_sampling, only: runParaDRAMD
31 use iso_c_binding, only: c_funloc, c_null_char
32 use logfunc, only: SK, IK, RKG, NDIM, getLogFunc
33 character(len(c_null_char), SK), parameter :: NUL = c_null_char
34 integer(IK) :: stat
35 stat = runParaDRAMD(c_funloc(getLogFunc), NDIM, SK_"&paradram outputFileName = './out/runParaDRAM', outputStatus = 'retry', outputChainSize = 30000 /"//NUL)
36 if (stat /= 0) error stop "sampler failed."
37end
This module contains procedures and generic interfaces for the ParaMonte library sampler routines.
Definition: pm_sampling.F90:33
integer(IK) function runParaDRAMD(getLogFunc, ndim, input)

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Postprocessing of the example output
1#!/usr/bin/env python
2
3import matplotlib.pyplot as plt
4#import paramonte as pm
5import pandas as pd
6import numpy as np
7import glob
8import sys
9import os
10funcname = os.path.basename(os.getcwd())
11
12linewidth = 2
13fontsize = 17
14files = glob.glob("./out/*_chain.txt")
15
16for file in files:
17
18 df = pd.read_csv(file, delimiter = ",")
19 sindex = list(df.columns).index("sampleLogFunc") + 1 # start column index.
20
21 # traceplot
22
23 #print(df.values)
24 fig = plt.figure(figsize = (8, 6))
25 ax = plt.subplot(1,1,1)
26 ax.plot ( range(len(df.values[:, 0]))
27 , df.values[:, sindex:]
28 , zorder = 1000
29 )
30 plt.minorticks_on()
31 ax.set_xscale("log")
32 ax.set_xlabel("MCMC Count", fontsize = 17)
33 ax.set_ylabel("MCMC State", fontsize = 17)
34 ax.tick_params(axis = "x", which = "minor")
35 ax.tick_params(axis = "y", which = "minor")
36 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
37 #ax.legend(df.columns.values, fontsize = fontsize)
38 plt.tight_layout()
39 plt.savefig(funcname + ".traceplot.png")
40
41 # scatterplot
42
43 if len(df.values[1, sindex:]) > 1:
44 #print(df.values)
45 fig = plt.figure(figsize = (8, 6))
46 ax = plt.subplot(1,1,1)
47 ax.scatter ( df.values[:, sindex]
48 , df.values[:, sindex + 1]
49 , zorder = 1000
50 , alpha = 0.5
51 , s = 1
52 )
53 plt.minorticks_on()
54 #ax.set_xscale("log")
55 ax.set_xlabel(df.columns.values[sindex], fontsize = 17)
56 ax.set_ylabel(df.columns.values[sindex + 1], fontsize = 17)
57 ax.tick_params(axis = "x", which = "minor")
58 ax.tick_params(axis = "y", which = "minor")
59 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
60 #ax.legend(df.columns.values, fontsize = fontsize)
61 plt.tight_layout()
62 plt.savefig(funcname + ".scatterplot.png")
63
64 # adaptation measure.
65
66 #print(df.values)
67 if "proposalAdaptation" in df.columns.values:
68 if any(df["proposalAdaptation"].values != 0):
69 fig = plt.figure(figsize = (8, 6))
70 ax = plt.subplot(1,1,1)
71 ax.scatter ( range(len(df.values[:, 0]))
72 , df["proposalAdaptation"].values
73 , zorder = 1000
74 , alpha = 0.5
75 , c = "red"
76 , s = 1
77 )
78 plt.minorticks_on()
79 #ax.set_xscale("log")
80 ax.set_yscale("log")
81 ax.set_xlabel("MCMC Count", fontsize = 17)
82 ax.set_ylabel("proposalAdaptation", fontsize = 17)
83 ax.tick_params(axis = "x", which = "minor")
84 ax.tick_params(axis = "y", which = "minor")
85 plt.grid(visible = True, which = "both", axis = "both", color = "0.85", linestyle = "-")
86 #ax.legend(df.columns.values, fontsize = fontsize)
87 plt.tight_layout()
88 plt.savefig(funcname + ".proposalAdaptation.png")
89
90 #sim = pm.ParaDRAM()
91 #sample = sim.readSample("./out/", renabled = True)
92 #sample[0].plot.grid()

Visualization of the example output
Test:
test_pm_sampling
Bug:

Status: Unresolved
Source: Intel LLVM Fortran Compiler ifx version 2024.0.0 20231017
Description: Intel LLVM Fortran Compiler ifx This interface yields a segmentation fault error for all real types supported when linked with the ParaMonte library built with Intel LLVM Fortran Compiler ifx.

Remedy (as of ParaMonte Library version 2.0.0): For now, only Intel Classic Fortran Compiler ifort will be used.
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.11.1 20231117
Description: Intel Classic Fortran Compiler ifort returns an *already allocated error with the statement call setResized(scaling, lenScaling) which persists in both release and debug modes.
The full debug message is the following:
forrtl: severe (151): allocatable array is already allocated
Image PC Routine Line Source
libparamonte.so 00007FCFBA641EB8 Unknown Unknown Unknown
libparamonte.so 00007FCFB6E6D4A3 pm_arrayresize_MP 170 pm_arrayResize@routines.inc.F90
libparamonte.so 00007FCFB72B5BD6 pm_parallelism_MP 77 pm_parallelism@routines.inc.F90
libparamonte.so 00007FCFB6263F43 pm_sampling_MP_ge 394 pm_sampling@routines.inc.F90
libparamonte.so 00007FCFB61F0194 runParaDRAMD 136 pm_sampling@routines.inc.F90

Note that the line numbers for this file in the message above have changed because of code change in this file.
This error does not occur when the library is compiled with GNU Fortran Compiler gfortran version 13.
This error does not occur when the library is compiled with Intel LLVM Fortran Compiler ifx version 2025.0.0 20241008.
It seems like this error occurs because of placing the following typed variable speedup in the forkjoin_parallelism_block below.

Remedy (as of ParaMonte Library version 2.0.0): For now, the type definition and the typed variable declaration are taken out of the block and placed below.
This must be checked with newer Intel compilers as Intel Classic Fortran Compiler ifort version 2021.11.1 20231117 is being phased out by Intel.

Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.11.1 20231117
Description: The following declarations belong to only parallel multichain modes within a nested block.
type :: probKS_type
real(RKG) , allocatable :: values(:)
real(RKG) :: minval
integeR(IK) :: minloc
integeR(IK) :: minpid
end type
type :: sampleLogFuncState_type
character(:, SK), allocatable :: filePath
real(RKG) , allocatable :: thisImage(:,:)
real(RKG) , allocatable :: thatImage(:,:)
type(probKS_type) :: probKS
end type

But they had to be taken out of their local scope because Intel Classic Fortran Compiler ifort version 2021.11.1 20231117 cannot compile them with an ICE message as below.

pm_sampling@routines.inc.F90(1267): catastrophic error:
**Internal compiler error: internal abort** Please report this error along with the circumstances in which it occurred in a Software Problem Report.
Note: File and line given may not be explicit cause of this error.
compilation aborted for /home/amir/git/paramonte/src/fortran/main/pm_sampling@routines.F90 (code 1)


Remedy (as of ParaMonte Library version 2.0.0): For now, the type definition and the typed variable declaration are taken out of the block and placed below.
This must be checked with newer Intel compilers as Intel Classic Fortran Compiler ifort version 2021.11.1 20231117 is being phased out by Intel.
This bug may have the same origins as the bug in the above.

Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.11.0 20231010
Description: The runParaDRAML interface for long double yields a segmentation fault error.

Remedy (as of ParaMonte Library version 2.0.0): None as of today.
Todo:
Critical Priority: The current tests for long double interface fail with Intel LLVM Fortran Compiler ifx and Intel LLVM C/C++ Compiler icx compilers.
Apparently, there is a mismatch between long double and c_long_double storage.
This issue does not exist with GNU compilers, although the GNU definition of long double appears to yield incorrect values in some calculations, e.g., in isFailedGeomCyclicFit() of the ParaMonte library.


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, May 16 2016, 9:03 AM, Oden Institute for Computational Engineering and Sciences (ICES), UT Austin
integer(IK) function runParaDRAML (getLogFunc, ndim, input)
 
integer(IK) function runParaDRAMD (getLogFunc, ndim, input)
 
integer(IK) function runParaDRAMF (getLogFunc, ndim, input)
 

Variables

character(*, SK), parameter MODULE_NAME = "@pm_sampling"
 

Detailed Description

This module contains procedures and generic interfaces for the ParaMonte library sampler routines.

See the documentation of getErrSampling for example usage.

Test:
test_pm_sampling


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Monday 00:01 AM, January 1, 2018, Institute for Computational Engineering and Sciences, University of Texas Austin

Function/Subroutine Documentation

◆ runParaDRAMD()

integer(IK) function pm_sampling::runParaDRAMD ( type(c_funptr)  getLogFunc,
integer(IK), intent(in), value  ndim,
character(1, SK), dimension(*), intent(in), optional  input 
)

Definition at line 1006 of file pm_sampling.F90.

◆ runParaDRAMF()

integer(IK) function pm_sampling::runParaDRAMF ( type(c_funptr)  getLogFunc,
integer(IK), intent(in), value  ndim,
character(1, SK), dimension(*), intent(in), optional  input 
)

Definition at line 1021 of file pm_sampling.F90.

References pm_kind::IK, pm_kind::RKALL, and pm_kind::SK.

◆ runParaDRAML()

integer(IK) function pm_sampling::runParaDRAML ( type(c_funptr)  getLogFunc,
integer(IK), intent(in), value  ndim,
character(1, SK), dimension(*), intent(in), optional  input 
)

Definition at line 991 of file pm_sampling.F90.

Variable Documentation

◆ MODULE_NAME

character(*, SK), parameter pm_sampling::MODULE_NAME = "@pm_sampling"

Definition at line 43 of file pm_sampling.F90.