The ParaMonte Documentation Website
Current view: top level - kernel - SpecDRAM_DelayedRejectionScaleFactorVec_mod.f90 (source / functions) Hit Total Coverage
Test: ParaMonte 1.5.1 :: Coarray Parallel Kernel - Code Coverage Report Lines: 33 33 100.0 %
Date: 2021-01-08 12:59:07 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 `delayedRejectionScaleFactorVec` 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_DelayedRejectionScaleFactorVec_mod
      49             : 
      50             :     use Constants_mod, only: RK
      51             :     implicit none
      52             : 
      53             :     character(*), parameter         :: MODULE_NAME = "@SpecDRAM_DelayedRejectionScaleFactorVec_mod"
      54             : 
      55             :     real(RK), allocatable           :: DelayedRejectionScaleFactorVec(:) ! namelist input
      56             : 
      57             :     type                            :: DelayedRejectionScaleFactorVec_type
      58             :         real(RK), allocatable       :: Val(:)
      59             :         real(RK)                    :: def
      60             :         real(RK)                    :: null
      61             :         character(:), allocatable   :: desc
      62             :     contains
      63             :         procedure, pass             :: set, checkForSanity, nullifyNameListVar
      64             :     end type DelayedRejectionScaleFactorVec_type
      65             : 
      66             :     interface DelayedRejectionScaleFactorVec_type
      67             :         module procedure            :: construct
      68             :     end interface DelayedRejectionScaleFactorVec_type
      69             : 
      70             :     private :: construct, set, checkForSanity, nullifyNameListVar
      71             : 
      72             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      73             : 
      74             : contains
      75             : 
      76             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      77             : 
      78        1047 :     function construct(nd,methodName) result(self)
      79             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
      80             :         !DEC$ ATTRIBUTES DLLEXPORT :: construct
      81             : #endif
      82             :         use Decoration_mod, only: TAB
      83             :         use Constants_mod, only: IK, RK, NULL_RK
      84             :         use String_mod, only: num2str
      85             :         implicit none
      86             :         integer(IK), intent(in)                     :: nd
      87             :         character(*), intent(in)                    :: methodName
      88             :         type(DelayedRejectionScaleFactorVec_type)   :: self
      89        1047 :         self%def  = 0.5_RK**(1._RK/real(nd,kind=RK))  ! This gives a half volume to the covariance
      90        1047 :         self%null = NULL_RK
      91             :         self%desc = &
      92             :         "delayedRejectionScaleFactorVec is a real-valued positive vector of length (1:delayedRejectionCount) by which &
      93             :         &the covariance matrix of the proposal distribution of "//methodName//" sampler is scaled when the Delayed Rejection (DR) &
      94             :         &scheme is activated (by setting delayedRejectionCount>0). At each ith stage of the DR process, &
      95             :         &the proposal distribution from the last stage is scaled by the factor delayedRejectionScaleFactorVec(i). &
      96             :         &Missing elements of the delayedRejectionScaleFactorVec in the input to "//methodName//" will be set to the default value. &
      97             :         &The default value at all stages is 0.5^(1/ndim) = "//num2str(self%def)//", which reduces the &
      98             :         &volume of the covariance matrix of the proposal from the last DR stage by one half. The variable ndim represents the &
      99        1047 :         &number of dimensions of the Domain of the objective function."
     100        1047 :     end function construct
     101             : 
     102             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     103             : 
     104        1047 :     subroutine nullifyNameListVar(self)
     105             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     106             :         !DEC$ ATTRIBUTES DLLEXPORT :: nullifyNameListVar
     107             : #endif
     108        1047 :         use SpecDRAM_DelayedRejectionCount_mod, only: MAX_DELAYED_REJECTION_COUNT
     109             :         use Constants_mod, only: IK
     110             :         implicit none
     111             :         class(DelayedRejectionScaleFactorVec_type), intent(in) :: self
     112          12 :         if (allocated(DelayedRejectionScaleFactorVec)) deallocate(DelayedRejectionScaleFactorVec)
     113     1048050 :         allocate(DelayedRejectionScaleFactorVec(MAX_DELAYED_REJECTION_COUNT), source = self%null)
     114        1047 :     end subroutine nullifyNameListVar
     115             : 
     116             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     117             : 
     118        3153 :     pure subroutine set(self,delayedRejectionCount, DelayedRejectionScaleFactorVec)
     119             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     120             :         !DEC$ ATTRIBUTES DLLEXPORT :: set
     121             : #endif
     122        1047 :         use Constants_mod, only: IK, RK
     123             :         implicit none
     124             :         class(DelayedRejectionScaleFactorVec_type), intent(inout)   :: self
     125             :         integer(IK) , intent(in)                                    :: delayedRejectionCount
     126             :         real(RK)    , intent(in), optional                          :: DelayedRejectionScaleFactorVec(:) ! This must remain colon
     127        2088 :         if (present(DelayedRejectionScaleFactorVec)) then
     128     1048240 :             self%Val = pack( DelayedRejectionScaleFactorVec, mask = DelayedRejectionScaleFactorVec /= self%null )
     129        1065 :             if (size(self%Val) == 0) then
     130        1041 :                 deallocate(self%Val)
     131        2178 :                 allocate(self%Val(delayedRejectionCount), source = self%def)
     132             :             end if
     133             :         else
     134        1023 :             if (allocated(self%Val)) then
     135        1023 :                 if (size(self%Val) == 0) then
     136        1005 :                     deallocate(self%Val)
     137        7107 :                     allocate(self%Val(delayedRejectionCount), source = self%def)
     138             :                 end if
     139             :             end if
     140             :         end if
     141        2088 :     end subroutine set
     142             : 
     143             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     144             : 
     145        1035 :     subroutine checkForSanity(self,Err,methodName,delayedRejectionCount)
     146             : #if INTEL_COMPILER_ENABLED && defined DLL_ENABLED && (OS_IS_WINDOWS || defined OS_IS_DARWIN)
     147             :         !DEC$ ATTRIBUTES DLLEXPORT :: checkForSanity
     148             : #endif
     149        2088 :         use Constants_mod, only: IK, RK
     150             :         use Err_mod, only: Err_type
     151             :         use String_mod, only: num2str
     152             :         implicit none
     153             :         class(DelayedRejectionScaleFactorVec_type), intent(in)  :: self
     154             :         type(Err_type), intent(inout)                           :: Err
     155             :         character(*), intent(in)                                :: methodName
     156             :         integer(IK), intent(in)                                 :: delayedRejectionCount
     157             :         character(*), parameter                                 :: PROCEDURE_NAME = "@checkForSanity()"
     158             :         integer(IK)                                             :: i
     159        1035 :         if ( size(self%Val)/=delayedRejectionCount ) then
     160          12 :             Err%occurred = .true.
     161             :             Err%msg =   Err%msg // &
     162             :                         MODULE_NAME//PROCEDURE_NAME//": Error occurred. The length of the vector delayedRejectionScaleFactorVec (" // &
     163             :                         num2str(size(self%Val)) // ") is not equal to delayedRejectionCount = " // &
     164             :                         num2str(delayedRejectionCount) // ". If you are not sure how to set the values of &
     165             :                         &delayedRejectionScaleFactorVec, drop it from the input. " &
     166          12 :                         // methodName // " will automatically set the appropriate value for delayedRejectionScaleFactorVec.\n\n"
     167             :         end if
     168        7233 :         do i = 1,size(self%Val)
     169        7233 :             if ( self%Val(i)<=0._RK ) then
     170          12 :                 Err%occurred = .true.
     171             :                 Err%msg =   Err%msg // &
     172             :                             MODULE_NAME // PROCEDURE_NAME // ": Error occurred. The input value for the element " // &
     173          12 :                             num2str(i) // " of the variable delayedRejectionScaleFactorVec cannot be smaller than or equal to 0.\n\n"
     174             :             end if
     175             :         end do
     176        1035 :         deallocate(DelayedRejectionScaleFactorVec)
     177        2070 :     end subroutine checkForSanity
     178             : 
     179             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     180             : 
     181             : end module SpecDRAM_DelayedRejectionScaleFactorVec_mod ! LCOV_EXCL_LINE

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