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 : ! \bug
18 : ! Avoid Intel ifort bug for too many `use` statements in a submodule by placing them all in the submodule header.
19 : #if CHECK_ENABLED
20 : use pm_err, only: setAsserted
21 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG)call setAsserted(ASSERTION,getFine(__FILE__,LINE)//MODULE_NAME//MSG);
22 : #else
23 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG) continue;
24 : #endif
25 : ! Define intel-specific shared property of files on Windows.
26 : #if INTEL_ENABLED && WINDOWS_ENABLED
27 : #define SHARED, shared
28 : #else
29 : #define SHARED
30 : #endif
31 :
32 : use pm_err, only: getFine
33 : use pm_err, only: err_type
34 : use pm_kind, only: SK, IK, LK, SKC => SK
35 : use pm_sysPath, only: PATHLEN => MAX_LEN_FILE_PATH
36 : use pm_strASCII, only: getStrLower, setStrLower
37 : use pm_arrayReplace, only: getReplaced
38 : use pm_io, only: setContentsFrom
39 : use pm_val2str, only: getStr
40 :
41 : implicit none
42 :
43 : character(*, SK), parameter :: MODULE_NAME = "@pm_sampling_scio"
44 :
45 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46 :
47 : ! The namelist specifications of this module are all internal ParaMonte libraries and not meant to be accessible to end users.
48 : ! This utter io mess is entirely due to shear inadequacy of the namelist feature and
49 : ! the standard committee inaction to expand and improve this valuable feature.
50 : ! Unfortunately, the following variables must remain public, because they are
51 : ! used by the following internal modules:
52 : ! pm_sampling_base
53 : ! pm_sampling_mcmc
54 : ! pm_sampling_dram
55 : ! pm_sampling_nest
56 :
57 : ! specbase
58 : !> \cond excluded
59 : character(8191,SKC) :: description ! roughly 66Kb of memory
60 : character(63,SKC) :: domain
61 : character(63,SKC) , allocatable :: domainAxisName(:)
62 : real(RKC) , allocatable :: domainBallCenter(:)
63 : real(RKC) , allocatable :: domainBallCorMat(:,:)
64 : real(RKC) , allocatable :: domainBallCovMat(:,:)
65 : real(RKC) , allocatable :: domainBallStdVec(:)
66 : real(RKC) , allocatable :: domainCubeLimitLower(:)
67 : real(RKC) , allocatable :: domainCubeLimitUpper(:)
68 : integer(IK) :: domainErrCount
69 : integer(IK) :: domainErrCountMax
70 : logical(LK) :: inputFileHasPriority
71 : character(15,SKC) :: outputChainFileFormat
72 : integer(IK) :: outputColumnWidth
73 : character(PATHLEN,SKC) :: outputFileName
74 : character(15,SKC) :: outputStatus
75 : integer(IK) :: outputPrecision
76 : integer(IK) :: outputReportPeriod
77 : character(15,SKC) :: outputRestartFileFormat
78 : integer(IK) :: outputSampleSize
79 : character(63,SKC) :: outputSeparator
80 : character(63,SKC) :: outputSplashMode
81 : character(63,SKC) :: parallelism
82 : logical(LK) :: parallelismMpiFinalizeEnabled
83 : integer(IK) :: parallelismNumThread
84 : !character(511,SKC) :: plang
85 : integer(IK) :: randomSeed
86 : character(PATHLEN,SKC) :: sysInfoFilePath
87 : real(RKC) :: targetAcceptanceRate(2)
88 : !> \endcond excluded
89 :
90 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
91 :
92 : ! specmcmc
93 : !> \cond excluded
94 : integer(IK) :: outputChainSize
95 : integer(IK) :: outputSampleRefinementCount
96 : character(63,SKC) :: outputSampleRefinementMethod
97 : character(63,SKC) :: proposal
98 : real(RKC) , allocatable :: proposalCorMat(:,:)
99 : real(RKC) , allocatable :: proposalCovMat(:,:)
100 : character(127,SKC) :: proposalScaleFactor
101 : real(RKC) , allocatable :: proposalStart(:)
102 : real(RKC) , allocatable :: proposalStartDomainCubeLimitLower(:)
103 : real(RKC) , allocatable :: proposalStartDomainCubeLimitUpper(:)
104 : logical(LK) :: proposalStartRandomized
105 : real(RKC) , allocatable :: proposalStdVec(:)
106 : !> \endcond excluded
107 :
108 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109 :
110 : ! specdram
111 : !> \cond excluded
112 : real(RKC) :: burninAdaptationMeasure
113 : integer(IK) :: proposalAdaptationCount
114 : integer(IK) :: proposalAdaptationCountGreedy
115 : integer(IK) :: proposalAdaptationPeriod
116 : integer(IK) :: proposalDelayedRejectionCount
117 : real(RKC) , allocatable :: proposalDelayedRejectionScaleFactor(:)
118 : !> \endcond excluded
119 :
120 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
121 :
122 : ! specnest
123 : !> \cond excluded
124 : integer(IK) :: domainPartitionAdaptationCount
125 : integer(IK) :: domainPartitionAdaptationPeriod
126 : logical(LK) :: domainPartitionBiasCorrectionEnabled
127 : integer(IK) :: domainPartitionCountMax
128 : real(RKC) :: domainPartitionFactorExpansion
129 : real(RKC) :: domainPartitionFactorShrinkage
130 : integer(IK) :: domainPartitionKmeansClusterCountMax
131 : integer(IK) :: domainPartitionKmeansClusterSizeMin
132 : logical(LK) :: domainPartitionKmeansNormalizationEnabled
133 : integer(IK) :: domainPartitionKmeansNumFailMax
134 : integer(IK) :: domainPartitionKmeansNumRecursionMax
135 : integer(IK) :: domainPartitionKmeansNumTry
136 : integer(IK) :: domainPartitionKvolumeNumRecursionMax
137 : real(RKC) :: domainPartitionKvolumeWeightExponent
138 : character(63,SKC) :: domainPartitionMethod
139 : character(63,SKC) :: domainPartitionObject
140 : logical(LK) :: domainPartitionOptimizationScaleEnabled
141 : logical(LK) :: domainPartitionOptimizationShapeEnabled
142 : logical(LK) :: domainPartitionOptimizationShapeScaleEnabled
143 : real(RKC) :: domainPartitionScaleFactor
144 : character(63,SKC) :: domainSampler
145 : integer(IK) :: liveSampleSize
146 : real(RKC) :: tolerance
147 : !> \endcond excluded
148 :
149 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
150 :
151 : #define namelist_ENABLED 1
152 :
153 : #define NAMELIST paramonte
154 : #include "pm_sampling_scio.inc.F90"
155 : #undef NAMELIST
156 :
157 : #define NAMELIST paradram
158 : #include "pm_sampling_scio.inc.F90"
159 : #undef NAMELIST
160 :
161 : #define NAMELIST paradise
162 : #include "pm_sampling_scio.inc.F90"
163 : #undef NAMELIST
164 :
165 : #define NAMELIST paranest
166 : #include "pm_sampling_scio.inc.F90"
167 : #undef NAMELIST
168 :
169 : #undef namelist_ENABLED
170 :
171 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
172 :
173 : type :: cfhbase_type
174 : character(31,SKC) :: processID = "processID"
175 : character(31,SKC) :: meanAcceptanceRate = "meanAcceptanceRate"
176 : character(31,SKC) :: sampleLogFunc = "sampleLogFunc"
177 : end type
178 : type(cfhbase_type) , parameter :: cfhbase = cfhbase_type()
179 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
180 : !DIR$ ATTRIBUTES DLLEXPORT :: cfhbase
181 : #endif
182 :
183 : !> \brief
184 : !> This is the `abstract` derived type containing the components for storing the basic properties of output chain files by [ParaMonte samplers](@ref pm_sampling).<br>
185 : !>
186 : !> \see
187 : !> [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)<br>
188 : !> [cfcbase_type](@ref pm_sampling_scio::cfcbase_type)<br>
189 : !>
190 : !> \test
191 : !> [test_pm_sampling_scio](@ref test_pm_sampling_scio)<br>
192 : !>
193 : !> \bug
194 : !> \status \unresolved
195 : !> \source \gfortran{13}
196 : !> \desc
197 : !> Unfortunately, the lack of a fully functional implementation of PDTs in \gfortran{13} prevents the declaration of this derived type as a PDT.<br>
198 : !> \remedy
199 : !> See the mess with kind-specific sampling modules.
200 : !>
201 : !> \todo
202 : !> \pvhigh
203 : !> This derived type must be converted to a PDT as soon as PDTs are supported bug-free in \gfortran.<br>
204 : !>
205 : !> \finmain{cfcbase_type}
206 : !>
207 : !> \author
208 : !> \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
209 : type, abstract :: cfcbase_type
210 : integer(IK) :: numDefCol !< \public The scalar of type `integer` of default kind \IK. containing the number of the set of default columns that always appear in the chain file before columns corresponding to `sampleState` components appear.
211 : integer(IK) :: nsam !< \public The scalar of type `integer` of default kind \IK containing the chain size which may be less than the size of the chain components below (the number of (potentially weighted) sampled states).<br>
212 : character(:,SKC) , allocatable :: sep !< \public The scalar of type `character` of default kind \SK, containing the field separator used in the chain file.<br>
213 : character(:,SKC) , allocatable :: header !< \public The scalar of type `character` of default kind \SK, containing the delimited chain file column names.<br>
214 : !type(css_type) , allocatable :: colname(:) !< \public The vector of type [css_type](@ref pm_container::css_type), containing all chain file column names, such that the column name `sampleLogFunc` correspond to index `0` of colnames.<br>
215 : integer(IK) , allocatable :: processID(:) !< \public The `allocatable` vector of type `integer` of default kind \IK, containing the IDs of the processes contributing the corresponding accepted state in the chain file.<br>
216 : real(RKC) , allocatable :: meanAcceptanceRate(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the vector of the average acceptance rates at the given point in the chain.<br>
217 : real(RKC) , allocatable :: sampleLogFunc(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the natural logarithm of the function value corresponding to the state in the same entry in the chain file.<br>
218 : real(RKC) , allocatable :: sampleState(:,:) !< \public The `allocatable` matrix of type `real` of kind \RKALL, each `(1 : ndim, i)` slice of which contains the accepted `ndim`-dimensional state in the chain at the entry `i` of the chain file.<br>
219 : end type
220 :
221 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
222 :
223 : type, extends(cfhbase_type) :: cfhmcmc_type
224 : character(31,SKC) :: burninLocation = "burninLocation"
225 : character(31,SKC) :: sampleWeight = "sampleWeight"
226 : end type
227 : type(cfhmcmc_type) , parameter :: cfhmcmc = cfhmcmc_type()
228 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
229 : !DIR$ ATTRIBUTES DLLEXPORT :: cfhmcmc
230 : #endif
231 :
232 : !type :: size_type
233 : ! integer(IK) :: compact, verbose
234 : !end type
235 :
236 : !> \brief
237 : !> This is the `abstract` derived type containing the components for storing the MCMC properties of output chain files by [ParaMonte samplers](@ref pm_sampling).<br>
238 : !>
239 : !> \see
240 : !> [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)<br>
241 : !> [cfcbase_type](@ref pm_sampling_scio::cfcbase_type)<br>
242 : !>
243 : !> \test
244 : !> [test_pm_sampling_scio](@ref test_pm_sampling_scio)<br>
245 : !>
246 : !> \bug
247 : !> \status \unresolved
248 : !> \source \gfortran{13}
249 : !> \desc
250 : !> Unfortunately, the lack of a fully functional implementation of PDTs in \gfortran{13} prevents the declaration of this derived type as a PDT.<br>
251 : !> \remedy
252 : !> See the mess with kind-specific sampling modules.
253 : !>
254 : !> \todo
255 : !> \pvhigh
256 : !> This derived type must be converted to a PDT as soon as PDTs are supported bug-free in \gfortran.<br>
257 : !>
258 : !> \finmain{cfcbase_type}
259 : !>
260 : !> \author
261 : !> \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
262 : type, extends(cfcbase_type) :: cfcmcmc_type
263 : integer(IK) , allocatable :: burninLocation(:) !< \public The `allocatable` vector of type `integer` of default kind \IK, containing the burnin locations at the given locations in the chains.<br>
264 : integer(IK) , allocatable :: sampleWeight(:) !< \public The `allocatable` vector of type `integer` of default kind \IK, containing the weights of the accepted states.<br>
265 : integer(IK) :: sumw !< \public The sum of all sample weights (corresponding to the verbose chain size).<br>
266 : end type
267 :
268 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
269 :
270 : type, extends(cfhmcmc_type) :: cfhdram_type
271 : character(31,SKC) :: delayedRejectionStage = "delayedRejectionStage"
272 : character(31,SKC) :: adaptationMeasure = "adaptationMeasure"
273 : end type
274 : type(cfhdram_type) , parameter :: cfhdram = cfhdram_type()
275 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
276 : !DIR$ ATTRIBUTES DLLEXPORT :: cfhdram
277 : #endif
278 :
279 : character(31,SKC) , parameter :: chainFileColNameDRAM(*) =[ cfhdram%processID &
280 : , cfhdram%delayedRejectionStage &
281 : , cfhdram%meanAcceptanceRate &
282 : , cfhdram%adaptationMeasure &
283 : , cfhdram%burninLocation &
284 : , cfhdram%sampleWeight &
285 : , cfhdram%sampleLogFunc &
286 : ]
287 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
288 : !DIR$ ATTRIBUTES DLLEXPORT :: chainFileColNameDRAM
289 : #endif
290 :
291 : !> \brief
292 : !> This is the `abstract` derived type containing the components for storing the DRAM properties of output chain files by [ParaMonte samplers](@ref pm_sampling).<br>
293 : !>
294 : !> \see
295 : !> [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)<br>
296 : !> [cfcbase_type](@ref pm_sampling_scio::cfcbase_type)<br>
297 : !>
298 : !> \test
299 : !> [test_pm_sampling_scio](@ref test_pm_sampling_scio)<br>
300 : !>
301 : !> \bug
302 : !> \status \unresolved
303 : !> \source \gfortran{13}
304 : !> \desc
305 : !> Unfortunately, the lack of a fully functional implementation of PDTs in \gfortran{13} prevents the declaration of this derived type as a PDT.<br>
306 : !> \remedy
307 : !> See the mess with kind-specific sampling modules.
308 : !>
309 : !> \todo
310 : !> \pvhigh
311 : !> This derived type must be converted to a PDT as soon as PDTs are supported bug-free in \gfortran.<br>
312 : !>
313 : !> \finmain{cfcbase_type}
314 : !>
315 : !> \author
316 : !> \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
317 : type, extends(cfcmcmc_type) :: cfcdram_type
318 : real(RKC) , allocatable :: adaptationMeasure(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the vector of the adaptation measures at the accepted states.<br>
319 : integer(IK) , allocatable :: delayedRejectionStage(:) !< \public The `allocatable` vector of type `integer` of default kind \IK, containing the delayed rejection stages at which the proposed states were accepted.<br>
320 : end type
321 :
322 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
323 :
324 : type, extends(cfhdram_type) :: cfhdise_type
325 : end type
326 : type(cfhdise_type) , parameter :: cfhdise = cfhdise_type()
327 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
328 : !DIR$ ATTRIBUTES DLLEXPORT :: cfhdise
329 : #endif
330 :
331 : character(31,SKC) , parameter :: chainFileColNameDISE(*) =[ cfhdise%processID &
332 : , cfhdise%delayedRejectionStage &
333 : , cfhdise%meanAcceptanceRate &
334 : , cfhdise%adaptationMeasure &
335 : , cfhdise%burninLocation &
336 : , cfhdise%sampleWeight &
337 : , cfhdise%sampleLogFunc &
338 : ]
339 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
340 : !DIR$ ATTRIBUTES DLLEXPORT :: chainFileColNameDISE
341 : #endif
342 :
343 : !> \brief
344 : !> This is the `abstract` derived type containing the components for storing the DISE properties of output chain files by [ParaMonte samplers](@ref pm_sampling).<br>
345 : !>
346 : !> \see
347 : !> [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)<br>
348 : !> [cfcbase_type](@ref pm_sampling_scio::cfcbase_type)<br>
349 : !>
350 : !> \test
351 : !> [test_pm_sampling_scio](@ref test_pm_sampling_scio)<br>
352 : !>
353 : !> \bug
354 : !> \status \unresolved
355 : !> \source \gfortran{13}
356 : !> \desc
357 : !> Unfortunately, the lack of a fully functional implementation of PDTs in \gfortran{13} prevents the declaration of this derived type as a PDT.<br>
358 : !> \remedy
359 : !> See the mess with kind-specific sampling modules.
360 : !>
361 : !> \todo
362 : !> \pvhigh
363 : !> This derived type must be converted to a PDT as soon as PDTs are supported bug-free in \gfortran.<br>
364 : !>
365 : !> \finmain{cfcbase_type}
366 : !>
367 : !> \author
368 : !> \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
369 : type, extends(cfcdram_type) :: cfcdise_type
370 : end type
371 :
372 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
373 :
374 : type, extends(cfhbase_type) :: cfhnest_type
375 : character(31,SKC) :: domainLogVol = "domainLogVol"
376 : character(31,SKC) :: funcLogIntegral = "funcLogIntegral"
377 : character(31,SKC) :: sampleLogWeight = "sampleLogWeight"
378 : end type
379 : type(cfhnest_type) , parameter :: cfhnest = cfhnest_type()
380 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
381 : !DIR$ ATTRIBUTES DLLEXPORT :: cfhnest
382 : #endif
383 :
384 : character(31,SKC) , parameter :: chainFileColNameNest(*) =[ cfhnest%processID &
385 : , cfhnest%meanAcceptanceRate &
386 : , cfhnest%domainLogVol &
387 : , cfhnest%funcLogIntegral &
388 : , cfhnest%sampleLogWeight &
389 : , cfhnest%sampleLogFunc &
390 : ]
391 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
392 : !DIR$ ATTRIBUTES DLLEXPORT :: chainFileColNameNest
393 : #endif
394 :
395 : !> \brief
396 : !> This is the derived type containing `allocatable` storage for columns of data in output chain files resulting from [ParaNest](@ref pm_sampling::paranest_type) simulations.
397 : !>
398 : !> \details
399 : !> Objects of this derived type can be passed to the generic interface [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)
400 : !> to read the contents of a given [ParaNest](@ref pm_sampling::paranest_type) simulation chain file.<br>
401 : !> For more details, see the documentation of the components of the derived type.<br>
402 : !>
403 : !> \note
404 : !> This derived type does not have a custom constructor as of the most recent version of the ParaMonte library.<br>
405 : !>
406 : !> \see
407 : !> [getErrChainRead](@ref pm_sampling_scio::getErrChainRead)<br>
408 : !> [cfcbase_type](@ref pm_sampling_scio::cfcbase_type)<br>
409 : !>
410 : !> \test
411 : !> [test_pm_sampling_scio](@ref test_pm_sampling_scio)<br>
412 : !>
413 : !> \bug
414 : !> \status \unresolved
415 : !> \source \gfortran{13}
416 : !> \desc
417 : !> Unfortunately, the lack of a fully functional implementation of PDTs in \gfortran{13} prevents the declaration of this derived type as a PDT.<br>
418 : !> \remedy
419 : !> See the mess with kind-specific sampling modules.
420 : !>
421 : !> \todo
422 : !> \pvhigh
423 : !> This derived type must be converted to a PDT as soon as PDTs are supported bug-free in \gfortran.<br>
424 : !>
425 : !> \finmain{chainParaNest_type}
426 : !>
427 : !> \author
428 : !> \AmirShahmoradi, September 1, 2012, 12:00 AM, National Institute for Fusion Studies, The University of Texas at Austin
429 : type, extends(cfcbase_type) :: cfcnest_type
430 : real(RKC) , allocatable :: domainLogVol(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the natural logarithm of the size of the domain of integration at the corresponding stage in the chain of sampling/integration.<br>
431 : real(RKC) , allocatable :: logMaxRelErr(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the natural logarithm of the size of the domain of integration at the corresponding stage in the chain of sampling/integration.<br>
432 : real(RKC) , allocatable :: funcLogIntegral(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the natural logarithm of the size of the integral of the objective function at the corresponding stage in the chain.<br>
433 : real(RKC) , allocatable :: sampleLogWeight(:) !< \public The `allocatable` vector of type `real` of kind \RKALL, containing the natural logarithm of the contribution of the corresponding state in chain to the integral of the objective function.<br>
434 : end type
435 :
436 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
437 :
438 : contains
439 :
440 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
441 :
442 : !> \brief
443 : !> Set the module namelist variables declared in the parent module of this generic interface from the specified `inputFile`.<br>
444 : !>
445 : !> \details
446 : !> This is a low-level ParaMonte library interface that is meant to be directly accessible to the end users.<br>
447 : !> This routines of this generic interface are used by various explorers and samplers of the ParaMonte library
448 : !> to parse the user-specified simulation specifications.<br>
449 : !>
450 : !> \param[in] method : The input scalar object that can be
451 : !> <ol>
452 : !> <li> The constant `"ParaDISE"` indicating the use of the ParaDISE sampler,<br>
453 : !> <li> The constant `"ParaDRAM"` indicating the use of the ParaDRAM sampler,<br>
454 : !> <li> The constant `"ParaNest"` indicating the use of the ParaNest sampler,<br>
455 : !> </ol>
456 : !> \param[in] inputFile : The input scalar of type `character` of default kind \SK, containing either
457 : !> <ol>
458 : !> <li> The path to the external input file containing a ParaMonte simulation specifications.<br>
459 : !> <li> A string containing the simulation specifications as if it contains the contents of
460 : !> an external input simulation specification file in its entirety.
461 : !> </ol>
462 : !> The simulation specifications within the input file or string
463 : !> must follow the rules of the Fortran language `namelist` IO.<br>
464 : !> Four namelist grouping names are currently supported:<br>
465 : !> <ol>
466 : !> <li> The grouping name `&ParaDISE` (case-<b>insensitive</b>) signifies
467 : !> a specification grouping for the ParaDISE sampler and optimizer.<br>
468 : !> When the user-specified `method` is `"ParaDISE"`,
469 : !> this will be the default grouping name that will be search within the input file.<br>
470 : !> <li> The grouping name `&ParaDRAM` (case-<b>insensitive</b>) signifies
471 : !> a specification grouping for the ParaDRAM sampler and integrator.<br>
472 : !> When the user-specified `method` is `"ParaDRAM"`,
473 : !> this will be the default grouping name that will be search within the input file.<br>
474 : !> <li> The grouping name `&ParaNest` (case-<b>insensitive</b>) signifies
475 : !> a specification grouping for the ParaNest sampler and integrator.<br>
476 : !> When the user-specified `method` is `"ParaNest"`,
477 : !> this will be the default grouping name that will be search within the input file.<br>
478 : !> <li> The grouping name `ParaMonte` (case-<b>insensitive</b>) specifies a specification grouping
479 : !> for any simulator/explorer supported by the ParaMonte library.<br>
480 : !> This group naming is used as the last resort to read the simulation specifications
481 : !> only if all of the above simulation-specific namelist groups are missing in the input file.<br>
482 : !> </ol>
483 : !> \param[inout] err : The input/output scalar of type [err_type](@ref pm_err::err_type) containing
484 : !> the relevant error information if any error occurs the IO within the routine.<br>
485 : !> On output, if any error has occurred, an error message will be appended to the
486 : !> `msg` component of `err` and the `occurred` component is set to `.true.`.<br>
487 : !> Otherwise, the input contents of `err` remain intact on output.<br>
488 : !>
489 : !> \interface{setSpecFromInput}
490 : !> \code{.F90}
491 : !>
492 : !> use pm_sampling_scio, only: setSpecFromInput
493 : !>
494 : !> call setSpecFromInput(method, inputFile, err)
495 : !>
496 : !> \endcode
497 : !>
498 : !> \warnpure
499 : !>
500 : !> \see
501 : !> [pm_sampling](@ref pm_sampling)<br>
502 : !> [pm_sampling_base](@ref pm_sampling_base)<br>
503 : !> [pm_sampling_mcmc](@ref pm_sampling_mcmc)<br>
504 : !> [pm_sampling_dram](@ref pm_sampling_dram)<br>
505 : !> [pm_sampling_nest](@ref pm_sampling_nest)<br>
506 : !>
507 : !> \finmain{setSpecFromInput}
508 : !>
509 : !> \author
510 : !> \FatemehBagheri, Monday 02:15 AM, September 27, 2021, Dallas, TX<br>
511 3 : subroutine setSpecFromInput(method, inputFile, err)
512 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
513 : !DEC$ ATTRIBUTES DLLEXPORT :: setSpecFromInput
514 : #endif
515 : character(*,SKC) , intent(in) :: inputFile
516 : character(*,SKC) , intent(in) :: method
517 : type(err_type) , intent(inout) :: err
518 : logical(LK) :: opened
519 : integer(IK) :: unit
520 : integer :: iostat
521 : !logical(LK) :: missing
522 : !logical(LK) :: specIsNameList
523 : logical(LK) :: exist
524 : character(255, SK) :: iomsg
525 : character(:, SK), allocatable :: lcmethod
526 3 : character(:, SK), allocatable :: contents
527 : character(*, SK), parameter :: PROCEDURE_NAME = MODULE_NAME//SKC_"@setSpecFromInput()"
528 : #define RETURN_IF_ERRED(LINE,FAILED,MSG) \
529 : if (FAILED) then; \
530 : err%msg = PROCEDURE_NAME//getFine(__FILE__, LINE)//SK_": "//trim(MSG); \
531 : err%occurred = .true._LK; \
532 : err%stat = iostat; \
533 : return; \
534 : end if;
535 : ! Determine if the file is external vs. internal.
536 :
537 3 : exist = .false.
538 3 : inquire(file = inputFile, exist = exist, iostat = iostat, opened = opened, number = unit)
539 3 : exist = exist .and. iostat == 0
540 : ! Intel compiler throws error for `inputFile` longer than 4096 characters.
541 : !return_if_erred(__LINE__, iostat /= 0, SK_": Failed to inquire the properties of the user-specified input file: '"//trim(inputFile)//SK_"'. "//trim(adjustl(iomsg)))
542 :
543 : ! Read the contents of the input file.
544 :
545 3 : if (exist) then
546 : ! Close the input file and read its contents.
547 2 : if (opened) then
548 0 : close(unit, iostat = iostat, iomsg = iomsg)
549 0 : RETURN_IF_ERRED(__LINE__, iostat /= 0, SK_"Failed to close the user-specified input file: '"//trim(inputFile)//SK_"'. "//trim(adjustl(iomsg)))
550 : end if
551 2 : call setContentsFrom(inputFile, contents, iostat, iomsg)
552 : !open(newunit = unit, file = inputFile, form = "formatted", access = "sequential", status = "old", action = "read", iostat = iostat, iomsg = iomsg SHARED)
553 2 : RETURN_IF_ERRED(__LINE__, iostat /= 0, SK_"Failed to open the user-specified input file: '"//trim(inputFile)//SK_"'. "//trim(adjustl(iomsg)))
554 : else
555 1 : contents = inputFile
556 : end if
557 :
558 3 : iostat = index(getReplaced(contents, SK_" ", SK_""), SK_"&"//lcmethod, kind = IK)
559 3 : RETURN_IF_ERRED(__LINE__, iostat == 0, SK_"Failed to find a `&"//lcmethod//SK_" namelist group in user-specified input file: '"//trim(inputFile)//SK_"'. ")
560 :
561 : ! Now write the contents to a temporary file for the processor to parse it.
562 : ! The rewrite is important to avoid gfortran end of file IO errors.
563 :
564 3 : open(newunit = unit, form = "formatted", access = "sequential", status = "scratch", action = "readwrite", iostat = iostat, iomsg = iomsg)
565 3 : RETURN_IF_ERRED(__LINE__, iostat /= 0, SK_"Failed to open a temporary input file to hold user-specified simulation specifications: "//trim(adjustl(iomsg)))
566 3 : write(unit, "(A)", iostat = iostat, iomsg = iomsg) contents//new_line(contents) ! The new line is essential for successful namelist content parsing by the GNU compiler.
567 3 : RETURN_IF_ERRED(__LINE__, iostat /= 0, SK_"Failed to write the user-specified simulation specifications to a temporary input file: "//trim(adjustl(iomsg)))
568 3 : rewind(unit, iostat = iostat, iomsg = iomsg)
569 3 : RETURN_IF_ERRED(__LINE__, iostat /= 0, SK_"Failed to rewind the temporary input file containing the user-specified simulation specifications: "//trim(adjustl(iomsg)))
570 3 : lcmethod = getStrLower(method)
571 3 : if (lcmethod == SK_"paradise") then
572 0 : read(unit, nml = paradise, iostat = iostat, iomsg = iomsg)
573 3 : elseif (lcmethod == SK_"paradram") then
574 3 : read(unit, nml = paradram, iostat = iostat, iomsg = iomsg)
575 0 : elseif (lcmethod == SK_"paranest") then
576 0 : read(unit, nml = paranest, iostat = iostat, iomsg = iomsg)
577 : else
578 0 : error stop PROCEDURE_NAME//getFine(__FILE__, __LINE__)//SK_": Internal error occurred. Unrecognized input `method` = "//method
579 : end if
580 3 : close(unit, status = "delete", iostat = iostat)
581 :
582 : !if (iostat /= 0) then
583 : ! err%msg = err%msg//PROCEDURE_NAME//getFine(__FILE__, __LINE__)//SK_": Failed to read the sampler-specific namelist group `&"//method//SK_"`: "//trim(adjustl(iomsg))
584 : ! ! Try the generic namelist group name: ParaMonte
585 : ! read(INPUT_FILE, nml = paramonte, iostat = iostat, iomsg = iomsg)
586 : !end if
587 :
588 3 : if (iostat /= 0) then
589 0 : err%occurred = .true._LK
590 : ! Check if the any namelist group exists at all.
591 0 : if (index(getStrLower(contents), lcmethod, kind = IK) == 0_IK) then! .and. index(contents, SK_"paramonte", kind = IK) == 0_IK) then
592 0 : err%msg = err%msg//PROCEDURE_NAME//getFine(__FILE__, __LINE__)//SK_": Failed to detect any namelist groups named `&"//method//SK_"` in the specified input namelist file: "//trim(adjustl(iomsg))//SK_": """//contents//SK_""""
593 : else
594 0 : err%msg = err%msg//PROCEDURE_NAME//getFine(__FILE__, __LINE__)//SK_": Failed to read the simulation specifications from the specified input string/file: "//trim(adjustl(iomsg))//SK_": """//contents//SK_""""
595 : end if
596 : end if
597 :
598 6 : end subroutine
599 :
600 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
601 :
602 : #define ParaDISE_ENABLED 1
603 : #include "pm_sampling_scio.inc.F90"
604 : #undef ParaDISE_ENABLED
605 :
606 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
607 :
608 : #define ParaDRAM_ENABLED 1
609 : #include "pm_sampling_scio.inc.F90"
610 : #undef ParaDRAM_ENABLED
611 :
612 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
613 :
614 : #define ParaNest_ENABLED 1
615 : #include "pm_sampling_scio.inc.F90"
616 : #undef ParaNest_ENABLED
617 :
618 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
619 :
620 : ! check if is end of record or file error condition.
621 0 : pure function iserf(stat) result(ended)
622 : integer, intent(in) :: stat
623 : logical :: ended
624 0 : ended = is_iostat_eor(stat) .or. is_iostat_end(stat)
625 0 : end function
626 :
627 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
628 :
629 : #undef RETURN_IF_ERRED
630 : #undef SHARED
|