Octopus
minimizer_tests.F90
Go to the documentation of this file.
1!! Copyright (C) 2024 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use debug_oct_m
24 use global_oct_m
25 use io_oct_m
26 use iso_c_binding
27 use, intrinsic :: iso_fortran_env
29 use parser_oct_m
33 use mpi_oct_m
35 use unit_oct_m
37
38 implicit none
39
40 private
41 public :: &
43
45 character(len=:), allocatable :: current_log
46
47contains
48
49 !---------------------------------------------------------------------------
53 subroutine test_optimizers(namespace)
54 type(namespace_t), intent(in) :: namespace
55
56 real(real64), dimension(2) :: x0, x, mass
57 real(real64) :: tol_grad = 1.0e-6_real64
58 real(real64) :: dt, energy
59 integer :: ierr, integrator
60 integer, parameter :: max_iter = 1000
61
62 push_sub(test_optimizers)
63
64 ! See src/main/geom_opt.F90 for the variable description
65 call parse_variable(namespace, 'GOFireIntegrator', option__gofireintegrator__verlet, integrator)
66
67 !=== Rosenbrock =========================================================
68 call set_log_file("rosenbrock.log")
69 call io_rm(trim(current_log), namespace)
70
71 ! Testing the FIRE algorithm with the Rosenbrock function in 2D
72 ! Starting point
73 x0 = [3.0_real64, 4.0_real64]
74
75 ! Initial parameters for the FIRE algorithm with Velocity Verlet
76 x = x0
77 ! This timestep is taken to give a reasonable result. Too large vaues lead to unstable results
78 ! This is set such that run this through valgrind+verrou gives reproducible results
79 dt = 0.001_real64
80 mass = m_one
81 call minimize_fire(2, 2, x, dt, tol_grad, max_iter, rosenbrock_gradient_2d, write_iter_info, energy, ierr, mass, integrator)
82
83 write(message(1), "(a,i6,a)") "FIRE (Rosenbrock) algorithm converged in ", -ierr, " iterations"
84 write(message(2), "(a,2(f8.4,a))") " The minimum is found to be at (", x(1), ", ", x(2), ")"
85 write(message(3), "(a)") " Analytical minimum is (1,1)."
86 call messages_info(3)
88
89 !=== Sphere =============================================================
90 call set_log_file("sphere.log")
91 call io_rm(trim(current_log), namespace)
92
93 x0 = [3.0_real64, -4.0_real64]
94 x = x0
95 dt = 0.01_real64
96
97 call minimize_fire(2, 2, x, dt, tol_grad, max_iter, &
98 sphere_gradient_2d, write_iter_info, energy, ierr, mass, integrator)
99
100 write(message(1), "(a,i6,a)") "FIRE (Sphere) converged in ", -ierr, " iterations"
101 write(message(2), "(a,2(f8.4,a))") " Minimum found at (", x(1), ", ", x(2), ")"
102 write(message(3), "(a)") " Analytical minimum is (0,0)."
103 call messages_info(3)
104 call messages_new_line()
105
106
107 !=== Rastrigin ==========================================================
108 call set_log_file("rastrigin.log")
109 call io_rm(trim(current_log), namespace)
110
111 x0 = [3.0_real64, 4.0_real64]
112 x = x0
113 dt = 0.001_real64
114 tol_grad = 1e-10_real64
115
116 call minimize_fire(2, 2, x, dt, tol_grad, max_iter, &
117 rastrigin_gradient_2d, write_iter_info, energy, ierr, mass, integrator)
118
119 write(message(1), "(a,i6,a)") "FIRE (Rastrigin) converged in ", -ierr, " iterations"
120 write(message(2), "(a,2(f8.4,a))") " Minimum found at (", x(1), ", ", x(2), ")"
121 write(message(3), "(a)") " Analytical minimum is (0,0)."
122 call messages_info(3)
123
124 deallocate(current_log)
125
126 pop_sub(test_optimizers)
127 end subroutine test_optimizers
128
129 !---------------------------------------------------------------------------
136 subroutine rosenbrock_gradient_2d(n, x, val, getgrad, grad) !< get the new gradients given a new x
137 integer, intent(in) :: n
138 real(real64), intent(in) :: x(n)
139 real(real64), intent(inout) :: val
140 integer, intent(in) :: getgrad
141 real(real64), intent(inout) :: grad(n)
142
143 real(real64), parameter :: a = m_one
144 real(real64), parameter :: b = 100.0_real64
145
146 push_sub(rosenbrock_gradient_2d)
147
148 ! Only the 2D version is implemented here
149 assert(n==2)
150
151 grad(1) = -m_two * (a-x(1)) - m_four * b * (x(2) -x(1)**2) * x(1)
152 grad(2) = m_two * b * (x(2) -x(1)**2)
153
155 end subroutine rosenbrock_gradient_2d
156
157
158 !---------------------------------------------------------------------------
160 subroutine sphere_gradient_2d(n, x, val, getgrad, grad)
161 integer, intent(in) :: n
162 real(real64), intent(in) :: x(n)
163 real(real64), intent(inout) :: val
164 integer, intent(in) :: getgrad
165 real(real64), intent(inout) :: grad(n)
166
167 push_sub(sphere_gradient_2d)
168 assert(n==2)
169
170 val = dot_product(x, x)
171 grad = m_two * x
172
173 pop_sub(sphere_gradient_2d)
174 end subroutine sphere_gradient_2d
175
176
177 !---------------------------------------------------------------------------
185 subroutine rastrigin_gradient_2d(n, x, val, getgrad, grad)
186 integer, intent(in) :: n
187 real(real64), intent(in) :: x(n)
188 real(real64), intent(inout) :: val
189 integer, intent(in) :: getgrad
190 real(real64), intent(inout) :: grad(n)
191
192 real(real64), parameter :: a = 10.0_real64
193
194 push_sub(rastrigin_gradient_2d)
195 assert(n==2)
196
197 val = m_two*a + &
198 (x(1)**2 - a*cos(m_two * m_pi * x(1))) + &
199 (x(2)**2 - a*cos(m_two * m_pi * x(2)))
200
201 grad = m_two*(x + m_pi * a * sin(m_two * m_pi * x))
202
203 pop_sub(rastrigin_gradient_2d)
204 end subroutine rastrigin_gradient_2d
205
206 !---------------------------------------------------------------------------
210 subroutine write_iter_info(iter, n, val, maxdr, maxgrad, x) !< output for each iteration step
211 integer, intent(in) :: iter
212 integer, intent(in) :: n
213 real(real64), intent(in) :: val
214 real(real64), intent(in) :: maxdr
215 real(real64), intent(in) :: maxgrad
216 real(real64), intent(in) :: x(n)
217
218 integer :: iunit
219
220 if (mpi_world%is_root()) then
221 iunit = io_open(current_log, global_namespace, action = 'write', position = 'append')
222 if (iter == 1) then
223 write(iunit, '(a10, 3a20)') '# iter','x','y', 'max_force'
224 end if
225 write(iunit, '(i10,3f20.10)') iter, x(1), x(2), maxgrad
226 call io_close(iunit)
227 end if
228
229 end subroutine write_iter_info
230
231 !---------------------------------------------------------------------------
233 subroutine set_log_file(name)
234 character(len=*), intent(in) :: name
235 if (allocated(current_log)) deallocate(current_log)
236 allocate(character(len=len_trim(name)) :: current_log)
237 current_log = trim(name)
238 end subroutine set_log_file
239
240end module minimizer_tests_oct_m
241
242!! Local Variables:
243!! mode: f90
244!! coding: utf-8
245!! End:
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:192
real(real64), parameter, public m_four
Definition: global.F90:194
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:188
real(real64), parameter, public m_one
Definition: global.F90:191
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:466
subroutine, public io_rm(fname, namespace)
Definition: io.F90:391
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:401
subroutine, public messages_new_line()
Definition: messages.F90:1118
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
This modules takes care of testing optimizers using standard test functions.
subroutine write_iter_info(iter, n, val, maxdr, maxgrad, x)
Helper function required by the minimizer.
subroutine rosenbrock_gradient_2d(n, x, val, getgrad, grad)
Gradient of the Rosenbrock function.
subroutine set_log_file(name)
Set the log‐file name before running each test.
subroutine rastrigin_gradient_2d(n, x, val, getgrad, grad)
Gradient of the Rastrigin function.
subroutine, public test_optimizers(namespace)
Unit tests for different optimizers.
subroutine sphere_gradient_2d(n, x, val, getgrad, grad)
Gradient of the sphere function.
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
type(namespace_t), public global_namespace
Definition: namespace.F90:134
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
This module defines the unit system, used for input and output.
for(l=0;l< nri;l++)
Definition: operate_inc.c:15