ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
Partition_mod.F90
Go to the documentation of this file.
1!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!!
4!!!! MIT License
5!!!!
6!!!! ParaMonte: plain powerful parallel Monte Carlo library.
7!!!!
8!!!! Copyright (C) 2012-present, The Computational Data Science Lab
9!!!!
10!!!! This file is part of the ParaMonte library.
11!!!!
12!!!! Permission is hereby granted, free of charge, to any person obtaining a
13!!!! copy of this software and associated documentation files (the "Software"),
14!!!! to deal in the Software without restriction, including without limitation
15!!!! the rights to use, copy, modify, merge, publish, distribute, sublicense,
16!!!! and/or sell copies of the Software, and to permit persons to whom the
17!!!! Software is furnished to do so, subject to the following conditions:
18!!!!
19!!!! The above copyright notice and this permission notice shall be
20!!!! included in all copies or substantial portions of the Software.
21!!!!
22!!!! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
23!!!! EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
24!!!! MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
25!!!! IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM,
26!!!! DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR
27!!!! OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE
28!!!! OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
29!!!!
30!!!! ACKNOWLEDGMENT
31!!!!
32!!!! ParaMonte is an honor-ware and its currency is acknowledgment and citations.
33!!!! As per the ParaMonte library license agreement terms, if you use any parts of
34!!!! this library for any purposes, kindly acknowledge the use of ParaMonte in your
35!!!! work (education/research/industry/development/...) by citing the ParaMonte
36!!!! library as described on this page:
37!!!!
38!!!! https://github.com/cdslaborg/paramonte/blob/main/ACKNOWLEDGMENT.md
39!!!!
40!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
41!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
42
48
49!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
50
52
53 use Except_mod, only: NEGBIG_RK
54 use Kind_mod, only: IK, LK, RK, SK
55 use Err_mod, only: Err_type
56
57 implicit none
58
59 private
60 public :: IK, LK, RK, SK
61 public :: MODULE_NAME
62 public :: PartitionMinVol_type
64
65 character(*, SK), parameter :: MODULE_NAME = "@Partition_mod"
66
70 type, abstract :: BaseProp_type
71 integer(IK) , allocatable :: Size(:)
72 integer(IK) , allocatable :: CumSumSize(:)
73 integer(IK) , allocatable :: Membership(:)
74 integer(IK) , allocatable :: EffectiveSize(:)
75 real(RK) , allocatable :: LogVolNormed(:)
76 real(RK) , allocatable :: ChoLowCovUpp(:,:,:)
77 real(RK) , allocatable :: InvCovMat(:,:,:)
78 real(RK) , allocatable :: ChoDia(:,:)
79 real(RK) , allocatable :: Center(:,:)
80 end type BaseProp_type
81
85 type, extends(BaseProp_type) :: ParProp_type
86 real(RK) , allocatable :: ScaleFactorSq(:)
87 real(RK) , allocatable :: ScaleFactor(:)
88 real(RK) , allocatable :: LogVolRatio(:)
89 real(RK) , allocatable :: MahalSq(:,:)
90 end type ParProp_type
91
95 type, extends(ParProp_type) :: ParPropTry_type
96 integer(IK) :: ic
97 integer(IK) :: nc
98 integer(IK) :: ne
99 integer(IK) :: np
100 integer(IK) :: nps
101 integer(IK) :: npe
102 integer(IK) :: itp
103 integer(IK) :: its
104 logical(LK) :: deallocated = .true._LK
105 logical(LK) , allocatable :: SubclusteringWarranted(:)
106 end type ParPropTry_type
107
112 real(RK) :: logSumVolNormed
113 logical(LK) :: subclusteringEnabled
114 !logical(LK) :: isBetterThanParent !< Logical flag that is `.true.` the sum of the volumes of the current-level
115 ! !! partitions is smaller than the shrinkage-corrected parent volume.
116 !real(RK) , allocatable :: LogVolNormedExpected(:) !< An array containing the expected normalized log-volumes of the partition at any given level.
117 end type ParPropTryVol_type
118
123 real(RK) , allocatable :: LogLikeFitness(:)
124 real(RK) , allocatable :: SumLogLikeFitness(:)
125 end type ParPropTryDen_type
126
130 type KmeansTry_type
131 real(RK) :: potential
132 real(RK) , allocatable :: Center(:,:)
133 integer(IK) , allocatable :: Membership(:)
134 integer(IK) , allocatable :: Size(:)
135 type(Err_type) :: Err
136 end type KmeansTry_type
137
140 type, abstract, extends(BaseProp_type) :: PartitionBase_type
141 integer(IK) :: ne
142 integer(IK) , allocatable :: Final(:,:)
143 real(RK) , allocatable :: VolNormedCDF(:)
144 end type PartitionBase_type
145
149 type, abstract, extends(PartitionBase_type) :: Partition_type
150 integer(IK) :: nd,np
151 integer(IK) :: nc,nt
152 integer(IK) :: nemax
153 integer(IK) :: minSize
154 integer(IK) :: kmeansNumFailMax
155 integer(IK) :: numRecursiveCall
156 integer(IK) :: kmeansNumRecursionMax
157 integer(IK) :: convergenceFailureCount
158 integer(IK) :: kvolumeNumRecursionMax
159 logical(LK) :: biasCorrectionEnabled
160 logical(LK) :: mahalSqWeightEnabled
161 logical(LK) :: stanEnabled
162 real(RK) :: nplog
163 real(RK) :: kmeansRelTol
164 real(RK) :: logExpansion
165 real(RK) :: logShrinkage
166 real(RK) :: pointLogVolNormed
167 real(RK) :: parentLogVolNormed
168 real(RK) :: mahalSqWeightExponent
169 real(RK) , allocatable :: BiasCorrectionScaleFactorSq(:)
170 integer(IK) , allocatable :: PointIndex(:)
171 type(Err_type) :: Err
172 contains
173 procedure(getPartition_proc), deferred :: get
174 procedure(write_proc) , deferred :: write
175 end type Partition_type
176
177 abstract interface
178 subroutine getPartition_proc(Partition, Point, parentLogVolNormed)
179 import :: Partition_type, RK
180 class(Partition_type) , intent(inout) :: Partition
181 real(RK) , intent(inout) , contiguous :: Point(:,:) ! (Partition%nd, Partition%np)
182 real(RK) , intent(in) , optional :: parentLogVolNormed
183 end subroutine getPartition_proc
184 end interface
185
186 interface
187 subroutine write_proc(Partition, fileUnit, Point)
188 import :: Partition_type, RK, IK
189 class(Partition_type) , intent(in) :: Partition
190 integer(IK) , intent(in) :: fileUnit
191 real(RK) , intent(in) :: Point(Partition%nd,Partition%np)
192 end subroutine write_proc
193 end interface
194
195!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
196!%%%% PartitionMaxDen specs
197!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
198
202 logical(LK) :: scaleOptimizationEnabled
203 logical(LK) :: shapeOptimizationEnabled
204 logical(LK) :: shapeAdaptationEnabled
205 logical(LK) :: optimizationEnabled
206 real(RK) , allocatable :: LogLikeFitness(:)
207 type(ParPropTryDen_type), allocatable :: Try(:)
208 contains
209 procedure, pass :: get => getPartitionMaxDen
210 procedure, pass :: write => writePartitionMaxDen
211 end type PartitionMaxDen_type
212
213 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
214
215 interface
216
254 module function constructPartitionMaxDen( Point & ! LCOV_EXCL_LINE
255 , nc,nt & ! LCOV_EXCL_LINE
256 , nemax & ! LCOV_EXCL_LINE
257 , minSize & ! LCOV_EXCL_LINE
258 , stanEnabled & ! LCOV_EXCL_LINE
259 , trimEnabled & ! LCOV_EXCL_LINE
260 , kmeansRelTol & ! LCOV_EXCL_LINE
261 , logExpansion & ! LCOV_EXCL_LINE
262 , logShrinkage & ! LCOV_EXCL_LINE
263 , parentLogVolNormed & ! LCOV_EXCL_LINE
264 , mahalSqWeightExponent & ! LCOV_EXCL_LINE
265 , kmeansNumFailMax & ! LCOV_EXCL_LINE
266 , kmeansNumRecursionMax & ! LCOV_EXCL_LINE
267 , kvolumeNumRecursionMax & ! LCOV_EXCL_LINE
268 , scaleOptimizationEnabled & ! LCOV_EXCL_LINE
269 , shapeOptimizationEnabled & ! LCOV_EXCL_LINE
270 , shapeAdaptationEnabled & ! LCOV_EXCL_LINE
271 , biasCorrectionEnabled & ! LCOV_EXCL_LINE
272 ) result(Partition)
273#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
274 !DEC$ ATTRIBUTES DLLEXPORT :: constructPartitionMaxDen
275#endif
276 implicit none
277 real(RK) , intent(inout) , contiguous :: Point(:,:)
278 real(RK) , intent(in) :: parentLogVolNormed
279 integer(IK) , intent(in) , optional :: nc,nt
280 integer(IK) , intent(in) , optional :: nemax
281 integer(IK) , intent(in) , optional :: minSize
282 integer(IK) , intent(in) , optional :: kmeansNumFailMax
283 integer(IK) , intent(in) , optional :: kmeansNumRecursionMax
284 integer(IK) , intent(in) , optional :: kvolumeNumRecursionMax
285 logical(LK) , intent(in) , optional :: scaleOptimizationEnabled
286 logical(LK) , intent(in) , optional :: shapeOptimizationEnabled
287 logical(LK) , intent(in) , optional :: shapeAdaptationEnabled
288 logical(LK) , intent(in) , optional :: biasCorrectionEnabled
289 logical(LK) , intent(in) , optional :: stanEnabled
290 logical(LK) , intent(in) , optional :: trimEnabled
291 real(RK) , intent(in) , optional :: kmeansRelTol
292 real(RK) , intent(in) , optional :: logExpansion
293 real(RK) , intent(in) , optional :: logShrinkage
294 real(RK) , intent(in) , optional :: mahalSqWeightExponent
295 type(PartitionMaxDen_type) :: Partition
296 end function constructPartitionMaxDen
297 end interface
298
299 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
300
301 interface
302 module subroutine getPartitionMaxDen(Partition, Point, parentLogVolNormed)
303#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
304 !DEC$ ATTRIBUTES DLLEXPORT :: getPartitionMaxDen
305#endif
306 class(PartitionMaxDen_type) , intent(inout) :: Partition
307 real(RK) , intent(inout) , contiguous :: Point(:,:) ! (Partition%nd, Partition%np)
308 real(RK) , intent(in) , optional :: parentLogVolNormed
309 end subroutine getPartitionMaxDen
310 end interface
311
312 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
313
314 interface
315 module subroutine writePartitionMaxDen ( Partition & ! LCOV_EXCL_LINE
316 , fileUnit & ! LCOV_EXCL_LINE
317 , Point & ! LCOV_EXCL_LINE
318 )
319#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
320 !DEC$ ATTRIBUTES DLLEXPORT :: writePartitionMaxDen
321#endif
322 class(PartitionMaxDen_type) , intent(in) :: Partition
323 integer(IK) , intent(in) :: fileUnit
324 real(RK) , intent(in) :: Point(Partition%nd,Partition%np)
325 character(*,SK) , parameter :: fileFormat = "(*(g0,:,','))"
326 integer(IK) :: i, j, ic
327 end subroutine writePartitionMaxDen
328 end interface
329
330 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
331
332 interface PartitionMaxDen_type
333 module procedure :: constructPartitionMaxDen
334 end interface PartitionMaxDen_type
335
336!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
337!%%%% PartitionMinVol specs
338!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339
343 logical(LK) :: isDynesty = .false._LK
344 logical(LK) :: isMultiNest = .false._LK
345 logical(LK) :: isParaMonte = .false._LK
346 character(:, SK), allocatable :: name
347 end type Method_type
348
352 logical(LK) :: inclusionFractionEnabled
353 real(RK) :: inclusionFraction
354 type(Method_type) :: Method
355 type(ParPropTryVol_type), allocatable :: Try(:)
356 contains
357 procedure, pass :: get => getPartitionMinVol
358 procedure, pass :: write => writePartitionMinVol
359 end type PartitionMinVol_type
360
361 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
362
363 interface
364
404 module function constructPartitionMinVol( Point & ! LCOV_EXCL_LINE
405 , nc,nt & ! LCOV_EXCL_LINE
406 , nemax & ! LCOV_EXCL_LINE
407 , minSize & ! LCOV_EXCL_LINE
408 , stanEnabled & ! LCOV_EXCL_LINE
409 , trimEnabled & ! LCOV_EXCL_LINE
410 , kmeansRelTol & ! LCOV_EXCL_LINE
411 , logExpansion & ! LCOV_EXCL_LINE
412 , logShrinkage & ! LCOV_EXCL_LINE
413 , inclusionFraction & ! LCOV_EXCL_LINE
414 , mahalSqWeightExponent & ! LCOV_EXCL_LINE
415 , parentLogVolNormed & ! LCOV_EXCL_LINE
416 , kmeansNumFailMax & ! LCOV_EXCL_LINE
417 , kmeansNumRecursionMax & ! LCOV_EXCL_LINE
418 , kvolumeNumRecursionMax & ! LCOV_EXCL_LINE
419 , scaleOptimizationEnabled & ! LCOV_EXCL_LINE
420 , shapeOptimizationEnabled & ! LCOV_EXCL_LINE
421 , shapeAdaptationEnabled & ! LCOV_EXCL_LINE
422 , biasCorrectionEnabled & ! LCOV_EXCL_LINE
423 , method & ! LCOV_EXCL_LINE
424 ) result(Partition)
425#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
426 !DEC$ ATTRIBUTES DLLEXPORT :: constructPartitionMinVol
427#endif
428 implicit none
429 real(RK) , intent(inout) , contiguous :: Point(:,:)
430 character(*, SK), intent(in) , optional :: method
431 integer(IK) , intent(in) , optional :: nc,nt
432 integer(IK) , intent(in) , optional :: nemax
433 integer(IK) , intent(in) , optional :: minSize
434 integer(IK) , intent(in) , optional :: kmeansNumFailMax
435 integer(IK) , intent(in) , optional :: kmeansNumRecursionMax
436 integer(IK) , intent(in) , optional :: kvolumeNumRecursionMax
437 logical(LK) , intent(in) , optional :: scaleOptimizationEnabled
438 logical(LK) , intent(in) , optional :: shapeOptimizationEnabled
439 logical(LK) , intent(in) , optional :: shapeAdaptationEnabled
440 logical(LK) , intent(in) , optional :: biasCorrectionEnabled
441 logical(LK) , intent(in) , optional :: stanEnabled
442 logical(LK) , intent(in) , optional :: trimEnabled
443 real(RK) , intent(in) , optional :: kmeansRelTol
444 real(RK) , intent(in) , optional :: logExpansion
445 real(RK) , intent(in) , optional :: logShrinkage
446 real(RK) , intent(in) , optional :: inclusionFraction
447 real(RK) , intent(in) , optional :: mahalSqWeightExponent
448 real(RK) , intent(in) , optional :: parentLogVolNormed
449 type(PartitionMinVol_type) :: Partition
450 end function constructPartitionMinVol
451 end interface
452
453 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
454
455 interface
456 module subroutine getPartitionMinVol(Partition, Point, parentLogVolNormed)
457#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
458 !DEC$ ATTRIBUTES DLLEXPORT :: getPartitionMinVol
459#endif
460 class(PartitionMinVol_type) , intent(inout) :: Partition
461 real(RK) , intent(inout) , contiguous :: Point(:,:) ! (Partition%nd, Partition%np)
462 real(RK) , intent(in) , optional :: parentLogVolNormed
463 end subroutine getPartitionMinVol
464 end interface
465
466 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
467
468 interface
469 module subroutine writePartitionMinVol ( Partition & ! LCOV_EXCL_LINE
470 , fileUnit & ! LCOV_EXCL_LINE
471 , Point & ! LCOV_EXCL_LINE
472 )
473#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
474 !DEC$ ATTRIBUTES DLLEXPORT :: writePartitionMinVol
475#endif
476 class(PartitionMinVol_type) , intent(in) :: Partition
477 integer(IK) , intent(in) :: fileUnit
478 real(RK) , intent(in) :: Point(Partition%nd,Partition%np)
479 character(*,SK) , parameter :: fileFormat = "(*(g0,:,','))"
480 integer(IK) :: i, j, ic
481 end subroutine writePartitionMinVol
482 end interface
483
484 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
485
486 interface PartitionMinVol_type
487 module procedure :: constructPartitionMinVol
488 end interface PartitionMinVol_type
489
490!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
491
492contains
493
494!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
495
498 function getBiasCorrectionScaleFactorSq(nd,npmin,npmax) result(BiasCorrectionScaleFactorSq)
499#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
500 !DEC$ ATTRIBUTES DLLEXPORT :: getBiasCorrectionScaleFactorSq
501#endif
502 use InterpLinear_mod, only: InterpLinear_type
503 use Kind_mod, only: IK, RK
504 implicit none
505
506 integer(IK) , intent(in) :: nd
507 integer(IK) , intent(in) :: npmin
508 integer(IK) , intent(in) :: npmax
509 real(RK) :: LogCount(npmin:npmax)
510 real(RK) :: BiasCorrectionScaleFactorSq(npmin:npmax)
511 integer(IK) :: ip, ndim !, npminDef, icase
512 type(InterpLinear_type) :: Interp
513
514 type :: Case_type
515 real(RK), allocatable :: LogNumPointInSpheroid(:)
516 real(RK), allocatable :: LogMeanLogScaleFactor(:)
517 end type Case_type
518#include "Partition_mod.inc.scaleFactor.inc.F90"
519
520 ndim = nd
521 if (npmin < nd + 1_IK) error stop
522 if (ndim < lbound(Case, dim = 1)) error stop
523 if (ndim > ubound(Case, dim = 1)) ndim = ubound(Case, dim = 1)
524
525 LogCount(npmin:npmax) = [( log(real(ip,RK)), ip = npmin, npmax )]
526 Interp = InterpLinear_type(SortedX = Case(ndim)%LogNumPointInSpheroid, SortedY = Case(ndim)%LogMeanLogScaleFactor)
527 BiasCorrectionScaleFactorSq(npmin:npmax) = exp(2_IK*exp(Interp%predict(QueryX = LogCount)))
528
529#if DEBUG_ENABLED || TESTING_ENABLED || CODECOV_ENABLED
530 ! WARNING: scaleFactorSq must be larger than 1 in all circumstances. This is crucial later on.
531 if (any(BiasCorrectionScaleFactorSq < 1._RK)) error stop
532#else
533 where (BiasCorrectionScaleFactorSq < 1._RK); BiasCorrectionScaleFactorSq = 1._RK; end where
534#endif
535
537
538!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
539
540 ! ATTN Developer: Do NOT lower the case of `PURE` below. This is a macro.
541 PURE function getLogLikeFitness(count, logVolRatio) result(logLikeFitness)
542#if INTEL_COMPILER_ENABLED && DLL_ENABLED && (WINDOWS_ENABLED || DARWIN_ENABLED)
543 !DEC$ ATTRIBUTES DLLEXPORT :: getLogLikeFitness
544#endif
545 use Except_mod, only: NEGBIG_RK
546 use Kind_mod, only: IK, RK
547 implicit none
548 integer(IK) , intent(in) :: count
549 real(RK) , intent(in) :: logVolRatio
550 real(RK) :: logLikeFitness
551 logLikeFitness = count * (1._RK + logVolRatio - exp(logVolRatio))
552#if CHECK_ENABLED
553 if (logLikeFitness < NEGBIG_RK .or. logLikeFitness > 0._RK) then
554 write(*,*) .or."Internal error occurred: logLikeFitness < NEGBIG_RK logLikeFitness > 0._RK"
555 write(*,*) "count, logVolRatio, logLikeFitness:", count, logVolRatio, logLikeFitness
556 error stop
557 end if
558#endif
559 end function getLogLikeFitness
560
561!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
562
563end module Partition_mod ! LCOV_EXCL_LINE
This module contains procedures and routines for ellipsoidal partitioning of a given set of points by...
PURE real(RK) function getLogLikeFitness(count, logVolRatio)
character(*, SK), parameter MODULE_NAME
real(RK) function, dimension(npmin:npmax) getBiasCorrectionScaleFactorSq(nd, npmin, npmax)
Abstract class containing the basic properties of partitions that are common between the exploration ...
The base class for generating a vector of Kmeans objects that hold the properties of different Kmeans...
The Method_type class for specifying which minimum-volume partitioning method to be used.
The subclass for generating objects that hold partition properties at any given level of partitioning...
The subclass for generating objects that hold partition properties at any given level of partitioning...
The subclass for generating objects that hold partition properties at any given level of partitioning...
The base class for generating objects that hold partition properties at any given level of partitioni...
The class for generating objects that hold the final output partition properties.
The PartitionMaxDen_type class.
The PartitionMinVol_type class.
The Partition_type abstract class for Partition derived types. Partitions an input array Point(nd,...