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)
166 ierr = loct_minimize_direct(method, dim, x(1), step, toldr, maxiter, f, write_iter_info, minimum)
169 pop_sub(minimize_multidim_nograd)
171 end subroutine minimize_multidim_nograd
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
183 subroutine f(val, n, x, grad, need_gradient, f_data)
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
193 real(real64),
intent(out) :: minimum
194 real(real64),
intent(in),
optional :: lb(:), ub(:)
195#if defined(HAVE_NLOPT)
198 subroutine nlo_create(opt, alg, n)
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)
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
212 subroutine nlo_set_upper_bounds(ret, opt, upper_bounds)
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
219 subroutine nlo_set_min_objective(ret, opt, f, f_data)
221 integer(c_int),
intent(out) :: ret
222 type(c_ptr),
intent(inout) :: opt
224 subroutine f(val, n, x, grad, need_gradient, f_data)
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
234 type(c_ptr),
intent(in) :: f_data
235 end subroutine nlo_set_min_objective
237 subroutine nlo_set_xtol_abs1(ret, opt, xtol_abs)
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
244 subroutine nlo_set_initial_step1(ret, opt, initial_step1)
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
251 subroutine nlo_set_maxeval(ret, opt, maxeval)
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
258 subroutine nlo_optimize(ret, opt, x, optf)
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
266 subroutine nlo_destroy(opt)
268 type(c_ptr),
intent(inout) :: opt
269 end subroutine nlo_destroy
275 integer,
parameter :: NLOPT_LD_LBFGS = 11
276 integer,
parameter :: NLOPT_LN_BOBYQA = 34
281 call nlo_create(opt, nlopt_ln_bobyqa, dim)
283 call nlo_create(opt, nlopt_ld_lbfgs, dim)
286 if (
present(lb))
then
287 call nlo_set_lower_bounds(ires, opt, lb)
289 if (
present(ub))
then
290 call nlo_set_upper_bounds(ires, opt, ub)
293 call nlo_set_min_objective(ires, opt,
f, c_null_ptr)
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)
301 call nlo_optimize(ires, opt, x, minimum)
303 call nlo_destroy(opt)
308 end subroutine minimize_multidim_nlopt
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
322 procedure(
info_i) :: write_iter_info
323 real(real64),
intent(out) :: minimum
324 integer,
intent(out) :: ierr
326 push_sub(minimize_multidim)
328 assert(ubound(x, dim = 1) >= dim)
332 call minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
335 ierr = loct_minimize(method, dim, x(1), step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum)
339 pop_sub(minimize_multidim)
341 end subroutine minimize_multidim
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
351 procedure(
info_i) :: write_iter_info
352 real(real64),
intent(out) :: minimum
353 integer,
intent(out) :: ierr
356 real(real64),
allocatable :: grad(:)
357 real(real64) :: step2, maxgrad
359 push_sub(minimize_sd)
361 safe_allocate(grad(1:dim))
363 step2 = step*10.0_real64
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)
371 x(1:dim) = x(1:dim) - step2*grad(1:dim)
373 step2 = step2*0.99_real64
378 end subroutine minimize_sd
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
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
408 real(real64),
allocatable :: grad(:)
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(:)
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
421 real(real64),
parameter :: tol_zero_forces = 1e-10_real64
423 push_sub(minimize_fire)
426 assert(mod(dim,space_dim) == 0)
427 assert(
size(x) == dim)
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))
442 dt_max = 10.0_real64 * dt
443 dt_min = 0.02_real64 * dt
446 call gradf(dim, x, en, 1, grad)
452 p_value = dot_product(grad, vel)
455 if (p_value >
m_zero)
then
456 p_times = p_times + 1
458 if (p_times > n_min)
then
459 dt = min(dt * f_inc , dt_max)
460 alpha = alpha * f_alpha
464 n_times = n_times + 1
466 if (n_times > n_decrease_max)
exit
468 if (iter >= n_min)
then
469 dt = max(dt * f_dec, dt_min)
485 if (mix > tol_zero_forces)
then
491 select case (integrator)
492 case (option__gofireintegrator__verlet)
494 vel = vel +
m_half * grad * dt / mass
497 if (p_value >
m_zero)
then
498 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
503 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
507 call gradf(dim, x, en, 1, grad)
511 vel = vel +
m_half * grad * dt / mass
513 case (option__gofireintegrator__euler)
515 if (p_value >
m_zero)
then
516 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
521 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
525 call gradf(dim, x, en, 1, grad)
529 vel = vel + grad * dt / mass
532 case (option__gofireintegrator__semi_implicit_euler)
534 vel = vel + grad * dt / mass
537 if (p_value >
m_zero)
then
538 vel = (
m_one - alpha) * vel + alpha * mix * grad *
lalg_nrm2(dim, vel)
543 call prevent_large_changes(dim, space_dim, dr_i, dr_atoms)
547 call gradf(dim, x, en, 1, grad)
552 do ia = 0, dim/space_dim - 1
553 grad_atoms(ia+1) = norm2(grad(space_dim*ia+1:space_dim*ia+space_dim))
557 call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
560 if (maxval(abs(grad_atoms(1:))) < tolgrad)
exit
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)
572 pop_sub(minimize_fire)
573 end subroutine minimize_fire
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(:)
583 integer :: i1, i2, ia
585 real(real64),
parameter :: maxmove = 0.2_real64 *
p_ang
588 do ia = 0, dim/space_dim - 1
590 i2 = space_dim*ia+space_dim
591 dr_atoms(ia+1) = norm2(dr_i(i1:i2))
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
598 end subroutine prevent_large_changes
Returns the euclidean norm of a vector.
real(real64), parameter, public m_huge
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)