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 :: &
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 select case (method)
165 ierr = loct_minimize_direct(method, dim, x(1), step, toldr, maxiter, f, write_iter_info, minimum)
166 end select
167
168 pop_sub(minimize_multidim_nograd)
169
170 end subroutine minimize_multidim_nograd
171
172
173 subroutine minimize_multidim_nlopt(ierr, method, dim, x, step, toldr, maxiter, f, minimum, lb, ub)
174 integer, intent(out) :: ierr
175 integer, intent(in) :: method
176 integer, intent(in) :: dim
177 real(real64), intent(inout) :: x(:)
178 real(real64), intent(in) :: step
179 real(real64), intent(in) :: toldr
180 integer, intent(in) :: maxiter
181 interface
182 subroutine f(val, n, x, grad, need_gradient, f_data)
183 use iso_c_binding
184 real(c_double), intent(out) :: val
185 integer(c_int), intent(in) :: n
186 real(c_double), intent(in) :: x(*)
187 real(c_double), intent(out) :: grad(*)
188 integer(c_int), intent(in) :: need_gradient
189 type(c_ptr), intent(in) :: f_data
190 end subroutine f
191 end interface
192 real(real64), intent(out) :: minimum
193 real(real64), intent(in), optional :: lb(:), ub(:)
194#if defined(HAVE_NLOPT)
195
196 interface
197 subroutine nlo_create(opt, alg, n)
198 use iso_c_binding
199 type(c_ptr), intent(out) :: opt
200 integer(c_int), intent(in) :: alg
201 integer(c_int), intent(in) :: n
202 end subroutine nlo_create
204 subroutine nlo_set_lower_bounds(ret, opt, lower_bounds)
205 use iso_c_binding
206 integer(c_int), intent(out) :: ret
207 type(c_ptr), intent(inout) :: opt
208 real(c_double), intent(in) :: lower_bounds(*)
209 end subroutine nlo_set_lower_bounds
210
211 subroutine nlo_set_upper_bounds(ret, opt, upper_bounds)
212 use iso_c_binding
213 integer(c_int), intent(out) :: ret
214 type(c_ptr), intent(inout) :: opt
215 real(c_double), intent(in) :: upper_bounds(*)
216 end subroutine nlo_set_upper_bounds
217
218 subroutine nlo_set_min_objective(ret, opt, f, f_data)
219 use iso_c_binding
220 integer(c_int), intent(out) :: ret
221 type(c_ptr), intent(inout) :: opt
222 interface
223 subroutine f(val, n, x, grad, need_gradient, f_data)
224 use iso_c_binding
225 real(c_double), intent(out) :: val
226 integer(c_int), intent(in) :: n
227 real(c_double), intent(in) :: x(*)
228 real(c_double), intent(out) :: grad(*)
229 integer(c_int), intent(in) :: need_gradient
230 type(c_ptr), intent(in) :: f_data
231 end subroutine f
232 end interface
233 type(c_ptr), intent(in) :: f_data
234 end subroutine nlo_set_min_objective
235
236 subroutine nlo_set_xtol_abs1(ret, opt, xtol_abs)
237 use iso_c_binding
238 integer(c_int), intent(out) :: ret
239 type(c_ptr), intent(inout) :: opt
240 real(c_double), intent(in) :: xtol_abs
241 end subroutine nlo_set_xtol_abs1
242
243 subroutine nlo_set_initial_step1(ret, opt, initial_step1)
244 use iso_c_binding
245 integer(c_int), intent(out) :: ret
246 type(c_ptr), intent(inout) :: opt
247 real(c_double), intent(in) :: initial_step1
248 end subroutine nlo_set_initial_step1
249
250 subroutine nlo_set_maxeval(ret, opt, maxeval)
251 use iso_c_binding
252 integer(c_int), intent(out) :: ret
253 type(c_ptr), intent(inout) :: opt
254 integer(c_int), intent(in) :: maxeval
255 end subroutine nlo_set_maxeval
256
257 subroutine nlo_optimize(ret, opt, x, optf)
258 use iso_c_binding
259 integer(c_int), intent(out) :: ret
260 type(c_ptr), intent(inout) :: opt
261 real(c_double), intent(inout) :: x(*)
262 real(c_double), intent(out) :: optf
263 end subroutine nlo_optimize
264
265 subroutine nlo_destroy(opt)
266 use iso_c_binding
267 type(c_ptr), intent(inout) :: opt
268 end subroutine nlo_destroy
269 end interface
270
271 type(c_ptr) :: opt
272 integer :: ires
273 ! The following values are taken from the ''nlopt.f'' file installed by NLopt.
274 integer, parameter :: NLOPT_LD_LBFGS = 11
275 integer, parameter :: NLOPT_LN_BOBYQA = 34
276
277 select case (method)
279 call nlo_create(opt, nlopt_ln_bobyqa, dim)
281 call nlo_create(opt, nlopt_ld_lbfgs, dim)
282 end select
283
284 if (present(lb)) then
285 call nlo_set_lower_bounds(ires, opt, lb)
286 end if
287 if (present(ub)) then
288 call nlo_set_upper_bounds(ires, opt, ub)
289 end if
290
291 call nlo_set_min_objective(ires, opt, f, c_null_ptr)
292 ! This would set an inequality constraint (TODO)
293 ! call nlo_add_inequality_constraint(ires, opt, myconstraint, d1, 1.0e-8_real64)
294
295 call nlo_set_xtol_abs1(ires, opt, toldr)
296 call nlo_set_initial_step1(ires, opt, step)
297 call nlo_set_maxeval(ires, opt, maxiter)
298
299 call nlo_optimize(ires, opt, x, minimum)
300 ierr = ires
301 call nlo_destroy(opt)
302#else
303 ierr = 0
304#endif
305 end subroutine minimize_multidim_nlopt
306
307
308 !----------------------------------------------
309 subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
310 integer, intent(in) :: method
311 integer, intent(in) :: dim
312 real(real64), intent(inout) :: x(:)
313 real(real64), intent(in) :: step
314 real(real64), intent(in) :: line_tol
315 real(real64), intent(in) :: tolgrad
316 real(real64), intent(in) :: toldr
317 integer, intent(in) :: maxiter
318 procedure(minimizer_with_grad_i) :: f
319 procedure(info_i) :: write_iter_info
320 real(real64), intent(out) :: minimum
321 integer, intent(out) :: ierr
322
323 push_sub(minimize_multidim)
324
325 assert(ubound(x, dim = 1) >= dim)
326
327 select case (method)
329 call minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
330
331 case default
332 ierr = loct_minimize(method, dim, x(1), step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum)
333
334 end select
335
337
338 end subroutine minimize_multidim
339
340 !----------------------------------------------
341
342 subroutine minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
343 integer, intent(in) :: dim
344 real(real64), intent(inout) :: x(:)
345 real(real64), intent(in) :: step
346 integer, intent(in) :: maxiter
347 procedure(minimizer_with_grad_i) :: f
348 procedure(info_i) :: write_iter_info
349 real(real64), intent(out) :: minimum
350 integer, intent(out) :: ierr
351
352 integer :: iter
353 real(real64), allocatable :: grad(:)
354 real(real64) :: step2, maxgrad
355
356 push_sub(minimize_sd)
357
358 safe_allocate(grad(1:dim))
359
360 step2 = step*10.0_real64
361 do iter = 1, maxiter
362 call f(dim, x, minimum, 1, grad)
363
364 maxgrad = maxval(abs(grad))
366 call write_iter_info(iter, dim, minimum, maxgrad*step2, maxgrad, x)
368 x(1:dim) = x(1:dim) - step2*grad(1:dim)
369
370 step2 = step2*0.99_real64
371 end do
372 ierr = 0
373
374 pop_sub(minimize_sd)
375 end subroutine minimize_sd
376
377 !----------------------------------------------
388 subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
389 integer, intent(in) :: dim
390 integer, intent(in) :: space_dim
391 real(real64), intent(inout) :: x(:)
392 real(real64), intent(in) :: step
393 real(real64), intent(in) :: tolgrad
394 integer, intent(in) :: maxiter
395 procedure(minimizer_with_grad_i) :: f
396 procedure(info_i) :: write_iter_info
397 real(real64), intent(out) :: en
398 integer, intent(out) :: ierr
399 real(real64), intent(in) :: mass(:)
400 integer, intent(in) :: integrator
401
402 real(real64), allocatable :: grad(:), old_grad(:)
403
404 integer :: p_times, iter, ia, offset
405 real(real64) :: dt, alpha, p_value, dt_max
406 real(real64), allocatable :: grad_atoms(:), vel(:), dr_i(:), x_old(:), dr_atoms(:)
407
408 ! Parameters from the original paper
409 integer, parameter :: n_min = 5
410 real(real64), parameter :: f_alpha = 0.99_real64
411 real(real64), parameter :: f_inc = 1.1_real64
412 real(real64), parameter :: f_dec = 0.5_real64
413 real(real64), parameter :: alpha_start = 0.1_real64
414
415 real(real64), parameter :: maxmove = 0.2_real64 * p_ang
416
417 push_sub(minimize_fire)
418
419 ! dim must be a multiple of space_dim
420 assert(mod(dim,space_dim) == 0)
421 assert(size(x) == dim)
422
423 safe_allocate(grad_atoms(1:dim/space_dim))
424 safe_allocate(grad(1:dim))
425 safe_allocate(old_grad(1:dim))
426 safe_allocate(vel(1:dim))
427 safe_allocate(dr_atoms(1:dim/space_dim))
428 safe_allocate(x_old(1:dim))
429 safe_allocate(dr_i(1:dim))
430
431 ! Initial values
432 ierr = 0
433 x_old = x
434 vel = m_zero
435 dt = step
436 alpha = alpha_start
437 dr_i = m_zero
438 dt_max = 10.0_real64 * dt
439 p_times = 0
440 call f(dim, x, en, 1, grad)
441 grad = -grad
442 old_grad = grad
443
444 do iter = 1, maxiter
445
446 ! Perform the MD step: get new gradients from the new positions
447 select case (integrator)
448 ! Velocity verlet - update the positions
449 case (option__gofireintegrator__verlet)
450 dr_i(1:dim) = vel(1:dim) * dt + m_half * grad(1:dim) / mass(1:dim) * dt **2
451 case (option__gofireintegrator__euler)
452 ! Euler method
453 dr_i(1:dim) = vel(1:dim)*dt
454 end select
455
456 ! Get the norms of the displacements for each atom
457 do ia = 1, dim/space_dim
458 offset = space_dim * (ia -1)
459 dr_atoms(ia) = norm2(dr_i(offset+1:offset+space_dim))
460 ! Rescale the displacement to avoid too large changes
461 if (dr_atoms(ia) > maxmove) then
462 dr_i(offset+1:offset+space_dim) = maxmove * dr_i(offset+1:offset+space_dim) / dr_atoms(ia)
463 dr_atoms(ia) = maxmove
464 end if
465 end do
466
467 ! Get the new position
468 x(1:dim) = x_old(1:dim) + dr_i(1:dim)
469
470 ! Get the new gradient
471 call f(dim, x, en, 1, grad)
472 grad = -grad
473
474 ! Molecular dynamics step for getting the new velocities
475 select case (integrator)
476 case (option__gofireintegrator__verlet)
477 ! Velocity Verlet - update velocities
478 vel(1:dim) = vel(1:dim) + m_half*(grad(1:dim) + old_grad(1:dim))*dt/mass(1:dim)
479 case (option__gofireintegrator__euler)
480 ! Euler method
481 vel(1:dim) = vel(1:dim) + grad(1:dim)*dt/mass(1:dim)
482 end select
483
484 ! Now that we performed the MD part, we correct the velocitites
485
486 ! Perform step F1: compute P = F.v
487 p_value = dot_product(grad, vel)
488
489 ! Step F2: v-> (1-\alpha)v + \alpha \hat{F}.|v|
490 if (iter > 1) then
491 vel(1:dim) = (m_one - alpha) * vel(1:dim) + alpha * grad(1:dim) * lalg_nrm2(dim, vel) / lalg_nrm2(dim,grad)
492 end if
493
494 ! Perform step F3
495 if (p_value > m_zero) then
496 p_times = p_times + 1
497 if (p_times > n_min) then
498 dt = min(dt * f_inc , dt_max)
499 alpha = alpha * f_alpha
500 end if
501
502 else ! or step F4
503 p_times = 0
504 dt = dt * f_dec
505 vel = m_zero
506 alpha = alpha_start
507 end if
508
509 ! Get the norms of the gradients for each atom
510 do ia = 0, dim/space_dim - 1
511 grad_atoms(ia+1) = norm2(grad(space_dim*ia+1:space_dim*ia+space_dim))
512 end do
513
514 ! Output for each iteration step
515 call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
516
517 ! Check convergence
518 if (maxval(abs(grad_atoms(1:))) < tolgrad) exit
519
520 x_old(1:dim) = x(1:dim)
521 old_grad = grad
522 end do
523
524 ierr = -iter
525
526 safe_deallocate_a(dr_atoms)
527 safe_deallocate_a(x_old)
528 safe_deallocate_a(dr_i)
529 safe_deallocate_a(vel)
530 safe_deallocate_a(grad)
531 safe_deallocate_a(old_grad)
532 safe_deallocate_a(grad_atoms)
533
534 pop_sub(minimize_fire)
535
536 end subroutine minimize_fire
537
538end module minimizer_oct_m
539
540!! Local Variables:
541!! mode: f90
542!! coding: utf-8
543!! End:
Returns the euclidean norm of a vector.
Definition: lalg_basic.F90:197
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:403
subroutine minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:436
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_zero
Definition: global.F90:187
real(real64), parameter, public p_ang
Definition: global.F90:219
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
integer minmethod_fr_cg
Definition: minimizer.F90:134
integer minmethod_fire
Definition: minimizer.F90:134
integer minmethod_pr_cg
Definition: minimizer.F90:134
integer minmethod_nmsimplex
Definition: minimizer.F90:134
integer minmethod_bfgs2
Definition: minimizer.F90:134
integer minmethod_nlopt_lbfgs
Definition: minimizer.F90:134
integer minmethod_nlopt_bobyqa
Definition: minimizer.F90:134
integer minmethod_bfgs
Definition: minimizer.F90:134
integer minmethod_sd_native
Definition: minimizer.F90:134