25 use,
intrinsic :: iso_fortran_env
37 minimize_multidim_nograd, &
38 minimize_multidim_nlopt
41 integer,
public,
parameter :: &
42 MINMETHOD_STEEPEST_DESCENT = 1, &
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)
70 subroutine info_i(iter, n, val, maxdr, maxgrad, x)
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)
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)
93 subroutine oct_1dminimize(a, b, m, f, status)
96 real(real64),
intent(inout) :: a, b, m
101 real(real64),
intent(in) :: x
102 real(real64),
intent(out) :: fx
105 integer,
intent(out) :: status
106 end subroutine oct_1dminimize
109 interface loct_minimize
110 integer function oct_minimize(method, dim, x, step, line_tol, &
111 tolgrad, toldr, maxiter, f, write_iter_info, minimum)
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
128 interface loct_minimize_direct
129 function oct_minimize_direct(method, dim, x, step, toldr, maxiter, f, write_iter_info, minimum)
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
141 real(real64),
intent(out) :: minimum
142 end function oct_minimize_direct
143 end interface loct_minimize_direct
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
156 real(real64),
intent(out) :: minimum
157 integer,
intent(out) :: ierr
159 push_sub(minimize_multidim_nograd)
161 assert(ubound(x, dim = 1) >= dim)
165 ierr = loct_minimize_direct(method, dim, x(1), step, toldr, maxiter, f, write_iter_info, minimum)
168 pop_sub(minimize_multidim_nograd)
170 end subroutine minimize_multidim_nograd
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
182 subroutine f(val, n, x, grad, need_gradient, f_data)
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
192 real(real64),
intent(out) :: minimum
193 real(real64),
intent(in),
optional :: lb(:), ub(:)
194#if defined(HAVE_NLOPT)
197 subroutine nlo_create(opt, alg, n)
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)
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
211 subroutine nlo_set_upper_bounds(ret, opt, upper_bounds)
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
218 subroutine nlo_set_min_objective(ret, opt, f, f_data)
220 integer(c_int),
intent(out) :: ret
221 type(c_ptr),
intent(inout) :: opt
223 subroutine f(val, n, x, grad, need_gradient, f_data)
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
233 type(c_ptr),
intent(in) :: f_data
234 end subroutine nlo_set_min_objective
236 subroutine nlo_set_xtol_abs1(ret, opt, xtol_abs)
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
243 subroutine nlo_set_initial_step1(ret, opt, initial_step1)
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
250 subroutine nlo_set_maxeval(ret, opt, maxeval)
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
257 subroutine nlo_optimize(ret, opt, x, optf)
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
265 subroutine nlo_destroy(opt)
267 type(c_ptr),
intent(inout) :: opt
268 end subroutine nlo_destroy
274 integer,
parameter :: NLOPT_LD_LBFGS = 11
275 integer,
parameter :: NLOPT_LN_BOBYQA = 34
279 call nlo_create(opt, nlopt_ln_bobyqa, dim)
281 call nlo_create(opt, nlopt_ld_lbfgs, dim)
284 if (
present(lb))
then
285 call nlo_set_lower_bounds(ires, opt, lb)
287 if (
present(ub))
then
288 call nlo_set_upper_bounds(ires, opt, ub)
291 call nlo_set_min_objective(ires, opt,
f, c_null_ptr)
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)
299 call nlo_optimize(ires, opt, x, minimum)
301 call nlo_destroy(opt)
305 end subroutine minimize_multidim_nlopt
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
319 procedure(
info_i) :: write_iter_info
320 real(real64),
intent(out) :: minimum
321 integer,
intent(out) :: ierr
323 push_sub(minimize_multidim)
325 assert(ubound(x, dim = 1) >= dim)
329 call minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
332 ierr = loct_minimize(method, dim, x(1), step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum)
336 pop_sub(minimize_multidim)
338 end subroutine minimize_multidim
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
348 procedure(
info_i) :: write_iter_info
349 real(real64),
intent(out) :: minimum
350 integer,
intent(out) :: ierr
353 real(real64),
allocatable :: grad(:)
354 real(real64) :: step2, maxgrad
356 push_sub(minimize_sd)
358 safe_allocate(grad(1:dim))
360 step2 = step*10.0_real64
362 call f(dim, x, minimum, 1, grad)
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)
370 step2 = step2*0.99_real64
375 end subroutine minimize_sd
391 subroutine minimize_fire(dim, space_dim, x, step, tolgrad, maxiter, gradf, 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
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
405 real(real64),
allocatable :: grad(:)
407 integer :: p_times, n_times, iter, ia
408 real(real64) :: dt, alpha, p_value, dt_max, dt_min, mix
409 real(real64),
allocatable :: grad_atoms(:), vel(:), dr_i(:), x_old(:), dr_atoms(:)
412 integer,
parameter :: n_min = 5, n_decrease_max = 500
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.25_real64
418 real(real64),
parameter :: tol_zero_forces = 1e-10_real64
420 push_sub(minimize_fire)
423 assert(mod(dim,space_dim) == 0)
424 assert(
size(x) == dim)
426 safe_allocate(grad_atoms(1:dim/space_dim))
427 safe_allocate(grad(1:dim))
428 safe_allocate(vel(1:dim))
429 safe_allocate(dr_atoms(1:dim/space_dim))
430 safe_allocate(x_old(1:dim))
431 safe_allocate(dr_i(1:dim))
439 dt_max = 10.0_real64 * dt
440 dt_min = 0.02_real64 * dt
443 call gradf(dim, x, en, 1, grad)
449 p_value = dot_product(grad, vel)
452 if (p_value >
m_zero)
then
453 p_times = p_times + 1
455 if (p_times > n_min)
then
456 dt = min(dt * f_inc , dt_max)
457 alpha = alpha * f_alpha
461 n_times = n_times + 1
463 if (n_times > n_decrease_max)
exit
465 if (iter >= n_min)
then
466 dt = max(dt * f_dec, dt_min)
482 if (mix > tol_zero_forces)
then
488 select case (integrator)
489 case (option__gofireintegrator__verlet)
491 vel = vel +
m_half * grad * dt / mass
494 if (p_value >
m_zero)
then
495 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
500 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
504 call gradf(dim, x, en, 1, grad)
508 vel = vel +
m_half * grad * dt / mass
510 case (option__gofireintegrator__euler)
512 if (p_value >
m_zero)
then
513 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
518 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
522 call gradf(dim, x, en, 1, grad)
526 vel = vel + grad * dt / mass
529 case (option__gofireintegrator__semi_implicit_euler)
531 vel = vel + grad * dt / mass
534 if (p_value >
m_zero)
then
535 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
540 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
544 call gradf(dim, x, en, 1, grad)
549 do ia = 0, dim/space_dim - 1
550 grad_atoms(ia+1) = norm2(grad(space_dim*ia+1:space_dim*ia+space_dim))
554 call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
557 if (maxval(abs(grad_atoms(1:))) < tolgrad)
exit
562 safe_deallocate_a(dr_atoms)
563 safe_deallocate_a(x_old)
564 safe_deallocate_a(dr_i)
565 safe_deallocate_a(vel)
566 safe_deallocate_a(grad)
567 safe_deallocate_a(grad_atoms)
569 pop_sub(minimize_fire)
570 end subroutine minimize_fire
575 subroutine prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
576 integer,
intent(in) :: dim
577 integer,
intent(in) :: space_dim
578 real(real64),
intent(inout) :: dr_i(:), dr_atoms(:)
580 integer :: i1, i2, ia
582 real(real64),
parameter :: maxmove = 0.2_real64 *
p_ang
585 do ia = 0, dim/space_dim - 1
587 i2 = space_dim*ia+space_dim
588 dr_atoms(ia+1) = norm2(dr_i(i1:i2))
590 if (dr_atoms(ia+1) > maxmove)
then
591 dr_i(i1:i2) = maxmove * dr_i(i1:i2) / dr_atoms(ia+1)
592 dr_atoms(ia+1) = maxmove
595 end subroutine prevent_large_changes
Returns the euclidean norm of a vector.
real(real64), parameter, public m_zero
real(real64), parameter, public p_ang
real(real64), parameter, public m_half
real(real64), parameter, public m_one
integer, parameter, public minmethod_fr_cg
integer, parameter, public minmethod_bfgs
integer, parameter, public minmethod_nmsimplex
integer, parameter, public minmethod_bfgs2
integer, parameter, public minmethod_nlopt_bobyqa
integer, parameter, public minmethod_nlopt_lbfgs
integer, parameter, public minmethod_pr_cg
integer, parameter, public minmethod_fire
integer, parameter, public minmethod_sd_native
static double f(double w, void *p)