The ParaMonte Documentation Website
Current view: top level - kernel/tests - Test_Optimization_mod.f90 (source / functions) Hit Total Coverage
Test: ParaMonte 1.5.1 :: MPI Parallel Kernel - Code Coverage Report Lines: 49 49 100.0 %
Date: 2021-01-08 13:07:16 Functions: 6 6 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 This module contains tests of the module [Optimization_mod](@ref optimization_mod).
      44             : !>  \author Amir Shahmoradi
      45             : 
      46             : module Test_Optimization_mod
      47             : 
      48             :     use Optimization_mod
      49             :     use Test_mod, only: Test_type
      50             :     implicit none
      51             : 
      52             :     private
      53             :     public :: test_Optimization
      54             : 
      55             :     type(Test_type) :: Test
      56             : 
      57             :     type :: TestFuncRosenBrock1D_type
      58             :         ! MATLAB reference code:
      59             :         ! function f = objectivefcn1(x)
      60             :         !     f = exp(-x^2) * cos(x) *sin(2*x);
      61             :         ! end
      62             :         ! CMD: x0 = [0.25];
      63             :         ! CMD: options = optimset('tolx',1.e-10)
      64             :         ! CMD: [x, f] = fminsearch(@objectivefcn1,x0,options)
      65             :         real(RK) :: XMIN_REF = -0.471354350447655_RK
      66             :         real(RK) :: FMIN_REF = -0.577293243101421_RK
      67             :     contains
      68             :         procedure, nopass :: get => getTestFuncRosenBrock1D
      69             :     end type TestFuncRosenBrock1D_type
      70             : 
      71             :     type :: TestFuncRosenBrock2D_type
      72             :         ! MATLAB reference code:
      73             :         ! function f = objectivefcn1(x)
      74             :         !     f = exp(-(x(1)-x(2))^2 - 2*x(1)^2)*cos(x(2))*sin(2*x(2));
      75             :         ! end
      76             :         ! CMD: x0 = [0.25,-0.25];
      77             :         ! CMD: options = optimset('tolx',1.e-10)
      78             :         ! CMD: [x, f] = fminsearch(@objectivefcn1,x0,options)
      79             :         real(RK) :: XMIN_REF(2) = [-0.169552635123540_RK, -0.508657910611664_RK]
      80             :         real(RK) :: FMIN_REF = -0.625285429865811_RK
      81             :     contains
      82             :         procedure, nopass :: get => getTestFuncRosenBrock2D
      83             :     end type TestFuncRosenBrock2D_type
      84             : 
      85             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      86             : 
      87             : contains
      88             : 
      89             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
      90             : 
      91           3 :     subroutine test_Optimization()
      92             :         implicit none
      93           3 :         Test = Test_type(moduleName=MODULE_NAME)
      94           3 :         call Test%run(test_BrentMinimum_type_1, "test_BrentMinimum_type_1")
      95           3 :         call Test%run(test_BrentMinimum_type_2, "test_BrentMinimum_type_2")
      96           3 :         call Test%run(test_PowellMinimum_type_1, "test_PowellMinimum_type_1") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
      97           3 :         call Test%finalize()
      98           3 :     end subroutine test_Optimization
      99             : 
     100             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     101             : 
     102           3 :     function test_BrentMinimum_type_1() result(assertion)
     103             : 
     104           3 :         use Constants_mod, only: RK, IK
     105             :         implicit none
     106             : 
     107             :         logical                         :: assertion
     108           3 :         type(BrentMinimum_type)         :: BrentMinimum
     109             :         type(TestFuncRosenBrock1D_type) :: TestFuncRosenBrock1D
     110             : 
     111           3 :         assertion = .true.
     112             : 
     113           3 :         BrentMinimum = minimize( getFunc = getTestFuncRosenBrock1D, xtol = 1.e-8_RK )
     114             : 
     115           3 :         assertion = .not. BrentMinimum%Err%occurred
     116           3 :         if (.not. assertion) then
     117             :         ! LCOV_EXCL_START
     118             :             if (Test%isDebugMode) then
     119             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     120             :                 write(Test%outputUnit,"(*(g0,:,', '))") BrentMinimum%Err%msg
     121             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     122             :             end if
     123             :             return
     124             :         end if
     125             :         ! LCOV_EXCL_STOP
     126             : 
     127           3 :         assertion = assertion .and. abs(BrentMinimum%xmin-TestFuncRosenBrock1D%XMIN_REF) / abs(TestFuncRosenBrock1D%XMIN_REF) < 1.e-6_RK
     128           3 :         assertion = assertion .and. abs(BrentMinimum%fmin-TestFuncRosenBrock1D%FMIN_REF) / abs(TestFuncRosenBrock1D%FMIN_REF) < 1.e-6_RK
     129             : 
     130           3 :         if (Test%isDebugMode .and. .not. assertion) then
     131             :         ! LCOV_EXCL_START
     132             :             write(Test%outputUnit,"(*(g0,:,', '))")
     133             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%xmin", BrentMinimum%xmin
     134             :             write(Test%outputUnit,"(*(g0,:,', '))") "XMIN_REF         ", TestFuncRosenBrock1D%XMIN_REF
     135             :             write(Test%outputUnit,"(*(g0,:,', '))")
     136             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%fmin", BrentMinimum%fmin
     137             :             write(Test%outputUnit,"(*(g0,:,', '))") "FMIN_REF         ", TestFuncRosenBrock1D%FMIN_REF
     138             :             write(Test%outputUnit,"(*(g0,:,', '))")
     139             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%niter", BrentMinimum%niter
     140             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%Bracket", BrentMinimum%Bracket
     141             :             write(Test%outputUnit,"(*(g0,:,', '))")
     142             :         end if
     143             :         ! LCOV_EXCL_STOP
     144             : 
     145           3 :     end function test_BrentMinimum_type_1
     146             : 
     147             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     148             : 
     149           3 :     function test_BrentMinimum_type_2() result(assertion)
     150             : 
     151           3 :         use Constants_mod, only: RK, IK
     152             :         implicit none
     153             : 
     154             :         logical                         :: assertion
     155           3 :         type(BrentMinimum_type)         :: BrentMinimum
     156             :         type(TestFuncRosenBrock1D_type) :: TestFuncRosenBrock1D
     157             : 
     158           3 :         assertion = .true.
     159             : 
     160             :         BrentMinimum = BrentMinimum_type( getFunc = getTestFuncRosenBrock1D &
     161             :                                         , xtol = 1.e-8_RK &
     162             :                                         , x0 = -1._RK &
     163             :                                         , x1 = 0._RK &
     164             :                                         , x2 = .5_RK &
     165           3 :                                         )
     166             : 
     167           3 :         assertion = .not. BrentMinimum%Err%occurred
     168           3 :         if (.not. assertion) then
     169             :             ! LCOV_EXCL_START
     170             :             if (Test%isDebugMode) then
     171             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     172             :                 write(Test%outputUnit,"(*(g0,:,', '))") BrentMinimum%Err%msg
     173             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     174             :             end if
     175             :             return
     176             :             ! LCOV_EXCL_STOP
     177             :         end if
     178             : 
     179           3 :         assertion = abs(BrentMinimum%xmin-TestFuncRosenBrock1D%XMIN_REF) / abs(TestFuncRosenBrock1D%XMIN_REF) < 1.e-6_RK
     180           3 :         assertion = abs(BrentMinimum%fmin-TestFuncRosenBrock1D%FMIN_REF) / abs(TestFuncRosenBrock1D%FMIN_REF) < 1.e-6_RK
     181             : 
     182           3 :         if (Test%isDebugMode .and. .not. assertion) then
     183             :         ! LCOV_EXCL_START
     184             :             write(Test%outputUnit,"(*(g0,:,', '))")
     185             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%xmin", BrentMinimum%xmin
     186             :             write(Test%outputUnit,"(*(g0,:,', '))") "XMIN_REF         ", TestFuncRosenBrock1D%XMIN_REF
     187             :             write(Test%outputUnit,"(*(g0,:,', '))")
     188             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%fmin", BrentMinimum%fmin
     189             :             write(Test%outputUnit,"(*(g0,:,', '))") "FMIN_REF         ", TestFuncRosenBrock1D%FMIN_REF
     190             :             write(Test%outputUnit,"(*(g0,:,', '))")
     191             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%niter", BrentMinimum%niter
     192             :             write(Test%outputUnit,"(*(g0,:,', '))") "BrentMinimum%Bracket", BrentMinimum%Bracket
     193             :             write(Test%outputUnit,"(*(g0,:,', '))")
     194             :         end if
     195             :         ! LCOV_EXCL_STOP
     196             : 
     197           3 :     end function test_BrentMinimum_type_2
     198             : 
     199             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     200             : 
     201           3 :     function test_PowellMinimum_type_1() result(assertion)
     202             : 
     203           3 :         use Constants_mod, only: RK, IK
     204             :         implicit none
     205             : 
     206             :         logical                         :: assertion
     207           3 :         type(PowellMinimum_type)        :: PowellMinimum
     208             :         type(TestFuncRosenBrock2D_type) :: TestFuncRosenBrock2D
     209             :         real(RK), allocatable           :: StartVec(:)
     210             : 
     211           3 :         assertion = .true.
     212             : 
     213           9 :         StartVec = [.25_RK,-.25_RK]
     214             : 
     215             :         PowellMinimum = PowellMinimum_type  ( ndim = 2_IK &
     216             :                                             , getFuncMD = getTestFuncRosenBrock2D &
     217             :                                             , StartVec = StartVec & ! LCOV_EXCL_LINE
     218             :                                             !, DirMat = reshape([1._RK, 0._RK, 0._RK, 1._RK], shape = [2,2])
     219             :                                             !, ftol = 1.e-8_RK &
     220           3 :                                             )
     221             : 
     222           3 :         assertion = .not. PowellMinimum%Err%occurred
     223           3 :         if (.not. assertion) then
     224             :         ! LCOV_EXCL_START
     225             :             if (Test%isDebugMode) then
     226             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     227             :                 write(Test%outputUnit,"(*(g0,:,', '))") PowellMinimum%Err%msg
     228             :                 write(Test%outputUnit,"(*(g0,:,', '))")
     229             :             end if
     230             :             return
     231             :         end if
     232             :         ! LCOV_EXCL_STOP
     233             : 
     234           3 :         assertion = assertion .and. any(abs(PowellMinimum%xmin-TestFuncRosenBrock2D%XMIN_REF) / abs(TestFuncRosenBrock2D%XMIN_REF) < 1.e-6_RK)
     235           3 :         assertion = assertion .and. abs(PowellMinimum%fmin-TestFuncRosenBrock2D%FMIN_REF) / abs(TestFuncRosenBrock2D%FMIN_REF) < 1.e-6_RK
     236             : 
     237           3 :         if (Test%isDebugMode .and. .not. assertion) then
     238             :         ! LCOV_EXCL_START
     239             :             write(Test%outputUnit,"(*(g0,:,', '))")
     240             :             write(Test%outputUnit,"(*(g0,:,', '))") "PowellMinimum%xmin", PowellMinimum%xmin
     241             :             write(Test%outputUnit,"(*(g0,:,', '))") "XMIN_REF          ", TestFuncRosenBrock2D%XMIN_REF
     242             :             write(Test%outputUnit,"(*(g0,:,', '))")
     243             :             write(Test%outputUnit,"(*(g0,:,', '))") "PowellMinimum%fmin", PowellMinimum%fmin
     244             :             write(Test%outputUnit,"(*(g0,:,', '))") "FMIN_REF          ", TestFuncRosenBrock2D%FMIN_REF
     245             :             write(Test%outputUnit,"(*(g0,:,', '))")
     246             :             write(Test%outputUnit,"(*(g0,:,', '))") "PowellMinimum%niter", PowellMinimum%niter
     247             :             write(Test%outputUnit,"(*(g0,:,', '))") "PowellMinimum%ftol", PowellMinimum%ftol
     248             :             write(Test%outputUnit,"(*(g0,:,', '))")
     249             :             write(Test%outputUnit,"(*(g0,:,', '))") "PowellMinimum%DirMat", PowellMinimum%DirMat
     250             :             write(Test%outputUnit,"(*(g0,:,', '))")
     251             :         end if
     252             :         ! LCOV_EXCL_STOP
     253             : 
     254          12 :     end function test_PowellMinimum_type_1
     255             : 
     256             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     257             : 
     258          81 :     pure function getTestFuncRosenBrock1D(x) result(testFuncVal)
     259           3 :         use Constants_mod, only: RK
     260             :         implicit none
     261             :         real(RK), intent(in)    :: x
     262             :         real(RK)                :: testFuncVal
     263          81 :         testFuncVal = exp(-x**2) * cos(x) * sin(2*x)
     264         162 :     end function getTestFuncRosenBrock1D
     265             : 
     266             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     267             : 
     268         519 :     pure function getTestFuncRosenBrock2D(ndim,Point) result(testFuncVal)
     269          81 :         use Constants_mod, only: RK, IK
     270             :         implicit none
     271             :         integer(IK) , intent(in)    :: ndim
     272             :         real(RK)    , intent(in)    :: Point(ndim)
     273             :         real(RK)                    :: testFuncVal
     274         519 :         testFuncVal = exp(-(Point(1)-Point(2))**2 - 2*Point(1)**2) * cos(Point(2)) * sin(2*Point(2))
     275         519 :     end function getTestFuncRosenBrock2D
     276             : 
     277             : !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
     278             : 
     279             : end module Test_Optimization_mod ! LCOV_EXCL_LINE

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