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
30
31!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
32
34
35 use pm_err, only: err_type
36 use pm_str, only: getCharSeq
37 use, intrinsic :: iso_c_binding, only: c_char, c_funptr, c_null_char, c_f_procpointer
38 use pm_kind, only: SK, IK, LK, RKL, RKD, RKH, RKALL, SKG => SK ! LCOV_EXCL_LINE
39 use pm_sysPath, only: PATHLEN => MAX_LEN_FILE_PATH
40
41 implicit none
42
43 character(*, SK), parameter :: MODULE_NAME = "@pm_sampling"
44
45!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46
81 character(:,SKG) , allocatable :: description
82 character(:,SKG) , allocatable :: domain
83 character(:,SKG) , allocatable :: domainAxisName(:)
84 real(RKH) , allocatable :: domainBallAvg(:)
85 real(RKH) , allocatable :: domainBallCor(:,:)
86 real(RKH) , allocatable :: domainBallCov(:,:)
87 real(RKH) , allocatable :: domainBallStd(:)
88 real(RKH) , allocatable :: domainCubeLimitLower(:)
89 real(RKH) , allocatable :: domainCubeLimitUpper(:)
90 integer(IK) , allocatable :: domainErrCount
91 integer(IK) , allocatable :: domainErrCountMax
92 character(:,SKG) , allocatable :: inputFile / / /
155 !logical(LK) , allocatable :: inputFileHasPriority !<
156 character(:,SKG) , allocatable :: outputChainFileFormat
157 integer(IK) , allocatable :: outputColumnWidth
158 character(:,SKG) , allocatable :: outputFileName
159 character(:,SKG) , allocatable :: outputStatus
160 integer(IK) , allocatable :: outputPrecision
161 integer(IK) , allocatable :: outputReportPeriod
162 character(:,SKG) , allocatable :: outputRestartFileFormat
163 integer(IK) , allocatable :: outputSampleSize
164 character(:,SKG) , allocatable :: outputSeparator
165 character(:,SKG) , allocatable :: outputSplashMode
166 character(:,SKG) , allocatable :: parallelism
167 logical(LK) , allocatable :: parallelismMpiFinalizeEnabled
168 integer(IK) , allocatable :: parallelismNumThread
169 integer(IK) , allocatable :: randomSeed
170 character(:,SKG) , allocatable :: sysInfoFilePath
171 real(RKH) , allocatable :: targetAcceptanceRate(:)
172 end type
173
209 type, extends(sampler_type) :: paramcmc_type
210 integer(IK) , allocatable :: outputChainSize
211 integer(IK) , allocatable :: outputSampleRefinementCount
212 character(:,SKG) , allocatable :: outputSampleRefinementMethod
213 character(:,SKG) , allocatable :: proposal
214 real(RKH) , allocatable :: proposalCor(:,:)
215 real(RKH) , allocatable :: proposalCov(:,:)
216 character(:,SKG) , allocatable :: proposalScale
217 real(RKH) , allocatable :: proposalStart(:)
218 real(RKH) , allocatable :: proposalStartDomainCubeLimitLower(:)
219 real(RKH) , allocatable :: proposalStartDomainCubeLimitUpper(:)
220 logical(LK) , allocatable :: proposalStartRandomized
221 real(RKH) , allocatable :: proposalStd(:)
222 end type
223
250 type, extends(paramcmc_type) :: paradram_type
251 real(RKH) , allocatable :: proposalAdaptationBurnin
252 integer(IK) , allocatable :: proposalAdaptationCount
253 integer(IK) , allocatable :: proposalAdaptationCountGreedy
254 integer(IK) , allocatable :: proposalAdaptationPeriod
255 integer(IK) , allocatable :: proposalDelayedRejectionCount
256 real(RKH) , allocatable :: proposalDelayedRejectionScale(:)
257 end type
258
285 type, extends(paradram_type) :: paradise_type
286 end type
287
315 type, extends(sampler_type) :: paranest_type
316 integer(IK) , allocatable :: domainPartitionAdaptationCount
317 integer(IK) , allocatable :: domainPartitionAdaptationPeriod
318 logical(LK) , allocatable :: domainPartitionBiasCorrectionEnabled
319 integer(IK) , allocatable :: domainPartitionCountMax
320 real(RKH) , allocatable :: domainPartitionFactorExpansion
321 real(RKH) , allocatable :: domainPartitionFactorShrinkage
322 integer(IK) , allocatable :: domainPartitionKmeansClusterCountMax
323 integer(IK) , allocatable :: domainPartitionKmeansClusterSizeMin
324 logical(LK) , allocatable :: domainPartitionKmeansNormalizationEnabled
325 integer(IK) , allocatable :: domainPartitionKmeansNumFailMax
326 integer(IK) , allocatable :: domainPartitionKmeansNumRecursionMax
327 integer(IK) , allocatable :: domainPartitionKmeansNumTry
328 integer(IK) , allocatable :: domainPartitionKvolumeNumRecursionMax
329 real(RKH) , allocatable :: domainPartitionKvolumeWeightExponent
330 character(:,SKG) , allocatable :: domainPartitionMethod
331 character(:,SKG) , allocatable :: domainPartitionObject
332 logical(LK) , allocatable :: domainPartitionOptimizationScaleEnabled
333 logical(LK) , allocatable :: domainPartitionOptimizationShapeEnabled
334 logical(LK) , allocatable :: domainPartitionOptimizationShapeScaleEnabled
335 real(RKH) , allocatable :: domainPartitionScaleFactor
336 character(:,SKG) , allocatable :: domainSampler
337 integer(IK) , allocatable :: liveSampleSize
338 real(RKH) , allocatable :: tolerance
339 end type
340
341!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
342
343 ! This abstract interface defines the set of objective function interfaces
344 ! that can be passed to the ParaMonte samplers (e.g., ParaDRAM, ParaNest, ...).
345 ! These interfaces only have internal library specific use cases.
347 abstract interface
348
349#if OMP_ENABLED && (MATLAB_ENABLED || PYTHON_ENABLED || R_ENABLED)
350
351 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
352
353#if RK5_ENABLED
354 function getLogFunc_proc_RK5(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
355#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
356 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK5
357#endif
358 use pm_kind, only: RKD, RKG => RK5
359 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
360 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
361 real(RKG) :: mold
362 end function
363#endif
364
365#if RK4_ENABLED
366 function getLogFunc_proc_RK4(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
367#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
368 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK4
369#endif
370 use pm_kind, only: RKD, RKG => RK4
371 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
372 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
373 real(RKG) :: mold
374 end function
375#endif
376
377#if RK3_ENABLED
378 function getLogFunc_proc_RK3(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
379#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
380 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK3
381#endif
382 use pm_kind, only: RKD, RKG => RK3
383 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
384 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
385 real(RKG) :: mold
386 end function
387#endif
388
389#if RK2_ENABLED
390 function getLogFunc_proc_RK2(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
391#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
392 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK2
393#endif
394 use pm_kind, only: RKD, RKG => RK2
395 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
396 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
397 real(RKG) :: mold
398 end function
399#endif
400
401#if RK1_ENABLED
402 function getLogFunc_proc_RK1(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
403#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
404 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK1
405#endif
406 use pm_kind, only: RKD, RKG => RK1
407 real(RKG), intent(inout), contiguous :: logFuncState(0:,:)
408 real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
409 real(RKG) :: mold
410 end function
411#endif
412
413 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
414
415#else
416
417 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
418
419#if RK5_ENABLED
420 function getLogFunc_proc_RK5(state) result(logFunc)
421#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
422 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK5
423#endif
424 use pm_kind, only: RKG => RK5
425 real(RKG), intent(in), contiguous :: state(:)
426 real(RKG) :: logFunc
427 end function
428#endif
429
430#if RK4_ENABLED
431 function getLogFunc_proc_RK4(state) result(logFunc)
432#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
433 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK4
434#endif
435 use pm_kind, only: RKG => RK4
436 real(RKG), intent(in), contiguous :: state(:)
437 real(RKG) :: logFunc
438 end function
439#endif
440
441#if RK3_ENABLED
442 function getLogFunc_proc_RK3(state) result(logFunc)
443#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
444 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK3
445#endif
446 use pm_kind, only: RKG => RK3
447 real(RKG), intent(in), contiguous :: state(:)
448 real(RKG) :: logFunc
449 end function
450#endif
451
452#if RK2_ENABLED
453 function getLogFunc_proc_RK2(state) result(logFunc)
454#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
455 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK2
456#endif
457 use pm_kind, only: RKG => RK2
458 real(RKG), intent(in), contiguous :: state(:)
459 real(RKG) :: logFunc
460 end function
461#endif
462
463#if RK1_ENABLED
464 function getLogFunc_proc_RK1(state) result(logFunc)
465#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
466 !DEC$ ATTRIBUTES DLLEXPORT :: getLogFunc_proc_RK1
467#endif
468 use pm_kind, only: RKG => RK1
469 real(RKG), intent(in), contiguous :: state(:)
470 real(RKG) :: logFunc
471 end function
472#endif
473
474 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
475
476#endif
477
478 end interface
480
481!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
482
675
676 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
677
678#if RK5_ENABLED
679 impure module function getErrParaDRAM_RK5(sampler, getLogFunc, ndim) result(err)
680#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
681 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK5
682#endif
683 use pm_kind, only: RKG => RK5
684 procedure(getLogFunc_proc_RK5) :: getLogFunc
685 type(paradram_type), intent(in) :: sampler
686 integer(IK), intent(in) :: ndim
687 type(err_type) :: err
688 end function
689#endif
690
691#if RK4_ENABLED
692 impure module function getErrParaDRAM_RK4(sampler, getLogFunc, ndim) result(err)
693#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
694 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK4
695#endif
696 use pm_kind, only: RKG => RK4
697 procedure(getLogFunc_proc_RK4) :: getLogFunc
698 type(paradram_type), intent(in) :: sampler
699 integer(IK), intent(in) :: ndim
700 type(err_type) :: err
701 end function
702#endif
703
704#if RK3_ENABLED
705 impure module function getErrParaDRAM_RK3(sampler, getLogFunc, ndim) result(err)
706#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
707 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK3
708#endif
709 use pm_kind, only: RKG => RK3
710 procedure(getLogFunc_proc_RK3) :: getLogFunc
711 type(paradram_type), intent(in) :: sampler
712 integer(IK), intent(in) :: ndim
713 type(err_type) :: err
714 end function
715#endif
716
717#if RK2_ENABLED
718 impure module function getErrParaDRAM_RK2(sampler, getLogFunc, ndim) result(err)
719#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
720 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK2
721#endif
722 use pm_kind, only: RKG => RK2
723 procedure(getLogFunc_proc_RK2) :: getLogFunc
724 type(paradram_type), intent(in) :: sampler
725 integer(IK), intent(in) :: ndim
726 type(err_type) :: err
727 end function
728#endif
729
730#if RK1_ENABLED
731 impure module function getErrParaDRAM_RK1(sampler, getLogFunc, ndim) result(err)
732#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
733 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDRAM_RK1
734#endif
735 use pm_kind, only: RKG => RK1
736 procedure(getLogFunc_proc_RK1) :: getLogFunc
737 type(paradram_type), intent(in) :: sampler
738 integer(IK), intent(in) :: ndim
739 type(err_type) :: err
740 end function
741#endif
742
743 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
744
745#if RK5_ENABLED
746 impure module function getErrParaDISE_RK5(sampler, getLogFunc, ndim) result(err)
747#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
748 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK5
749#endif
750 use pm_kind, only: RKG => RK5
751 procedure(getLogFunc_proc_RK5) :: getLogFunc
752 type(paradise_type), intent(in) :: sampler
753 integer(IK), intent(in) :: ndim
754 type(err_type) :: err
755 end function
756#endif
757
758#if RK4_ENABLED
759 impure module function getErrParaDISE_RK4(sampler, getLogFunc, ndim) result(err)
760#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
761 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK4
762#endif
763 use pm_kind, only: RKG => RK4
764 procedure(getLogFunc_proc_RK4) :: getLogFunc
765 type(paradise_type), intent(in) :: sampler
766 integer(IK), intent(in) :: ndim
767 type(err_type) :: err
768 end function
769#endif
770
771#if RK3_ENABLED
772 impure module function getErrParaDISE_RK3(sampler, getLogFunc, ndim) result(err)
773#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
774 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK3
775#endif
776 use pm_kind, only: RKG => RK3
777 procedure(getLogFunc_proc_RK3) :: getLogFunc
778 type(paradise_type), intent(in) :: sampler
779 integer(IK), intent(in) :: ndim
780 type(err_type) :: err
781 end function
782#endif
783
784#if RK2_ENABLED
785 impure module function getErrParaDISE_RK2(sampler, getLogFunc, ndim) result(err)
786#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
787 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK2
788#endif
789 use pm_kind, only: RKG => RK2
790 procedure(getLogFunc_proc_RK2) :: getLogFunc
791 type(paradise_type), intent(in) :: sampler
792 integer(IK), intent(in) :: ndim
793 type(err_type) :: err
794 end function
795#endif
796
797#if RK1_ENABLED
798 impure module function getErrParaDISE_RK1(sampler, getLogFunc, ndim) result(err)
799#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
800 !DEC$ ATTRIBUTES DLLEXPORT :: getErrParaDISE_RK1
801#endif
802 use pm_kind, only: RKG => RK1
803 procedure(getLogFunc_proc_RK1) :: getLogFunc
804 type(paradise_type), intent(in) :: sampler
805 integer(IK), intent(in) :: ndim
806 type(err_type) :: err
807 end function
808#endif
809
810 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
811
812 end interface
813
814!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
815
816contains
817
818!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
819
991 function runParaDRAML(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAML")
992#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
993 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAML
994#endif
995 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_long_double
996 integer, parameter :: RKG = merge(c_long_double, RKH, any(c_long_double == RKALL))
997 character(1, SK), intent(in), optional :: input(*)
998 type(c_funptr), intent(in), value :: getLogFunc
999 integer(IK), intent(in), value :: ndim
1000 integer(IK) :: stat
1001#define runParaDRAM_ENABLED 1
1002#include "pm_sampling@routines.inc.F90"
1003#undef runParaDRAM_ENABLED
1004 end function
1005
1006 function runParaDRAMD(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMD")
1007#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1008 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMD
1009#endif
1010 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_double
1011 integer, parameter :: RKG = merge(c_double, RKD, any(c_double == RKALL))
1012 character(1, SK), intent(in), optional :: input(*)
1013 type(c_funptr), intent(in), value :: getLogFunc
1014 integer(IK), intent(in), value :: ndim
1015 integer(IK) :: stat
1016#define runParaDRAM_ENABLED 1
1017#include "pm_sampling@routines.inc.F90"
1018#undef runParaDRAM_ENABLED
1019 end function
1020
1021 function runParaDRAMF(getLogFunc, ndim, input) result(stat) bind(C, name = "runParaDRAMF")
1022#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
1023 !DEC$ ATTRIBUTES DLLEXPORT :: runParaDRAMF
1024#endif
1025 use pm_kind, only: RKALL
1026 use, intrinsic :: iso_c_binding, only: SK => c_char, IK => c_int32_t, c_float
1027 integer, parameter :: RKG = merge(c_float, RKL, any(c_float == RKALL))
1028 character(1, SK), intent(in), optional :: input(*)
1029 type(c_funptr), intent(in), value :: getLogFunc
1030 integer(IK), intent(in), value :: ndim
1031 integer(IK) :: stat
1032#define runParaDRAM_ENABLED 1
1033#include "pm_sampling@routines.inc.F90"
1034#undef runParaDRAM_ENABLED
1035 end function
1037
1038!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1039
1040end 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:33
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:80