Octopus
minimizer.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
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
21module minimizer_oct_m
22 use debug_oct_m
23 use global_oct_m
24 use iso_c_binding
25 use, intrinsic :: iso_fortran_env
29
30 implicit none
31
32 private
33 public :: &
35 minimize_fire, &
36 minimize_multidim, &
37 minimize_multidim_nograd, &
38 minimize_multidim_nlopt
39
40
41 integer, public, parameter :: &
42 MINMETHOD_STEEPEST_DESCENT = 1, &
43 minmethod_fr_cg = 2, &
44 minmethod_pr_cg = 3, &
45 minmethod_bfgs = 4, &
46 minmethod_bfgs2 = 5, &
50 minmethod_fire = 8, &
52
53 abstract interface
54 subroutine minimizer_function_i(n, x, val)
55 import real64
56 implicit none
57 integer :: n
58 real(real64) :: x(n)
59 real(real64) :: val
60 end subroutine minimizer_function_i
61 subroutine minimizer_with_grad_i(n, x, val, getgrad, grad)
62 import real64
63 implicit none
64 integer, intent(in) :: n
65 real(real64), intent(in) :: x(n)
66 real(real64), intent(inout) :: val
67 integer, intent(in) :: getgrad
68 real(real64), intent(inout) :: grad(n)
69 end subroutine minimizer_with_grad_i
70 subroutine info_i(iter, n, val, maxdr, maxgrad, x)
71 import real64
72 implicit none
73 integer, intent(in) :: iter
74 integer, intent(in) :: n
75 real(real64), intent(in) :: val
76 real(real64), intent(in) :: maxdr
77 real(real64), intent(in) :: maxgrad
78 real(real64), intent(in) :: x(n)
79 end subroutine info_i
80 subroutine info_no_grad_i(iter, n, val, maxdr, x)
81 import real64
82 implicit none
83 integer, intent(in) :: iter
84 integer, intent(in) :: n
85 real(real64), intent(in) :: val
86 real(real64), intent(in) :: maxdr
87 real(real64), intent(in) :: x(n)
88 end subroutine info_no_grad_i
89 end interface
90
91
92 interface loct_1dminimize
93 subroutine oct_1dminimize(a, b, m, f, status)
94 import real64
95 implicit none
96 real(real64), intent(inout) :: a, b, m
97 interface
98 subroutine f(x, fx)
99 import real64
100 implicit none
101 real(real64), intent(in) :: x
102 real(real64), intent(out) :: fx
103 end subroutine f
104 end interface
105 integer, intent(out) :: status
106 end subroutine oct_1dminimize
107 end interface loct_1dminimize
108
109 interface loct_minimize
110 integer function oct_minimize(method, dim, x, step, line_tol, &
111 tolgrad, toldr, maxiter, f, write_iter_info, minimum)
112 import minimizer_with_grad_i, info_i, real64
113 implicit none
114 integer, intent(in) :: method
115 integer, intent(in) :: dim
116 real(real64), intent(inout) :: x
117 real(real64), intent(in) :: step
118 real(real64), intent(in) :: line_tol
119 real(real64), intent(in) :: tolgrad
120 real(real64), intent(in) :: toldr
121 integer, intent(in) :: maxiter
122 procedure(minimizer_with_grad_i) :: f
123 procedure(info_i) :: write_iter_info
124 real(real64), intent(out) :: minimum
125 end function oct_minimize
126 end interface loct_minimize
127
128 interface loct_minimize_direct
129 function oct_minimize_direct(method, dim, x, step, toldr, maxiter, f, write_iter_info, minimum)
131 implicit none
132 integer :: oct_minimize_direct
133 integer, intent(in) :: method
134 integer, intent(in) :: dim
135 real(real64), intent(inout) :: x
136 real(real64), intent(in) :: step
137 real(real64), intent(in) :: toldr
138 integer, intent(in) :: maxiter
139 procedure(minimizer_function_i) :: f
140 procedure(info_no_grad_i) :: write_iter_info
141 real(real64), intent(out) :: minimum
142 end function oct_minimize_direct
143 end interface loct_minimize_direct
144
145contains
146
147 subroutine minimize_multidim_nograd(method, dim, x, step, toldr, maxiter, f, write_iter_info, minimum, ierr)
148 integer, intent(in) :: method
149 integer, intent(in) :: dim
150 real(real64), intent(inout) :: x(:)
151 real(real64), intent(in) :: step
152 real(real64), intent(in) :: toldr
153 integer, intent(in) :: maxiter
154 procedure(minimizer_function_i) :: f
155 procedure(info_no_grad_i) :: write_iter_info
156 real(real64), intent(out) :: minimum
157 integer, intent(out) :: ierr
158
159 push_sub(minimize_multidim_nograd)
160
161 assert(ubound(x, dim = 1) >= dim)
162
163 ierr = 0
164 select case (method)
166 ierr = loct_minimize_direct(method, dim, x(1), step, toldr, maxiter, f, write_iter_info, minimum)
167 end select
168
169 pop_sub(minimize_multidim_nograd)
170
171 end subroutine minimize_multidim_nograd
172
173
174 subroutine minimize_multidim_nlopt(ierr, method, dim, x, step, toldr, maxiter, f, minimum, lb, ub)
175 integer, intent(out) :: ierr
176 integer, intent(in) :: method
177 integer, intent(in) :: dim
178 real(real64), contiguous, intent(inout) :: x(:)
179 real(real64), intent(in) :: step
180 real(real64), intent(in) :: toldr
181 integer, intent(in) :: maxiter
182 interface
183 subroutine f(val, n, x, grad, need_gradient, f_data)
184 use iso_c_binding
185 real(c_double), intent(out) :: val
186 integer(c_int), intent(in) :: n
187 real(c_double), intent(in) :: x(*)
188 real(c_double), intent(out) :: grad(*)
189 integer(c_int), intent(in) :: need_gradient
190 type(c_ptr), intent(in) :: f_data
191 end subroutine f
192 end interface
193 real(real64), intent(out) :: minimum
194 real(real64), intent(in), optional :: lb(:), ub(:)
195#if defined(HAVE_NLOPT)
196
197 interface
198 subroutine nlo_create(opt, alg, n)
199 use iso_c_binding
200 type(c_ptr), intent(out) :: opt
201 integer(c_int), intent(in) :: alg
202 integer(c_int), intent(in) :: n
203 end subroutine nlo_create
205 subroutine nlo_set_lower_bounds(ret, opt, lower_bounds)
206 use iso_c_binding
207 integer(c_int), intent(out) :: ret
208 type(c_ptr), intent(inout) :: opt
209 real(c_double), intent(in) :: lower_bounds(*)
210 end subroutine nlo_set_lower_bounds
211
212 subroutine nlo_set_upper_bounds(ret, opt, upper_bounds)
213 use iso_c_binding
214 integer(c_int), intent(out) :: ret
215 type(c_ptr), intent(inout) :: opt
216 real(c_double), intent(in) :: upper_bounds(*)
217 end subroutine nlo_set_upper_bounds
218
219 subroutine nlo_set_min_objective(ret, opt, f, f_data)
220 use iso_c_binding
221 integer(c_int), intent(out) :: ret
222 type(c_ptr), intent(inout) :: opt
223 interface
224 subroutine f(val, n, x, grad, need_gradient, f_data)
225 use iso_c_binding
226 real(c_double), intent(out) :: val
227 integer(c_int), intent(in) :: n
228 real(c_double), intent(in) :: x(*)
229 real(c_double), intent(out) :: grad(*)
230 integer(c_int), intent(in) :: need_gradient
231 type(c_ptr), intent(in) :: f_data
232 end subroutine f
233 end interface
234 type(c_ptr), intent(in) :: f_data
235 end subroutine nlo_set_min_objective
236
237 subroutine nlo_set_xtol_abs1(ret, opt, xtol_abs)
238 use iso_c_binding
239 integer(c_int), intent(out) :: ret
240 type(c_ptr), intent(inout) :: opt
241 real(c_double), intent(in) :: xtol_abs
242 end subroutine nlo_set_xtol_abs1
243
244 subroutine nlo_set_initial_step1(ret, opt, initial_step1)
245 use iso_c_binding
246 integer(c_int), intent(out) :: ret
247 type(c_ptr), intent(inout) :: opt
248 real(c_double), intent(in) :: initial_step1
249 end subroutine nlo_set_initial_step1
250
251 subroutine nlo_set_maxeval(ret, opt, maxeval)
252 use iso_c_binding
253 integer(c_int), intent(out) :: ret
254 type(c_ptr), intent(inout) :: opt
255 integer(c_int), intent(in) :: maxeval
256 end subroutine nlo_set_maxeval
257
258 subroutine nlo_optimize(ret, opt, x, optf)
259 use iso_c_binding
260 integer(c_int), intent(out) :: ret
261 type(c_ptr), intent(inout) :: opt
262 real(c_double), intent(inout) :: x(*)
263 real(c_double), intent(out) :: optf
264 end subroutine nlo_optimize
265
266 subroutine nlo_destroy(opt)
267 use iso_c_binding
268 type(c_ptr), intent(inout) :: opt
269 end subroutine nlo_destroy
270 end interface
271
272 type(c_ptr) :: opt
273 integer :: ires
274 ! The following values are taken from the ''nlopt.f'' file installed by NLopt.
275 integer, parameter :: NLOPT_LD_LBFGS = 11
276 integer, parameter :: NLOPT_LN_BOBYQA = 34
277
278 opt = c_null_ptr
279 select case (method)
281 call nlo_create(opt, nlopt_ln_bobyqa, dim)
283 call nlo_create(opt, nlopt_ld_lbfgs, dim)
284 end select
285
286 if (present(lb)) then
287 call nlo_set_lower_bounds(ires, opt, lb)
288 end if
289 if (present(ub)) then
290 call nlo_set_upper_bounds(ires, opt, ub)
291 end if
292
293 call nlo_set_min_objective(ires, opt, f, c_null_ptr)
294 ! This would set an inequality constraint (TODO)
295 ! call nlo_add_inequality_constraint(ires, opt, myconstraint, d1, 1.0e-8_real64)
296
297 call nlo_set_xtol_abs1(ires, opt, toldr)
298 call nlo_set_initial_step1(ires, opt, step)
299 call nlo_set_maxeval(ires, opt, maxiter)
300
301 call nlo_optimize(ires, opt, x, minimum)
302 ierr = ires
303 call nlo_destroy(opt)
304#else
305 ierr = 0
306 minimum = -m_huge
307#endif
308 end subroutine minimize_multidim_nlopt
309
310
311 !----------------------------------------------
312 subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
313 integer, intent(in) :: method
314 integer, intent(in) :: dim
315 real(real64), intent(inout) :: x(:)
316 real(real64), intent(in) :: step
317 real(real64), intent(in) :: line_tol
318 real(real64), intent(in) :: tolgrad
319 real(real64), intent(in) :: toldr
320 integer, intent(in) :: maxiter
321 procedure(minimizer_with_grad_i) :: f
322 procedure(info_i) :: write_iter_info
323 real(real64), intent(out) :: minimum
324 integer, intent(out) :: ierr
325
326 push_sub(minimize_multidim)
327
328 assert(ubound(x, dim = 1) >= dim)
329
330 select case (method)
332 call minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
333
334 case default
335 ierr = loct_minimize(method, dim, x(1), step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum)
336
337 end select
338
339 pop_sub(minimize_multidim)
340
341 end subroutine minimize_multidim
342
343 !----------------------------------------------
344
345 subroutine minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
346 integer, intent(in) :: dim
347 real(real64), intent(inout) :: x(:)
348 real(real64), intent(in) :: step
349 integer, intent(in) :: maxiter
350 procedure(minimizer_with_grad_i) :: f
351 procedure(info_i) :: write_iter_info
352 real(real64), intent(out) :: minimum
353 integer, intent(out) :: ierr
354
355 integer :: iter
356 real(real64), allocatable :: grad(:)
357 real(real64) :: step2, maxgrad
358
359 push_sub(minimize_sd)
360
361 safe_allocate(grad(1:dim))
362
363 step2 = step*10.0_real64
364 do iter = 1, maxiter
365 call f(dim, x, minimum, 1, grad)
366
367 maxgrad = maxval(abs(grad))
368
369 call write_iter_info(iter, dim, minimum, maxgrad*step2, maxgrad, x)
370
371 x(1:dim) = x(1:dim) - step2*grad(1:dim)
372
373 step2 = step2*0.99_real64
374 end do
375 ierr = 0
376
377 pop_sub(minimize_sd)
378 end subroutine minimize_sd
379
380 !----------------------------------------------
394 subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, gradf, write_iter_info, en, ierr, mass, integrator)
395 integer, intent(in) :: dim
396 integer, intent(in) :: space_dim
397 real(real64), intent(inout) :: x(:)
398 real(real64), intent(in) :: step
399 real(real64), intent(in) :: tolgrad
400 integer, intent(in) :: maxiter
401 procedure(minimizer_with_grad_i) :: gradf
402 procedure(info_i) :: write_iter_info
403 real(real64), intent(out) :: en
404 integer, intent(out) :: ierr
405 real(real64), intent(in) :: mass(:)
406 integer, intent(in) :: integrator
407
408 real(real64), allocatable :: grad(:)
409
410 integer :: p_times, n_times, iter, ia
411 real(real64) :: dt, alpha, p_value, dt_max, dt_min, mix
412 real(real64), allocatable :: grad_atoms(:), vel(:), dr_i(:), x_old(:), dr_atoms(:)
413
414 ! Parameters from the original paper
415 integer, parameter :: n_min = 5, n_decrease_max = 500
416 real(real64), parameter :: f_alpha = 0.99_real64
417 real(real64), parameter :: f_inc = 1.1_real64
418 real(real64), parameter :: f_dec = 0.5_real64
419 real(real64), parameter :: alpha_start = 0.25_real64
420
421 real(real64), parameter :: tol_zero_forces = 1e-10_real64
422
423 push_sub(minimize_fire)
424
425 ! dim must be a multiple of space_dim
426 assert(mod(dim,space_dim) == 0)
427 assert(size(x) == dim)
428
429 safe_allocate(grad_atoms(1:dim/space_dim))
430 safe_allocate(grad(1:dim))
431 safe_allocate(vel(1:dim))
432 safe_allocate(dr_atoms(1:dim/space_dim))
433 safe_allocate(x_old(1:dim))
434 safe_allocate(dr_i(1:dim))
435
436 ! Initial values
437 ierr = 0
438 vel = m_zero
439 dt = step
440 alpha = alpha_start
441 dr_i = m_zero
442 dt_max = 10.0_real64 * dt
443 dt_min = 0.02_real64 * dt
444 p_times = 0
445 n_times = 0
446 call gradf(dim, x, en, 1, grad)
447 grad = -grad
448
449 do iter = 1, maxiter
450
451 ! Perform step F1: compute P = F.v
452 p_value = dot_product(grad, vel)
453
454 ! Perform step F3
455 if (p_value > m_zero) then
456 p_times = p_times + 1
457 n_times = 0
458 if (p_times > n_min) then
459 dt = min(dt * f_inc , dt_max)
460 alpha = alpha * f_alpha
461 end if
462 else ! or step F4
463 p_times = 0
464 n_times = n_times + 1
465
466 if (n_times > n_decrease_max) exit ! If we are stuck, we exit
467
468 if (iter >= n_min) then ! No decrease for the first few steps
469 dt = max(dt * f_dec, dt_min)
470 alpha = alpha_start
471 end if
472 x = x - m_half * vel * dt ! Correct uphill motion
473 vel = m_zero
474 end if
475
476 ! Store old values
477 x_old = x
478
479 ! Perform the MD step: get new gradients from the new positions
480 ! Step F2: v-> (1-\alpha)v + \alpha \hat{F}.|v| is done inside the MD step
481 ! We follow the integration schemes of the FIRE 2.0 paper
482
483 ! Prevent division by zero in the mixing
484 mix = lalg_nrm2(dim,grad)
485 if (mix > tol_zero_forces) then
486 mix = m_one / mix
487 else
488 mix = m_zero
489 end if
490
491 select case (integrator)
492 case (option__gofireintegrator__verlet) ! Velocity Verlet, see Algorithm 6
493 ! Velocity Verlet - update velocities 1
494 vel = vel + m_half * grad * dt / mass
495
496 ! Mixing
497 if (p_value > m_zero) then
498 vel = (m_one - alpha) * vel + alpha * mix * grad * lalg_nrm2(dim, vel)
499 end if
500
501 ! Get the new position
502 dr_i = vel * dt
503 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
504 x = x_old + dr_i
505
506 ! Get the new gradient
507 call gradf(dim, x, en, 1, grad)
508 grad = -grad
509
510 ! Velocity Verlet - update velocities 2
511 vel = vel + m_half * grad * dt / mass
512
513 case (option__gofireintegrator__euler) ! Explicit Euler method, see Algorithm 3 in Fire2.0 paper
514 ! Mixing
515 if (p_value > m_zero) then
516 vel = (m_one - alpha) * vel + alpha * mix * grad * lalg_nrm2(dim, vel)
517 end if
518
519 ! Get the new position
520 dr_i = vel * dt
521 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
522 x = x_old + dr_i
523
524 ! Get the new gradient
525 call gradf(dim, x, en, 1, grad)
526 grad = -grad
527
528 ! Velocity Verlet - update velocities
529 vel = vel + grad * dt / mass
530
531
532 case (option__gofireintegrator__semi_implicit_euler) ! Semi-implicit Euler method, see Algorithm 4
533 ! Velocity Verlet - update velocities
534 vel = vel + grad * dt / mass
535
536 ! Mixing
537 if (p_value > m_zero) then
538 vel = (m_one - alpha) * vel + alpha * mix * grad * lalg_nrm2(dim, vel)
539 end if
540
541 ! Get the new position
542 dr_i = vel * dt
543 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
544 x = x_old + dr_i
545
546 ! Get the new gradient
547 call gradf(dim, x, en, 1, grad)
548 grad = -grad
549 end select
550
551 ! Get the norms of the gradients for each atom
552 do ia = 0, dim/space_dim - 1
553 grad_atoms(ia+1) = norm2(grad(space_dim*ia+1:space_dim*ia+space_dim))
554 end do
555
556 ! Output for each iteration step
557 call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
558
559 ! Check convergence
560 if (maxval(abs(grad_atoms(1:))) < tolgrad) exit
561 end do
562
563 ierr = -iter
565 safe_deallocate_a(dr_atoms)
566 safe_deallocate_a(x_old)
567 safe_deallocate_a(dr_i)
568 safe_deallocate_a(vel)
569 safe_deallocate_a(grad)
570 safe_deallocate_a(grad_atoms)
571
572 pop_sub(minimize_fire)
573 end subroutine minimize_fire
574
576 ! TODO: we should instead change the velocities,
577 ! see S. Echeverri Restrepo and P. Andric, Computational Materials Science 218 (2023) 111978
578 subroutine prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
579 integer, intent(in) :: dim
580 integer, intent(in) :: space_dim
581 real(real64), intent(inout) :: dr_i(:), dr_atoms(:)
582
583 integer :: i1, i2, ia
584
585 real(real64), parameter :: maxmove = 0.2_real64 * p_ang
586
587 ! Get the norms of the displaments for each atoms
588 do ia = 0, dim/space_dim - 1
589 i1 = space_dim*ia+1
590 i2 = space_dim*ia+space_dim
591 dr_atoms(ia+1) = norm2(dr_i(i1:i2))
592 ! Rescale the displacement to avoid too large changes
593 if (dr_atoms(ia+1) > maxmove) then
594 dr_i(i1:i2) = maxmove * dr_i(i1:i2) / dr_atoms(ia+1)
595 dr_atoms(ia+1) = maxmove
596 end if
597 end do
598 end subroutine prevent_large_changes
599
600end module minimizer_oct_m
601
602!! Local Variables:
603!! mode: f90
604!! coding: utf-8
605!! End:
Returns the euclidean norm of a vector.
Definition: lalg_basic.F90:200
real(real64), parameter, public m_huge
Definition: global.F90:218
real(real64), parameter, public m_zero
Definition: global.F90:200
real(real64), parameter, public p_ang
Definition: global.F90:238
real(real64), parameter, public m_half
Definition: global.F90:206
real(real64), parameter, public m_one
Definition: global.F90:201
integer, parameter, public minmethod_fr_cg
Definition: minimizer.F90:136
integer, parameter, public minmethod_bfgs
Definition: minimizer.F90:136
integer, parameter, public minmethod_nmsimplex
Definition: minimizer.F90:136
integer, parameter, public minmethod_bfgs2
Definition: minimizer.F90:136
integer, parameter, public minmethod_nlopt_bobyqa
Definition: minimizer.F90:136
integer, parameter, public minmethod_nlopt_lbfgs
Definition: minimizer.F90:136
integer, parameter, public minmethod_pr_cg
Definition: minimizer.F90:136
integer, parameter, public minmethod_fire
Definition: minimizer.F90:136
integer, parameter, public minmethod_sd_native
Definition: minimizer.F90:136
static double f(double w, void *p)