ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_sampling.F90
Go to the documentation of this file.
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
27
28!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
29
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, SKG => 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
77 character(:,SKG) , allocatable :: description
78 character(:,SKG) , allocatable :: domain
79 character(:,SKG) , allocatable :: domainAxisName(:)
80 real(RKH) , allocatable :: domainBallAvg(:)
81 real(RKH) , allocatable :: domainBallCor(:,:)
82 real(RKH) , allocatable :: domainBallCov(:,:)
83 real(RKH) , allocatable :: domainBallStd(:)
84 real(RKH) , allocatable :: domainCubeLimitLower(:)
85 real(RKH) , allocatable :: domainCubeLimitUpper(:)
86 integer(IK) , allocatable :: domainErrCount
87 integer(IK) , allocatable :: domainErrCountMax
88 character(:,SKG) , allocatable :: inputFile / / /
151 !logical(LK) , allocatable :: inputFileHasPriority !<
152 character(:,SKG) , allocatable :: outputChainFileFormat
153 integer(IK) , allocatable :: outputColumnWidth
154 character(:,SKG) , allocatable :: outputFileName
155 character(:,SKG) , allocatable :: outputStatus
156 integer(IK) , allocatable :: outputPrecision
157 integer(IK) , allocatable :: outputReportPeriod
158 character(:,SKG) , allocatable :: outputRestartFileFormat
159 integer(IK) , allocatable :: outputSampleSize
160 character(:,SKG) , allocatable :: outputSeparator
161 character(:,SKG) , allocatable :: outputSplashMode
162 character(:,SKG) , allocatable :: parallelism
163 logical(LK) , allocatable :: parallelismMpiFinalizeEnabled
164 integer(IK) , allocatable :: parallelismNumThread
165 integer(IK) , allocatable :: randomSeed
166 character(:,SKG) , allocatable :: sysInfoFilePath
167 real(RKH) , allocatable :: targetAcceptanceRate(:)
168 end type
169
204 type, extends(sampler_type) :: paramcmc_type
205 integer(IK) , allocatable :: outputChainSize
206 integer(IK) , allocatable :: outputSampleRefinementCount
207 character(:,SKG) , allocatable :: outputSampleRefinementMethod
208 character(:,SKG) , allocatable :: proposal
209 real(RKH) , allocatable :: proposalCor(:,:)
210 real(RKH) , allocatable :: proposalCov(:,:)
211 character(:,SKG) , allocatable :: proposalScale
212 real(RKH) , allocatable :: proposalStart(:)
213 real(RKH) , allocatable :: proposalStartDomainCubeLimitLower(:)
214 real(RKH) , allocatable :: proposalStartDomainCubeLimitUpper(:)
215 logical(LK) , allocatable :: proposalStartRandomized
216 real(RKH) , allocatable :: proposalStd(:)
217 end type
218
244 type, extends(paramcmc_type) :: paradram_type
245 real(RKH) , allocatable :: proposalAdaptationBurnin
246 integer(IK) , allocatable :: proposalAdaptationCount
247 integer(IK) , allocatable :: proposalAdaptationCountGreedy
248 integer(IK) , allocatable :: proposalAdaptationPeriod
249 integer(IK) , allocatable :: proposalDelayedRejectionCount
250 real(RKH) , allocatable :: proposalDelayedRejectionScale(:)
251 end type
252
278 type, extends(paradram_type) :: paradise_type
279 end type
280
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(:,SKG) , allocatable :: domainPartitionMethod
323 character(:,SKG) , allocatable :: domainPartitionObject
324 logical(LK) , allocatable :: domainPartitionOptimizationScaleEnabled
325 logical(LK) , allocatable :: domainPartitionOptimizationShapeEnabled
326 logical(LK) , allocatable :: domainPartitionOptimizationShapeScaleEnabled
327 real(RKH) , allocatable :: domainPartitionScaleFactor
328 character(:,SKG) , 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.
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, RKG => RK5
351 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
352 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
353 real(RKG) :: 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, RKG => RK4
363 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
364 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
365 real(RKG) :: 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, RKG => RK3
375 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
376 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
377 real(RKG) :: 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, RKG => RK2
387 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
388 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
389 real(RKG) :: 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, RKG => RK1
399 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
400 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
401 real(RKG) :: 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: RKG => RK5
417 real(RKG), intent(in), contiguous :: state(:)
418 real(RKG) :: 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: RKG => RK4
428 real(RKG), intent(in), contiguous :: state(:)
429 real(RKG) :: 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: RKG => RK3
439 real(RKG), intent(in), contiguous :: state(:)
440 real(RKG) :: 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: RKG => RK2
450 real(RKG), intent(in), contiguous :: state(:)
451 real(RKG) :: 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: RKG => RK1
461 real(RKG), intent(in), contiguous :: state(:)
462 real(RKG) :: logFunc
463 end function
464#endif
465
466 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
467
468#endif
469
470 end interface
472
473!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
474
653
654 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
655
656#if RK5_ENABLED
657 impure module function getErrParaDRAM_RK5(sampler, getLogFunc, ndim) result(err)
658#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
659 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK5
660#endif
661 use pm_kind, only: RKG => RK5
662 procedure(getLogFunc_proc_RK5) :: getLogFunc
663 type(paradram_type), intent(in) :: sampler
664 integer(IK), intent(in) :: ndim
665 type(err_type) :: err
666 end function
667#endif
668
669#if RK4_ENABLED
670 impure module function getErrParaDRAM_RK4(sampler, getLogFunc, ndim) result(err)
671#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
672 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK4
673#endif
674 use pm_kind, only: RKG => RK4
675 procedure(getLogFunc_proc_RK4) :: getLogFunc
676 type(paradram_type), intent(in) :: sampler
677 integer(IK), intent(in) :: ndim
678 type(err_type) :: err
679 end function
680#endif
681
682#if RK3_ENABLED
683 impure module function getErrParaDRAM_RK3(sampler, getLogFunc, ndim) result(err)
684#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
685 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK3
686#endif
687 use pm_kind, only: RKG => RK3
688 procedure(getLogFunc_proc_RK3) :: getLogFunc
689 type(paradram_type), intent(in) :: sampler
690 integer(IK), intent(in) :: ndim
691 type(err_type) :: err
692 end function
693#endif
694
695#if RK2_ENABLED
696 impure module function getErrParaDRAM_RK2(sampler, getLogFunc, ndim) result(err)
697#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
698 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK2
699#endif
700 use pm_kind, only: RKG => RK2
701 procedure(getLogFunc_proc_RK2) :: getLogFunc
702 type(paradram_type), intent(in) :: sampler
703 integer(IK), intent(in) :: ndim
704 type(err_type) :: err
705 end function
706#endif
707
708#if RK1_ENABLED
709 impure module function getErrParaDRAM_RK1(sampler, getLogFunc, ndim) result(err)
710#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
711 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK1
712#endif
713 use pm_kind, only: RKG => RK1
714 procedure(getLogFunc_proc_RK1) :: getLogFunc
715 type(paradram_type), intent(in) :: sampler
716 integer(IK), intent(in) :: ndim
717 type(err_type) :: err
718 end function
719#endif
720
721 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
722
723#if RK5_ENABLED
724 impure module function getErrParaDISE_RK5(sampler, getLogFunc, ndim) result(err)
725#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
726 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK5
727#endif
728 use pm_kind, only: RKG => RK5
729 procedure(getLogFunc_proc_RK5) :: getLogFunc
730 type(paradise_type), intent(in) :: sampler
731 integer(IK), intent(in) :: ndim
732 type(err_type) :: err
733 end function
734#endif
735
736#if RK4_ENABLED
737 impure module function getErrParaDISE_RK4(sampler, getLogFunc, ndim) result(err)
738#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
739 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK4
740#endif
741 use pm_kind, only: RKG => RK4
742 procedure(getLogFunc_proc_RK4) :: getLogFunc
743 type(paradise_type), intent(in) :: sampler
744 integer(IK), intent(in) :: ndim
745 type(err_type) :: err
746 end function
747#endif
748
749#if RK3_ENABLED
750 impure module function getErrParaDISE_RK3(sampler, getLogFunc, ndim) result(err)
751#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
752 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK3
753#endif
754 use pm_kind, only: RKG => RK3
755 procedure(getLogFunc_proc_RK3) :: getLogFunc
756 type(paradise_type), intent(in) :: sampler
757 integer(IK), intent(in) :: ndim
758 type(err_type) :: err
759 end function
760#endif
761
762#if RK2_ENABLED
763 impure module function getErrParaDISE_RK2(sampler, getLogFunc, ndim) result(err)
764#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
765 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK2
766#endif
767 use pm_kind, only: RKG => RK2
768 procedure(getLogFunc_proc_RK2) :: getLogFunc
769 type(paradise_type), intent(in) :: sampler
770 integer(IK), intent(in) :: ndim
771 type(err_type) :: err
772 end function
773#endif
774
775#if RK1_ENABLED
776 impure module function getErrParaDISE_RK1(sampler, getLogFunc, ndim) result(err)
777#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
778 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK1
779#endif
780 use pm_kind, only: RKG => RK1
781 procedure(getLogFunc_proc_RK1) :: getLogFunc
782 type(paradise_type), intent(in) :: sampler
783 integer(IK), intent(in) :: ndim
784 type(err_type) :: err
785 end function
786#endif
787
788 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
789
790 end interface
791
792!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
793
794contains
795
796!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
797
905 function runParaDRAML(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAML")
906#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
907 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAML
908#endif
909 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_long_double
910 integer, parameter :: RKG = merge(c_long_double, RKH, any(c_long_double == RKALL))
911 character(1, SK), intent(in), optional :: input(*)
912 type(c_funptr), intent(in), value :: getLogFunc
913 integer(IK), intent(in), value :: ndim
914 integer(IK) :: stat
915#define runParaDRAM_ENABLED 1
916#include "pm_sampling@routines.inc.F90"
917#undef runParaDRAM_ENABLED
918 end function
919
920 function runParaDRAMD(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMD")
921#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
922 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMD
923#endif
924 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_double
925 integer, parameter :: RKG = merge(c_double, RKD, any(c_double == RKALL))
926 character(1, SK), intent(in), optional :: input(*)
927 type(c_funptr), intent(in), value :: getLogFunc
928 integer(IK), intent(in), value :: ndim
929 integer(IK) :: stat
930#define runParaDRAM_ENABLED 1
931#include "pm_sampling@routines.inc.F90"
932#undef runParaDRAM_ENABLED
933 end function
934
935 function runParaDRAMF(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMF")
936#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
937 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMF
938#endif
939 use pm_kind, only: RKALL
940 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_float
941 integer, parameter :: RKG = merge(c_float, RKL, any(c_float == RKALL))
942 character(1, SK), intent(in), optional :: input(*)
943 type(c_funptr), intent(in), value :: getLogFunc
944 integer(IK), intent(in), value :: ndim
945 integer(IK) :: stat
946#define runParaDRAM_ENABLED 1
947#include "pm_sampling@routines.inc.F90"
948#undef runParaDRAM_ENABLED
949 end function
951
952!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
953
954end module pm_sampling ! LCOV_EXCL_LINE
Generate and return .true. if the procedure fails to fully accomplish the task of Monte Carlo samplin...
Generate and return a scalar string resulting from the concatenation of the individual elements (sing...
Definition: pm_str.F90:1347
This module contains classes and procedures for reporting and handling errors.
Definition: pm_err.F90:52
character(*, SK), parameter MODULE_NAME
Definition: pm_err.F90:58
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK5
Definition: pm_kind.F90:478
integer, dimension(*), parameter RKALL
The integer vector containing all defined real kinds within the ParaMonte library.
Definition: pm_kind.F90:636
integer, parameter RK4
Definition: pm_kind.F90:489
integer, parameter RK2
Definition: pm_kind.F90:511
integer, parameter RK3
Definition: pm_kind.F90:500
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter RKL
The scalar integer constant of intrinsic default kind, representing the lowest-precision real kind ty...
Definition: pm_kind.F90:779
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKH
The scalar integer constant of intrinsic default kind, representing the highest-precision real kind t...
Definition: pm_kind.F90:858
integer, parameter RK1
Definition: pm_kind.F90:522
This module contains procedures and generic interfaces for the ParaMonte library sampler routines.
Definition: pm_sampling.F90:30
integer(IK) function runParaDRAML(getLogFunc, ndim, input)
integer(IK) function runParaDRAMF(getLogFunc, ndim, input)
integer(IK) function runParaDRAMD(getLogFunc, ndim, input)
This module contains classes and procedures for various string manipulations and inquiries.
Definition: pm_str.F90:49
This module contains classes and procedures for manipulating system file/folder paths.
Definition: pm_sysPath.F90:274
integer(IK), parameter MAX_LEN_FILE_PATH
The scalar integer constant of default kind IK, representing the maximum system-allowed path length (...
Definition: pm_sysPath.F90:294
This is the derived type for generating objects to gracefully and verbosely handle runtime unexpected...
Definition: pm_err.F90:157
This is a derived type for constructing objects containing the optional simulation properties of the ...
This is a derived type for constructing objects containing the optional simulation properties of the ...
This is a derived type for constructing objects containing the optional simulation properties of the ...
This is a derived type for constructing objects containing the optional simulation properties of the ...
This is a derived type for constructing objects containing the optional simulation properties of the ...
Definition: pm_sampling.F90:76