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 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
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
204
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 !----------------------------------------------
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)
367 maxgrad = maxval(abs(grad))
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 !----------------------------------------------
391 subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
392 integer, intent(in) :: dim
393 integer, intent(in) :: space_dim
394 real(real64), intent(inout) :: x(:)
395 real(real64), intent(in) :: step
396 real(real64), intent(in) :: tolgrad
397 integer, intent(in) :: maxiter
398 procedure(minimizer_with_grad_i) :: f
399 procedure(info_i) :: write_iter_info
400 real(real64), intent(out) :: en
401 integer, intent(out) :: ierr
402 real(real64), intent(in) :: mass(:)
403 integer, intent(in) :: integrator
404
405 real(real64), allocatable :: grad(:), old_grad(:)
406
407 integer :: p_times, iter, ia, offset
408 real(real64) :: dt, alpha, p_value, dt_max
409 real(real64), allocatable :: grad_atoms(:), vel(:), dr_i(:), x_old(:), dr_atoms(:)
410
411 ! Parameters from the original paper
412 integer, parameter :: n_min = 5
413 real(real64), parameter :: f_alpha = 0.99_real64
414 real(real64), parameter :: f_inc = 1.1_real64
415 real(real64), parameter :: f_dec = 0.5_real64
416 real(real64), parameter :: alpha_start = 0.1_real64
417
418 real(real64), parameter :: maxmove = 0.2_real64 * p_ang
419
420 push_sub(minimize_fire)
421
422 ! dim must be a multiple of space_dim
423 assert(mod(dim,space_dim) == 0)
424 assert(size(x) == dim)
425
426 safe_allocate(grad_atoms(1:dim/space_dim))
427 safe_allocate(grad(1:dim))
428 safe_allocate(old_grad(1:dim))
429 safe_allocate(vel(1:dim))
430 safe_allocate(dr_atoms(1:dim/space_dim))
431 safe_allocate(x_old(1:dim))
432 safe_allocate(dr_i(1:dim))
433
434 ! Initial values
435 ierr = 0
436 x_old(:) = x(:)
437 vel = m_zero
438 dt = step
439 alpha = alpha_start
440 dr_i = m_zero
441 dt_max = 10.0_real64 * dt
442 p_times = 0
443 call f(dim, x, en, 1, grad)
444 grad = -grad
445 old_grad(:) = grad(:)
446
447 do iter = 1, maxiter
448
449 ! Perform the MD step: get new gradients from the new positions
450 select case (integrator)
451 ! Velocity verlet - update the positions
452 case (option__gofireintegrator__verlet)
453 dr_i(1:dim) = vel(1:dim) * dt + m_half * grad(1:dim) / mass(1:dim) * dt **2
454 case (option__gofireintegrator__euler)
455 ! Euler method
456 dr_i(1:dim) = vel(1:dim)*dt
457 end select
458
459 ! Get the norms of the displacements for each atom
460 do ia = 1, dim/space_dim
461 offset = space_dim * (ia -1)
462 dr_atoms(ia) = norm2(dr_i(offset+1:offset+space_dim))
463 ! Rescale the displacement to avoid too large changes
464 if (dr_atoms(ia) > maxmove) then
465 dr_i(offset+1:offset+space_dim) = maxmove * dr_i(offset+1:offset+space_dim) / dr_atoms(ia)
466 dr_atoms(ia) = maxmove
467 end if
468 end do
469
470 ! Get the new position
471 x(1:dim) = x_old(1:dim) + dr_i(1:dim)
472
473 ! Get the new gradient
474 call f(dim, x, en, 1, grad)
475 grad = -grad
476
477 ! Molecular dynamics step for getting the new velocities
478 select case (integrator)
479 case (option__gofireintegrator__verlet)
480 ! Velocity Verlet - update velocities
481 vel(1:dim) = vel(1:dim) + m_half*(grad(1:dim) + old_grad(1:dim))*dt/mass(1:dim)
482 case (option__gofireintegrator__euler)
483 ! Euler method
484 vel(1:dim) = vel(1:dim) + grad(1:dim)*dt/mass(1:dim)
485 end select
486
487 ! Now that we performed the MD part, we correct the velocitites
488
489 ! Perform step F1: compute P = F.v
490 p_value = dot_product(grad, vel)
491
492 ! Step F2: v-> (1-\alpha)v + \alpha \hat{F}.|v|
493 if (iter > 1) then
494 vel(1:dim) = (m_one - alpha) * vel(1:dim) + alpha * grad(1:dim) * lalg_nrm2(dim, vel) / lalg_nrm2(dim,grad)
495 end if
496
497 ! Perform step F3
498 if (p_value > m_zero) then
499 p_times = p_times + 1
500 if (p_times > n_min) then
501 dt = min(dt * f_inc , dt_max)
502 alpha = alpha * f_alpha
503 end if
504
505 else ! or step F4
506 p_times = 0
507 dt = dt * f_dec
508 vel = m_zero
509 alpha = alpha_start
510 end if
511
512 ! Get the norms of the gradients for each atom
513 do ia = 0, dim/space_dim - 1
514 grad_atoms(ia+1) = norm2(grad(space_dim*ia+1:space_dim*ia+space_dim))
515 end do
516
517 ! Output for each iteration step
518 call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
519
520 ! Check convergence
521 if (maxval(abs(grad_atoms(1:))) < tolgrad) exit
522
523 x_old(1:dim) = x(1:dim)
524 old_grad(:) = grad(:)
525 end do
526
527 ierr = -iter
528
529 safe_deallocate_a(dr_atoms)
530 safe_deallocate_a(x_old)
531 safe_deallocate_a(dr_i)
532 safe_deallocate_a(vel)
533 safe_deallocate_a(grad)
534 safe_deallocate_a(old_grad)
535 safe_deallocate_a(grad_atoms)
536
537 pop_sub(minimize_fire)
538
539 end subroutine minimize_fire
540
541end module minimizer_oct_m
542
543!! Local Variables:
544!! mode: f90
545!! coding: utf-8
546!! End:
Returns the euclidean norm of a vector.
Definition: lalg_basic.F90:198
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:406
subroutine minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
Definition: minimizer.F90:439
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:485
real(real64), parameter, public m_huge
Definition: global.F90:206
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public p_ang
Definition: global.F90:220
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
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
static double f(double w, void *p)