The ParaMonte Documentation Website
Current view: top level - kernel - SpecDRAM_AdaptiveUpdatePeriod_mod.f90 (source / functions) Hit Total Coverage
Test: ParaMonte 1.5.1 :: MPI Parallel Kernel - Code Coverage Report Lines: 19 19 100.0 %
Date: 2021-01-08 13:07:16 Functions: 4 4 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       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             : 
      43             : !> \brief
      44             : !> This module contains the classes and procedures for setting up the `adaptiveUpdatePeriod` attribute of samplers of class [ParaDRAM_type](@ref paradram_mod::paradram_type).
      45             : !> For more information, see the description of this attribute in the body of the module.
      46             : !> \author Amir Shahmoradi
      47             : 
      48             : module SpecDRAM_AdaptiveUpdatePeriod_mod
      49             : 
      50             :     use Constants_mod, only: IK
      51             :     implicit none
      52             : 
      53             :     character(*), parameter         :: MODULE_NAME = "@SpecDRAM_AdaptiveUpdatePeriod_mod"
      54             : 
      55             :     integer(IK)                     :: adaptiveUpdatePeriod ! namelist input
      56             : 
      57             :     type                            :: AdaptiveUpdatePeriod_type
      58             :         integer(IK)                 :: val
      59             :         integer(IK)                 :: def
      60             :         integer(IK)                 :: null
      61             :         character(:), allocatable   :: desc
      62             :     contains
      63             :         procedure, pass             :: set => setAdaptiveUpdatePeriod, checkForSanity, nullifyNameListVar
      64             :     end type AdaptiveUpdatePeriod_type
      65             : 
      66             :     interface AdaptiveUpdatePeriod_type
      67             :         module procedure            :: constructAdaptiveUpdatePeriod
      68             :     end interface AdaptiveUpdatePeriod_type
      69             : 
      70             :     private :: constructAdaptiveUpdatePeriod, setAdaptiveUpdatePeriod, checkForSanity, nullifyNameListVar
      71             : 
      72             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      73             : 
      74             : contains
      75             : 
      76             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      77             : 
      78        1047 :     function constructAdaptiveUpdatePeriod(ndim,methodName) result(AdaptiveUpdatePeriodObj)
      79             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
      80             :         !DEC$ ATTRIBUTES DLLEXPORT :: constructAdaptiveUpdatePeriod
      81             : #endif
      82             :         use Constants_mod, only: IK, NULL_IK
      83             :         use String_mod, only: num2str
      84             :         implicit none
      85             :         integer(IK), intent(in)         :: ndim
      86             :         character(*), intent(in)        :: methodName
      87             :         type(AdaptiveUpdatePeriod_type) :: AdaptiveUpdatePeriodObj
      88        1047 :         AdaptiveUpdatePeriodObj%def     = ndim * 4_IK !+ 1_IK ! max(ndim+1_IK,100_IK)
      89        1047 :         AdaptiveUpdatePeriodObj%null    = NULL_IK
      90             :         AdaptiveUpdatePeriodObj%desc    = &
      91             :         "Every adaptiveUpdatePeriod calls to the objective function, the parameters of the proposal distribution will be updated. &
      92             :         &The variable adaptiveUpdatePeriod must be a positive integer (>0). The smaller the value of adaptiveUpdatePeriod, &
      93             :         &the easier it will be for the " // methodName // " kernel to adapt the proposal distribution to the covariance structure &
      94             :         &of the objective function. However, this will happen at the expense of slower simulation runtime as the adaptation &
      95             :         &process can become computationally expensive, in particular, for very high dimensional objective functions (ndim>>1). &
      96             :         &The larger the value of adaptiveUpdatePeriod, the easier it will be &
      97             :         &for the " // methodName // " kernel to keep the sampling efficiency close to the requested target acceptance rate range &
      98             :         &(if specified via the input variable targetAcceptanceRate). &
      99             :         &However, too large values for adaptiveUpdatePeriod will only delay the adaptation of the proposal distribution to &
     100             :         &the global structure of the objective function that is being sampled. &
     101             :         &If adaptiveUpdatePeriod>=chainSize, then no adaptive updates to the proposal distribution will be made. &
     102             :         &The default value is 4 * ndim, where ndim is the dimension of the domain of the objective function to be sampled. &
     103             :         &In this particular " // methodName // " simulation, this corresponds to the value " // &
     104        1047 :         num2str(AdaptiveUpdatePeriodObj%def) // "."
     105        1047 :     end function constructAdaptiveUpdatePeriod
     106             : 
     107             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     108             : 
     109        1047 :     subroutine nullifyNameListVar(AdaptiveUpdatePeriodObj)
     110             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     111             :         !DEC$ ATTRIBUTES DLLEXPORT :: nullifyNameListVar
     112             : #endif
     113             :         implicit none
     114             :         class(AdaptiveUpdatePeriod_type), intent(in)  :: AdaptiveUpdatePeriodObj
     115        1047 :         adaptiveUpdatePeriod = AdaptiveUpdatePeriodObj%null
     116        1047 :     end subroutine nullifyNameListVar
     117             : 
     118             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     119             : 
     120        1107 :     subroutine setAdaptiveUpdatePeriod(AdaptiveUpdatePeriodObj,adaptiveUpdatePeriod)
     121             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     122             :         !DEC$ ATTRIBUTES DLLEXPORT :: setAdaptiveUpdatePeriod
     123             : #endif
     124        1047 :         use Constants_mod, only: IK
     125             :         implicit none
     126             :         class(AdaptiveUpdatePeriod_type), intent(inout) :: AdaptiveUpdatePeriodObj
     127             :         integer(IK), intent(in)                         :: adaptiveUpdatePeriod
     128        1107 :         AdaptiveUpdatePeriodObj%val = adaptiveUpdatePeriod
     129        1029 :         if ( AdaptiveUpdatePeriodObj%val==AdaptiveUpdatePeriodObj%null ) AdaptiveUpdatePeriodObj%val = AdaptiveUpdatePeriodObj%def
     130        1107 :     end subroutine setAdaptiveUpdatePeriod
     131             : 
     132             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     133             : 
     134        1035 :     subroutine checkForSanity(AdaptiveUpdatePeriodObj,Err,methodName)
     135             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     136             :         !DEC$ ATTRIBUTES DLLEXPORT :: checkForSanity
     137             : #endif
     138        1107 :         use Constants_mod, only: IK
     139             :         use Err_mod, only: Err_type
     140             :         use String_mod, only: num2str
     141             :         implicit none
     142             :         class(AdaptiveUpdatePeriod_type), intent(in)    :: AdaptiveUpdatePeriodObj
     143             :         character(*), intent(in)                        :: methodName
     144             :         type(Err_type), intent(inout)                   :: Err
     145             :         character(*), parameter                         :: PROCEDURE_NAME = "@checkForSanity()"
     146        1035 :         if ( AdaptiveUpdatePeriodObj%val<1_IK) then
     147           6 :             Err%occurred = .true.
     148             :             Err%msg =   Err%msg // &
     149             :                         MODULE_NAME // PROCEDURE_NAME // ": Error occurred. &
     150             :                         &Invalid requested value for adaptiveUpdatePeriod. &
     151             :                         &The input requested value for adaptiveUpdatePeriod (" // num2str(AdaptiveUpdatePeriodObj%val) // &
     152             :                         ") cannot be less than 1. If you are not sure of the appropriate value for adaptiveUpdatePeriod, drop it &
     153           6 :                         &from the input list. " // methodName // " will automatically assign an appropriate value to it.\n\n"
     154             :         end if
     155        2070 :     end subroutine checkForSanity
     156             : 
     157             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     158             : 
     159             : end module SpecDRAM_AdaptiveUpdatePeriod_mod ! LCOV_EXCL_LINE

ParaMonte: Plain Powerful Parallel Monte Carlo Library 
The Computational Data Science Lab
© Copyright 2012 - 2021