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
32 use mpi_oct_m
34 use unit_oct_m
36
37 implicit none
38
39 private
40 public :: &
42
43contains
44
45 !---------------------------------------------------------------------------
49 subroutine test_optimizers(namespace)
50 type(namespace_t), intent(in) :: namespace
51
52 real(real64), dimension(2) :: x0, x, mass
53 real(real64), parameter :: tol_grad = 1.0e-6_real64
54 real(real64) :: dt, energy
55 integer :: ierr
56 integer, parameter :: max_iter = 1000
57
58 push_sub(test_optimizers)
59
60 call io_rm("optimization.log", namespace)
61
62 ! Testing the FIRE algorithm with the Rosenbrock function in 2D
63 ! Starting point
64 x0 = [3.0_real64, 4.0_real64]
65
66 ! Initial parameters for the FIRE algorithm with Velocity Verlet
67 x = x0
68 ! This timestep is taken to give a reasonable result. Too large vaues lead to unstable results
69 dt = 0.002
70 mass = m_one
71 call minimize_fire(2, 2, x, dt, tol_grad, max_iter, rosenbrock_gradient_2d, write_iter_info, energy, ierr, mass, 1)
72
73 write(message(1), "(a,i6,a)") "FIRE algorithm converged in ", -ierr, " iterations"
74 write(message(2), "(a,2(f8.4,a))") " The minimum is found to be at (", x(1), ", ", x(2), ")"
75 write(message(3), "(a)") " Analytical minimum is (1,1)."
76 call messages_info(3)
77
78 pop_sub(test_optimizers)
79 end subroutine test_optimizers
80
81 !---------------------------------------------------------------------------
88 subroutine rosenbrock_gradient_2d(n, x, val, getgrad, grad) !< get the new gradients given a new x
89 integer, intent(in) :: n
90 real(real64), intent(in) :: x(n)
91 real(real64), intent(inout) :: val
92 integer, intent(in) :: getgrad
93 real(real64), intent(inout) :: grad(n)
94
95 real(real64), parameter :: a = m_one
96 real(real64), parameter :: b = 100.0_real64
97
98 push_sub(f_rosenbrok)
99
100 ! Only the 2D version is implemented here
101 assert(n==2)
102
103 grad(1) = -m_two * (a-x(1)) - m_four * b * (x(2) -x(1)**2) * x(1)
104 grad(2) = m_two * b * (x(2) -x(1)**2)
105
106 pop_sub(f_rosenbrok)
107 end subroutine rosenbrock_gradient_2d
108
109 !---------------------------------------------------------------------------
113 subroutine write_iter_info(iter, n, val, maxdr, maxgrad, x) !< output for each iteration step
114 integer, intent(in) :: iter
115 integer, intent(in) :: n
116 real(real64), intent(in) :: val
117 real(real64), intent(in) :: maxdr
118 real(real64), intent(in) :: maxgrad
119 real(real64), intent(in) :: x(n)
120
121 integer :: iunit
122
123 if (mpi_grp_is_root(mpi_world)) then
124 iunit = io_open(trim('optimization.log'), global_namespace, action = 'write', position = 'append')
125 if (iter == 1) then
126 write(iunit, '(a10, 3a20)') '# iter','x','y', 'max_force'
127 end if
128 write(iunit, '(i10,3f20.10)') iter, x(1), x(2), maxgrad
129 call io_close(iunit)
130 end if
131
132 end subroutine write_iter_info
133
134end module minimizer_tests_oct_m
135
136!! Local Variables:
137!! mode: f90
138!! coding: utf-8
139!! End:
subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
Implementation of the Fast Inertial Relaxation Engine (FIRE)
Definition: minimizer.F90:482
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_four
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:188
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_rm(fname, namespace)
Definition: io.F90:385
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
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, public test_optimizers(namespace)
Unit tests for different optimizers.
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
type(namespace_t), public global_namespace
Definition: namespace.F90:132
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
for(l=0;l< nri;l++)
Definition: operate_inc.c:15