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 : #if TESTING_ENABLED && (CAF_ENABLED || MPI_ENABLED)
18 : #define FAILED_IMAGE(FAILED)isFailedImage(FAILED)
19 : #else
20 : #define FAILED_IMAGE(FAILED)FAILED
21 : #endif
22 : #if ParaDISE_ENABLED
23 : #define SET_ParaDISE(X)X
24 : #else
25 : #define SET_ParaDISE(X)
26 : #endif
27 : #if CAF_ENABLED || MPI_ENABLED
28 : #define SET_CAFMPI(X)X
29 : #define SET_SERIAL(X)
30 : #else
31 : #define SET_CAFMPI(X)
32 : #define SET_SERIAL(X)X
33 : #endif
34 : #if CAF_ENABLED
35 : #define SET_CAF(X)X
36 : #else
37 : #define SET_CAF(X)
38 : #endif
39 : #if MPI_ENABLED
40 : use mpi !mpi_f08, only: mpi_bcast, mpi_gather, mpi_allreduce, mpi_sum, mpi_integer, mpi_double_precision, mpi_comm_world
41 : #define SET_MPI(X)X
42 : #else
43 : #define SET_MPI(X)
44 : #endif
45 : #if OMP_ENABLED
46 : #define SET_OMP(X)X
47 : #else
48 : #define SET_OMP(X)
49 : #endif
50 : #if DEBUG_ENABLED
51 : #define SET_DEBUG(X)X
52 : #else
53 : #define SET_DEBUG(X)
54 : #endif
55 : ! \bug
56 : ! Avoid Intel ifort bug for too many `use` statements in a submodule by placing them all in the submodule header.
57 : #if CHECK_ENABLED
58 : use pm_err, only: setAsserted
59 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG)call setAsserted(ASSERTION,getFine(__FILE__,LINE)//MODULE_NAME//MSG);
60 : #else
61 : #define CHECK_ASSERTION(LINE,ASSERTION,MSG) continue;
62 : #endif
63 : ! Define error handling.
64 : #define RETURN_IF_FAILED(LINE,FAILED,MSG) \
65 : if (FAILED) then; \
66 : err%occurred = .true._LK; \
67 : err%msg = PROCEDURE_NAME//getFine(__FILE__, LINE)//SKC_": "//trim(MSG); \
68 : call spec%disp%stop%show(err%msg); \
69 : return; \
70 : end if;
71 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
72 : #if ParaDISE_ENABLED || ParaDRAM_ENABLED
73 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
74 :
75 : use pm_err, only: getFine
76 : use pm_err, only: err_type
77 : use pm_val2str, only: getStr
78 : use pm_timer, only: timer_type
79 : use pm_kind, only: SKC => SK, SK, IK, LK
80 : use pm_paramonte, only: PARAMONTE_WEB_ISSUES
81 : use pm_parallelism, only: isFailedImage
82 : use pm_sampling_proposal, only: NL1, NL2
83 : use pm_sampling_proposal, only: spec_type
84 : use pm_sampling_proposal, only: stat_type
85 : use pm_sampling_proposal, only: getErrChainRead
86 : use pm_sampling_proposal, only: getErrChainWrite
87 : use pm_sampling_proposal, only: chainFileColName
88 : use pm_sampling_proposal, only: isFailedChainResize
89 : use pm_sampling_proposal, only: proposal_type, setProposalAdapted, readRestart, writeRestart, setProposalStateNew
90 : #if CAF_ENABLED || MPI_ENABLED
91 : use pm_sampling_proposal, only: bcastProposalAdaptation
92 : #endif
93 : ! used by the kernel routines.
94 : use pm_distUnif, only: setUnifRand
95 : use pm_arrayResize, only: setResized
96 : use pm_sysPath, only: isFailedRemove
97 : use pm_mathLogSubExp, only: getLogSubExp
98 : use pm_io, only: setContentsFrom, setContentsTo
99 : use, intrinsic :: iso_c_binding, only: c_sizeof
100 : use pm_kind, only: RKD
101 :
102 : implicit none
103 : #if ParaDISE_ENABLED
104 : character(*,SKC), parameter :: MODULE_NAME = SK_"@pm_sampling_kernel_dise"
105 : #elif ParaDRAM_ENABLED
106 : character(*,SKC), parameter :: MODULE_NAME = SK_"@pm_sampling_kernel_dram"
107 : #else
108 : #error "Unrecognized interface."
109 : #endif
110 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
111 :
112 : abstract interface
113 : #if OMP_ENABLED && (MATLAB_ENABLED || PYTHON_ENABLED || R_ENABLED)
114 : function getLogFunc_proc(logFuncState, avgTimePerFunCall, avgCommPerFunCall) result(mold)
115 : import :: RKC, RKD
116 : real(RKC), intent(inout), contiguous :: logFuncState(0:, :)
117 : real(RKD), intent(inout) :: avgTimePerFunCall, avgCommPerFunCall
118 : real(RKC) :: mold
119 : end function
120 : #else
121 : function getLogFunc_proc(state) result(logFunc)
122 : import :: RKC
123 : real(RKC), intent(in), contiguous :: state(:)
124 : real(RKC) :: logFunc
125 : end function
126 : #endif
127 : end interface
128 :
129 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
130 :
131 : contains
132 :
133 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
134 :
135 : #if ParaDISE_ENABLED
136 0 : subroutine killMeAlreadyCMake_RK5(); use pm_sampling_proposal_dise_RK5, only: RKC; end subroutine
137 0 : subroutine killMeAlreadyCMake_RK4(); use pm_sampling_proposal_dise_RK4, only: RKC; end subroutine
138 0 : subroutine killMeAlreadyCMake_RK3(); use pm_sampling_proposal_dise_RK3, only: RKC; end subroutine
139 0 : subroutine killMeAlreadyCMake_RK2(); use pm_sampling_proposal_dise_RK2, only: RKC; end subroutine
140 0 : subroutine killMeAlreadyCMake_RK1(); use pm_sampling_proposal_dise_RK1, only: RKC; end subroutine
141 : #elif ParaDRAM_ENABLED
142 0 : subroutine killMeAlreadyCMake_RK5(); use pm_sampling_proposal_dram_RK5, only: RKC; end subroutine
143 0 : subroutine killMeAlreadyCMake_RK4(); use pm_sampling_proposal_dram_RK4, only: RKC; end subroutine
144 0 : subroutine killMeAlreadyCMake_RK3(); use pm_sampling_proposal_dram_RK3, only: RKC; end subroutine
145 0 : subroutine killMeAlreadyCMake_RK2(); use pm_sampling_proposal_dram_RK2, only: RKC; end subroutine
146 0 : subroutine killMeAlreadyCMake_RK1(); use pm_sampling_proposal_dram_RK1, only: RKC; end subroutine
147 : #elif ParaNest_ENABLED
148 : subroutine killMeAlreadyCMake_RK5(); use pm_sampling_proposal_nest_RK5, only: RKC; end subroutine
149 : subroutine killMeAlreadyCMake_RK4(); use pm_sampling_proposal_nest_RK4, only: RKC; end subroutine
150 : subroutine killMeAlreadyCMake_RK3(); use pm_sampling_proposal_nest_RK3, only: RKC; end subroutine
151 : subroutine killMeAlreadyCMake_RK2(); use pm_sampling_proposal_nest_RK2, only: RKC; end subroutine
152 : subroutine killMeAlreadyCMake_RK1(); use pm_sampling_proposal_nest_RK1, only: RKC; end subroutine
153 : #else
154 : #error "Unrecognized interface."
155 : #endif
156 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
157 :
158 : !> \brief
159 : !> Return an error if the sampler kernel run failed.
160 : !>
161 : !> \warning
162 : !> This procedure changes the components of the global-scope variables `spec` and `proposal`. Otherwise, it is pure.
163 13 : function getErrKernelRun(getLogFunc, spec, stat, proposal) result(err)
164 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
165 : !DEC$ ATTRIBUTES DLLEXPORT :: getErrKernelRun
166 : #endif
167 : procedure(getLogFunc_proc) :: getLogFunc
168 : type(proposal_type) , intent(inout) :: proposal
169 : type(spec_type) , intent(inout) :: spec
170 : type(stat_type) , intent(inout) :: stat
171 : type(err_type) :: err
172 : #include "pm_sampling_kernel.inc.F90"
173 26 : end function getErrKernelRun
174 :
175 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
176 :
177 : !> \brief
178 : !> Write the latest simulated state to the output chain file contents.
179 285063 : subroutine writeCFC(spec, stat, adaptationMeasure)
180 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
181 : !DEC$ ATTRIBUTES DLLEXPORT :: writeCFC
182 : #endif
183 : type(spec_type), intent(in) :: spec
184 : type(stat_type), intent(in) :: stat
185 : real(RKC), intent(in), contiguous :: adaptationMeasure(:)
186 : integer(IK) :: j
187 : ! if new point has been sampled, write the previous sampled point to output file.
188 285063 : blockOutputWrite: if (0_IK < stat%numFunCallAccepted) then
189 285050 : if (spec%outputChainFileFormat%isCompact) then
190 285050 : write(spec%chainFile%unit,spec%chainFile%format%rows) stat%cfc%processID (stat%numFunCallAccepted) &
191 285050 : , stat%cfc%delayedRejectionStage(stat%numFunCallAccepted) &
192 285050 : , stat%cfc%meanAcceptanceRate (stat%numFunCallAccepted) &
193 285050 : , stat%cfc%adaptationMeasure (stat%numFunCallAccepted) &
194 285050 : , stat%cfc%burninLocation (stat%numFunCallAccepted) &
195 285050 : , stat%cfc%sampleWeight (stat%numFunCallAccepted) &
196 285050 : , stat%cfc%sampleLogFunc (stat%numFunCallAccepted) &
197 570100 : , stat%cfc%sampleState (1 : spec%ndim%val, stat%numFunCallAccepted)
198 0 : elseif (spec%outputChainFileFormat%isBinary) then
199 0 : write(spec%chainFile%unit ) stat%cfc%processID (stat%numFunCallAccepted) &
200 0 : , stat%cfc%delayedRejectionStage(stat%numFunCallAccepted) &
201 0 : , stat%cfc%meanAcceptanceRate (stat%numFunCallAccepted) &
202 0 : , stat%cfc%adaptationMeasure (stat%numFunCallAccepted) &
203 0 : , stat%cfc%burninLocation (stat%numFunCallAccepted) &
204 0 : , stat%cfc%sampleWeight (stat%numFunCallAccepted) &
205 0 : , stat%cfc%sampleLogFunc (stat%numFunCallAccepted) &
206 0 : , stat%cfc%sampleState (1 : spec%ndim%val, stat%numFunCallAccepted)
207 0 : elseif (spec%outputChainFileFormat%isVerbose) then
208 0 : do j = 1, stat%cfc%sampleWeight(stat%numFunCallAccepted)
209 0 : write(spec%chainFile%unit,spec%chainFile%format%rows) stat%cfc%processID (stat%numFunCallAccepted) &
210 0 : , stat%cfc%delayedRejectionStage(stat%numFunCallAccepted) &
211 0 : , stat%cfc%meanAcceptanceRate (stat%numFunCallAccepted) &
212 : , adaptationMeasure(j) &
213 : !, stat%cfc%adaptationMeasure(stat%numFunCallAccepted) &
214 0 : , stat%cfc%burninLocation (stat%numFunCallAccepted) &
215 0 : , 1_IK &
216 0 : , stat%cfc%sampleLogFunc (stat%numFunCallAccepted) &
217 0 : , stat%cfc%sampleState (1 : spec%ndim%val, stat%numFunCallAccepted)
218 : end do
219 : end if
220 : end if blockOutputWrite
221 285063 : end subroutine writeCFC
222 :
223 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
224 :
225 : !> \brief
226 : !> Update the user with the simulation progress information.
227 : !>
228 : !> \remark
229 : !> Objects that change state in this subroutine are: `timer`, `stat%progress%timeElapsedSinceStartInSeconds`, `stat%numFunCallAcceptedLastReport`.
230 839 : subroutine reportProgress(spec, stat, timeLeft)!, sumAccrAccRej
231 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
232 : !DEC$ ATTRIBUTES DLLEXPORT :: reportProgress
233 : #endif
234 : use pm_strASCII, only: CR
235 : use pm_arraySplit, only: setSplit
236 : use iso_fortran_env, only: output_unit
237 : class(spec_type), intent(inout) :: spec
238 : class(stat_type), intent(inout) :: stat
239 : real(RKD), intent(in), optional :: timeLeft
240 : real(RKD) :: timeElapsedSinceLastReportInSeconds
241 : real(RKD) :: timeRemainedToFinishInSeconds
242 : real(RKC) :: meanAccRateSinceLastReport
243 : real(RKC) :: meanAccRateSinceStart
244 : character(1 + len(spec%reportFile%indent) + 3 * (25 + 3), SK) :: status
245 839 : if (spec%run%is%new) then
246 839 : stat%timer%clock = stat%timer%time()
247 839 : if (stat%progress%clock < 0._RKC) stat%progress%clock = stat%timer%start
248 839 : timeElapsedSinceLastReportInSeconds = stat%timer%clock - stat%progress%clock
249 839 : stat%progress%timeElapsedSinceStartInSeconds = stat%progress%timeElapsedSinceStartInSeconds + timeElapsedSinceLastReportInSeconds
250 839 : stat%progress%clock = stat%timer%clock
251 839 : if (present(timeLeft)) then
252 13 : timeRemainedToFinishInSeconds = timeLeft
253 : else
254 826 : timeRemainedToFinishInSeconds = stat%progress%timeElapsedSinceStartInSeconds * real(spec%outputChainSize%val - stat%numFunCallAccepted, RKD) / stat%numFunCallAccepted
255 : end if
256 2517 : CHECK_ASSERTION(__LINE__, stat%progress%counterPRP == spec%outputReportPeriod%val .or. spec%outputReportPeriod%val == 0_IK .or. timeRemainedToFinishInSeconds == 0._RKC, \
257 : SK_"The condition `counterPRP == outputReportPeriod .or. outputReportPeriod == 0 .or. timeLeft == 0` must hold. counterPRP, outputReportPeriod, timeLeft = "//\
258 : getStr([stat%progress%counterPRP, spec%outputReportPeriod%val])//SK_", "//getStr(timeRemainedToFinishInSeconds))
259 839 : meanAccRateSinceStart = real(stat%numFunCallAccepted, RKC) / real(stat%numFunCallAcceptedRejected, RKC)
260 839 : meanAccRateSinceLastReport = (stat%numFunCallAccepted - stat%numFunCallAcceptedLastReport) * spec%outputReportPeriod%inv
261 839 : stat%numFunCallAcceptedLastReport = stat%numFunCallAccepted
262 839 : write(spec%progressFile%unit, spec%progressFile%format%rows ) stat%numFunCallAcceptedRejected &
263 839 : , stat%numFunCallAccepted &
264 : , meanAccRateSinceStart &
265 839 : , meanAccRateSinceLastReport &
266 839 : , timeElapsedSinceLastReportInSeconds &
267 839 : , stat%progress%timeElapsedSinceStartInSeconds &
268 1678 : , timeRemainedToFinishInSeconds
269 839 : flush(spec%progressFile%unit)
270 : else
271 : block
272 : integer(IK) :: lenrec
273 : character(500, SK) :: record
274 : integer(IK), allocatable :: index(:,:)
275 : integer(IK) :: numFunCallAcceptedRejectedLastReport
276 0 : read(spec%progressFile%unit, "(A)") record
277 0 : lenrec = len_trim(record, IK)
278 0 : call setSplit(index, record(1:lenrec), spec%outputSeparator%val)
279 0 : read(record(index(1,1) : index(2,1)), *) numFunCallAcceptedRejectedLastReport
280 0 : read(record(index(1,2) : index(2,2)), *) stat%numFunCallAcceptedLastReport
281 0 : read(record(index(1,3) : index(2,3)), *) meanAccRateSinceStart
282 0 : read(record(index(1,4) : index(2,4)), *) meanAccRateSinceLastReport
283 0 : read(record(index(1,5) : index(2,5)), *) timeElapsedSinceLastReportInSeconds
284 0 : read(record(index(1,6) : index(2,6)), *) stat%progress%timeElapsedSinceStartInSeconds
285 0 : read(record(index(1,7) : index(2,7)), *) timeRemainedToFinishInSeconds ! estimatedTimeToFinishInSeconds
286 : end block
287 0 : stat%progress%clock = stat%timer%time()
288 : end if
289 : ! report progress in the standard output
290 839 : if (spec%image%is%first) then
291 : write(status, fmt = "(A1,A,*(A25,3X))" ) CR, spec%reportFile%indent, & ! LCOV_EXCL_LINE
292 : getStr(getStr(stat%numFunCallAccepted)//SK_" / "//getStr(stat%numFunCallAcceptedRejected), SK_"(125A)"), & ! LCOV_EXCL_LINE
293 : getStr(trim(adjustl(getStr(meanAccRateSinceLastReport, SK_"(1F11.4)")))//SK_" / "//trim(adjustl(getStr(stat%numFunCallAccepted / real(stat%numFunCallAcceptedRejected), SK_"(1F11.4)"))), SK_"(125A)"), & ! LCOV_EXCL_LINE
294 1678 : getStr(trim(adjustl(getStr(stat%progress%timeElapsedSinceStartInSeconds, SK_"(1F11.4)")))//SK_" / "//trim(adjustl(getStr(min(999999._RKD, timeRemainedToFinishInSeconds), SK_"(1F11.4)"))), SK_"(125A)")
295 839 : call spec%disp%show(status, unit = output_unit, advance = SK_"NO", tmsize = 0_IK, bmsize = 0_IK)
296 839 : if (timeRemainedToFinishInSeconds == 0._RKC) call spec%disp%skip(unit = output_unit)
297 : !call execute_command_line(" ")
298 839 : flush(output_unit)
299 : ! LCOV_EXCL_STOP
300 : end if
301 : !numFunCallAcceptedRejectedLastReport = stat%numFunCallAcceptedRejected
302 839 : stat%progress%counterPRP = 0_IK
303 839 : end subroutine reportProgress
304 :
305 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
306 :
307 13 : subroutine initProgressReport(spec)
308 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
309 : !DEC$ ATTRIBUTES DLLEXPORT :: initProgressReport
310 : #endif
311 : use iso_fortran_env, only: output_unit
312 : type(spec_type), intent(inout) :: spec
313 : character(:, SK), allocatable :: header
314 : integer(IK), parameter :: SEGLEN = 25
315 13 : if (spec%image%is%first) then
316 : ! LCOV_EXCL_START
317 : header = spec%reportFile%indent//"Accepted/Total Func. Call "//"Dynamic/Overall Acc. Rate "//"Elapsed/Remained Time [s] "//NL1//&
318 : repeat(" ",len(spec%reportFile%indent))//repeat("=",SEGLEN)//" "//repeat("=",SEGLEN)//" "//repeat("=",SEGLEN)//" "
319 : call spec%disp%show(header, unit = output_unit, format = SK_"(A)", tmsize = 0_IK, bmsize = 0_IK)
320 : !call execute_command_line(" ")
321 : deallocate(header)
322 : flush(output_unit)
323 : ! LCOV_EXCL_STOP
324 : end if
325 13 : end subroutine initProgressReport
326 :
327 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
328 :
329 : !> \brief
330 : !> Return the predicted burnin location within the currently sampled chain.
331 : !>
332 : !> \param[in] lenLogFunc : The length of the input `lenLogFunc`.
333 : !> \param[in] refLogFunc : The reference logFunc value with respect to which the relative importance of the points is gauged.
334 : !> \param[in] logFunc : The vector of length `lenLogFunc` of log(function) values.
335 : !>
336 : !> \devnote
337 : !> This is a `simple` function.
338 : !>
339 : !> \return
340 : !> `burninLoc` : The location of burnin in the input `logFunc` vector.
341 285050 : pure function getBurninLoc(lenLogFunc, refLogFunc, logFunc) result(burninLoc)
342 : #if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
343 : !DEC$ ATTRIBUTES DLLEXPORT :: getBurninLoc
344 : #endif
345 : integer(IK), intent(in) :: lenLogFunc
346 : real(RKC), intent(in) :: refLogFunc, logFunc(lenLogFunc)
347 : real(RKC) :: negLogIncidenceProb
348 : integer(IK) :: burninLoc
349 285050 : negLogIncidenceProb = log(real(lenLogFunc, RKC))
350 : burninLoc = 0_IK
351 : do
352 2619480 : burninLoc = burninLoc + 1_IK
353 2619480 : if (burninLoc < lenLogFunc .and. refLogFunc - logFunc(burninLoc) > negLogIncidenceProb) cycle
354 : !burninLoc = burninLoc - 1
355 2334430 : exit
356 : end do
357 285050 : end function getBurninLoc
358 :
359 : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
360 :
361 : #undef RETURN_IF_FAILED
362 : #undef SET_ParaDISE
363 :
364 : #else
365 : !%%%%%%%%%%%%%%%%%%%%%%%%
366 : #error "Unrecognized interface."
367 : !%%%%%%%%%%%%%%%%%%%%%%%%
368 : #endif
369 : #undef RETURN_IF_FAILED
370 : #undef SET_CAFMPI
371 : #undef SET_SERIAL
372 : #undef SET_DEBUG
373 : #undef SET_CAF
374 : #undef SET_MPI
375 : #undef SET_OMP
|