ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
test_pm_optimization.F90
Go to the documentation of this file.
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
19
21
23 use pm_test, only: test_type, LK
24 implicit none
25
26 private
27 public :: setTest
28 type(test_type) :: test
29
31 ! MATLAB reference code:
32 ! function f = objectivefcn1(x)
33 ! f = exp(-x^2) * cos(x) *sin(2*x);
34 ! end
35 ! CMD: x0 = [0.25];
36 ! CMD: options = optimset('tolx',1.e-10)
37 ! CMD: [x, f] = fminsearch(@objectivefcn1,x0,options)
38 real(RK) :: XMIN_REF = -0.471354350447655_RK
39 real(RK) :: FMIN_REF = -0.577293243101421_RK
40 contains
41 procedure, nopass :: get => getTestFuncRosenBrock1D
43
45 ! MATLAB reference code:
46 ! function f = objectivefcn1(x)
47 ! f = exp(-(x(1)-x(2))^2 - 2*x(1)^2)*cos(x(2))*sin(2*x(2));
48 ! end
49 ! CMD: x0 = [0.25,-0.25];
50 ! CMD: options = optimset('tolx',1.e-10)
51 ! CMD: [x, f] = fminsearch(@objectivefcn1,x0,options)
52 real(RK) :: XMIN_REF(2) = [-0.169552635123540_RK, -0.508657910611664_RK]
53 real(RK) :: FMIN_REF = -0.625285429865811_RK
54 contains
55 procedure, nopass :: get => getTestFuncRosenBrock2D
57
58!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
59
60contains
61
62!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
63
64 subroutine setTest()
65
67 call test%run(test_BrentMinimum_type_1, SK_"test_BrentMinimum_type_1")
68 call test%run(test_BrentMinimum_type_2, SK_"test_BrentMinimum_type_2")
69 call test%run(test_powellMin_type_1, SK_"test_powellMin_type_1") ! The internal function passing as actual argument causes segfault with Gfortran (any version) on Windows subsystem for Linux.
70 call test%summarize()
71
72 end subroutine setTest
73
74!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
75
76 function test_BrentMinimum_type_1() result(assertion)
77
78 use pm_kind, only: RK, IK
79 implicit none
80
81 logical(LK) :: assertion
82 type(brentMinimum_type) :: BrentMinimum
83 type(TestFuncRosenBrock1D_type) :: TestFuncRosenBrock1D
84
85 assertion = .true._LK
86
87 BrentMinimum = minimize( getFunc = getTestFuncRosenBrock1D, xtol = 1.e-8_RK )
88
89 assertion = .not. BrentMinimum%err%occurred
90 if (.not. assertion) then
91 ! LCOV_EXCL_START
92 if (test%traceable) then
93 write(test%disp%unit,"(*(g0,:,', '))")
94 write(test%disp%unit,"(*(g0,:,', '))") BrentMinimum%err%msg
95 write(test%disp%unit,"(*(g0,:,', '))")
96 end if
97 end if
98 ! LCOV_EXCL_STOP
99 call test%assert(assertion)
100
101 assertion = assertion .and. abs(BrentMinimum%xmin-TestFuncRosenBrock1D%XMIN_REF) / abs(TestFuncRosenBrock1D%XMIN_REF) < 1.e-6_RK
102 call test%assert(assertion)
103 assertion = assertion .and. abs(BrentMinimum%fmin-TestFuncRosenBrock1D%FMIN_REF) / abs(TestFuncRosenBrock1D%FMIN_REF) < 1.e-6_RK
104 call test%assert(assertion)
105
106 if (test%traceable .and. .not. assertion) then
107 ! LCOV_EXCL_START
108 write(test%disp%unit,"(*(g0,:,', '))")
109 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%xmin", BrentMinimum%xmin
110 write(test%disp%unit,"(*(g0,:,', '))") "XMIN_REF ", TestFuncRosenBrock1D%XMIN_REF
111 write(test%disp%unit,"(*(g0,:,', '))")
112 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%fmin", BrentMinimum%fmin
113 write(test%disp%unit,"(*(g0,:,', '))") "FMIN_REF ", TestFuncRosenBrock1D%FMIN_REF
114 write(test%disp%unit,"(*(g0,:,', '))")
115 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%niter", BrentMinimum%niter
116 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%Bracket", BrentMinimum%Bracket
117 write(test%disp%unit,"(*(g0,:,', '))")
118 end if
119 ! LCOV_EXCL_STOP
120
121 end function test_BrentMinimum_type_1
122
123!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
124
125 function test_BrentMinimum_type_2() result(assertion)
126
127 use pm_kind, only: RK, IK
128 implicit none
129
130 logical(LK) :: assertion
131 type(brentMinimum_type) :: BrentMinimum
132 type(TestFuncRosenBrock1D_type) :: TestFuncRosenBrock1D
133
134 assertion = .true._LK
135
136 BrentMinimum = brentMinimum_type( getFunc = getTestFuncRosenBrock1D &
137 , xtol = 1.e-8_RK &
138 , x0 = -1._RK &
139 , x1 = 0._RK &
140 , x2 = .5_RK &
141 )
142
143 assertion = .not. BrentMinimum%err%occurred
144 if (.not. assertion) then
145 ! LCOV_EXCL_START
146 if (test%traceable) then
147 write(test%disp%unit,"(*(g0,:,', '))")
148 write(test%disp%unit,"(*(g0,:,', '))") BrentMinimum%err%msg
149 write(test%disp%unit,"(*(g0,:,', '))")
150 end if
151 ! LCOV_EXCL_STOP
152 end if
153 call test%assert(assertion)
154
155 assertion = abs(BrentMinimum%xmin-TestFuncRosenBrock1D%XMIN_REF) / abs(TestFuncRosenBrock1D%XMIN_REF) < 1.e-6_RK
156 call test%assert(assertion)
157 assertion = abs(BrentMinimum%fmin-TestFuncRosenBrock1D%FMIN_REF) / abs(TestFuncRosenBrock1D%FMIN_REF) < 1.e-6_RK
158 call test%assert(assertion)
159
160 if (test%traceable .and. .not. assertion) then
161 ! LCOV_EXCL_START
162 write(test%disp%unit,"(*(g0,:,', '))")
163 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%xmin", BrentMinimum%xmin
164 write(test%disp%unit,"(*(g0,:,', '))") "XMIN_REF ", TestFuncRosenBrock1D%XMIN_REF
165 write(test%disp%unit,"(*(g0,:,', '))")
166 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%fmin", BrentMinimum%fmin
167 write(test%disp%unit,"(*(g0,:,', '))") "FMIN_REF ", TestFuncRosenBrock1D%FMIN_REF
168 write(test%disp%unit,"(*(g0,:,', '))")
169 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%niter", BrentMinimum%niter
170 write(test%disp%unit,"(*(g0,:,', '))") "BrentMinimum%Bracket", BrentMinimum%Bracket
171 write(test%disp%unit,"(*(g0,:,', '))")
172 end if
173 ! LCOV_EXCL_STOP
174
175 end function test_BrentMinimum_type_2
176
177!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
178
179 function test_powellMin_type_1() result(assertion)
180
181 use pm_kind, only: RK, IK
182 implicit none
183
184 logical(LK) :: assertion
185 type(powellMin_type) :: powellMin
186 type(TestFuncRosenBrock2D_type) :: TestFuncRosenBrock2D
187 real(RK), allocatable :: StartVec(:)
188
189 assertion = .true._LK
190
191 StartVec = [.25_RK,-.25_RK]
192
193 powellMin = powellMin_type( ndim = 2_IK &
194 , getFuncMD = getTestFuncRosenBrock2D &
195 , StartVec = StartVec & ! LCOV_EXCL_LINE
196 !, DirMat = reshape([1._RK, 0._RK, 0._RK, 1._RK], shape = [2,2])
197 !, ftol = 1.e-8_RK &
198 )
199
200 assertion = .not. powellMin%err%occurred
201 if (.not. assertion) then
202 ! LCOV_EXCL_START
203 if (test%traceable) then
204 write(test%disp%unit,"(*(g0,:,', '))")
205 write(test%disp%unit,"(*(g0,:,', '))") powellMin%err%msg
206 write(test%disp%unit,"(*(g0,:,', '))")
207 end if
208 end if
209 ! LCOV_EXCL_STOP
210 call test%assert(assertion)
211
212 assertion = assertion .and. any(abs(powellMin%xmin-TestFuncRosenBrock2D%XMIN_REF) / abs(TestFuncRosenBrock2D%XMIN_REF) < 1.e-6_RK)
213 call test%assert(assertion)
214 assertion = assertion .and. abs(powellMin%fmin-TestFuncRosenBrock2D%FMIN_REF) / abs(TestFuncRosenBrock2D%FMIN_REF) < 1.e-6_RK
215 call test%assert(assertion)
216
217 if (test%traceable .and. .not. assertion) then
218 ! LCOV_EXCL_START
219 write(test%disp%unit,"(*(g0,:,', '))")
220 write(test%disp%unit,"(*(g0,:,', '))") "powellMin%xmin", powellMin%xmin
221 write(test%disp%unit,"(*(g0,:,', '))") "XMIN_REF ", TestFuncRosenBrock2D%XMIN_REF
222 write(test%disp%unit,"(*(g0,:,', '))")
223 write(test%disp%unit,"(*(g0,:,', '))") "powellMin%fmin", powellMin%fmin
224 write(test%disp%unit,"(*(g0,:,', '))") "FMIN_REF ", TestFuncRosenBrock2D%FMIN_REF
225 write(test%disp%unit,"(*(g0,:,', '))")
226 write(test%disp%unit,"(*(g0,:,', '))") "powellMin%niter", powellMin%niter
227 write(test%disp%unit,"(*(g0,:,', '))") "powellMin%ftol", powellMin%ftol
228 write(test%disp%unit,"(*(g0,:,', '))")
229 write(test%disp%unit,"(*(g0,:,', '))") "powellMin%DirMat", powellMin%DirMat
230 write(test%disp%unit,"(*(g0,:,', '))")
231 end if
232 ! LCOV_EXCL_STOP
233
234 end function test_powellMin_type_1
235
236!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
237
238 pure function getTestFuncRosenBrock1D(x) result(testFuncVal)
239 use pm_kind, only: RK
240 implicit none
241 real(RK), intent(in) :: x
242 real(RK) :: testFuncVal
243 testFuncVal = exp(-x**2) * cos(x) * sin(2*x)
244 end function getTestFuncRosenBrock1D
245
246!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
247
248 pure function getTestFuncRosenBrock2D(ndim,Point) result(testFuncVal)
249 use pm_kind, only: RK, IK
250 implicit none
251 integer(IK) , intent(in) :: ndim
252 real(RK) , intent(in) :: Point(ndim)
253 real(RK) :: testFuncVal
254 testFuncVal = exp(-(Point(1)-Point(2))**2 - 2*Point(1)**2) * cos(Point(2)) * sin(2*Point(2))
255 end function getTestFuncRosenBrock2D
256
257!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
258
259end module test_pm_optimization ! LCOV_EXCL_LINE
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
This module contains procedures, generic interfaces, and types for numerical optimizations of mathema...
character(*, SK), parameter MODULE_NAME
This module contains a simple unit-testing framework for the Fortran libraries, including the ParaMon...
Definition: pm_test.F90:42
This module contains tests of the module pm_optimization.
logical(LK) function test_powellMin_type_1()
pure real(RK) function getTestFuncRosenBrock1D(x)
pure real(RK) function getTestFuncRosenBrock2D(ndim, Point)
logical(LK) function test_BrentMinimum_type_1()
logical(LK) function test_BrentMinimum_type_2()
This is the derived type test_type for generating objects that facilitate testing of a series of proc...
Definition: pm_test.F90:209