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
 
  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
 
  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)
 
  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
 
  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)
 
  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
 
  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
 
  388  subroutine minimize_fire(dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
 
  389    integer,      
intent(in)    :: dim
 
  390    real(real64), 
intent(inout) :: x(:)
 
  391    real(real64), 
intent(in)    :: step
 
  392    real(real64), 
intent(in)    :: tolgrad
 
  393    integer,      
intent(in)    :: maxiter
 
  395    procedure(
info_i)           :: write_iter_info
 
  396    real(real64), 
intent(out)   :: en
 
  397    integer,      
intent(out)   :: ierr
 
  398    real(real64), 
intent(in)    :: mass(:)
 
  399    integer,      
intent(in)    :: integrator
 
  401    real(real64), 
allocatable :: grad(:), old_grad(:)
 
  403    integer :: p_times, iter, ia
 
  404    real(real64) :: dt, alpha, p_value, dt_max
 
  405    real(real64), 
allocatable :: grad_atoms(:), vel(:), dr_i(:), x_old(:), dr_atoms(:)
 
  407    integer, 
parameter :: n_min = 5
 
  408    real(real64), 
parameter :: f_alpha = 0.99_real64
 
  409    real(real64), 
parameter :: f_inc = 1.1_real64
 
  410    real(real64), 
parameter :: f_dec = 0.5_real64
 
  411    real(real64), 
parameter :: maxmove = 0.2_real64 * 
p_ang 
  412    real(real64), 
parameter :: alpha_start = 0.1_real64
 
  417    assert(mod(dim,3) == 0)
 
  419    safe_allocate(grad_atoms(1:dim/3))
 
  420    safe_allocate(grad(1:dim))
 
  421    safe_allocate(old_grad(1:dim))
 
  422    safe_allocate(vel(1:dim))
 
  423    safe_allocate(dr_atoms(1:dim/3))
 
  424    safe_allocate(x_old(1:dim))
 
  425    safe_allocate(dr_i(1:dim))
 
  434    dt_max = 10.0_real64 * dt
 
  436    call f(dim, x, en, 1, grad)
 
  443      select case (integrator)
 
  445      case (option__gofireintegrator__verlet)
 
  446        dr_i(1:dim) = vel(1:dim) * dt + 
m_half * grad(1:dim) / mass(1:dim) * dt **2
 
  447      case (option__gofireintegrator__euler)
 
  449        dr_i(1:dim) = vel(1:dim)*dt
 
  454        dr_atoms(ia+1) = norm2(dr_i(3*ia+1:3*ia+3))
 
  456        if (dr_atoms(ia+1) > maxmove) 
then 
  457          dr_i(3*ia+1:3*ia+3) = maxmove * dr_i(3*ia+1:3*ia+3) / dr_atoms(ia+1)
 
  458          dr_atoms(ia+1) = maxmove
 
  463      x(1:dim) = x_old(1:dim) + dr_i(1:dim)
 
  466      call f(dim, x, en, 1, grad)
 
  470      select case (integrator)
 
  471      case (option__gofireintegrator__verlet)
 
  473        vel(1:dim) = vel(1:dim) + 
m_half*(grad(1:dim) + old_grad(1:dim))*dt/mass(1:dim)
 
  474      case (option__gofireintegrator__euler)
 
  476        vel(1:dim) = vel(1:dim) + grad(1:dim)*dt/mass(1:dim)
 
  482      p_value = dot_product(grad, vel)
 
  490      if (p_value > 
m_zero) 
then 
  491        p_times = p_times + 1
 
  492        if (p_times > n_min) 
then 
  493          dt = min(dt * f_inc , dt_max)
 
  494          alpha = alpha * f_alpha
 
  506        grad_atoms(ia+1) = norm2(grad(3*ia+1:3*ia+3))
 
  510      call write_iter_info(iter, dim, en, maxval(dr_atoms), maxval(abs(grad_atoms)), x)
 
  513      if (maxval(abs(grad_atoms(1:))) < tolgrad) 
exit 
  515      x_old(1:dim) = x(1:dim)
 
  519    safe_deallocate_a(dr_atoms)
 
  520    safe_deallocate_a(x_old)
 
  521    safe_deallocate_a(dr_i)
 
  522    safe_deallocate_a(vel)
 
  523    safe_deallocate_a(grad)
 
  524    safe_deallocate_a(old_grad)
 
  525    safe_deallocate_a(grad_atoms)
 
Returns the euclidean norm of a vector.
 
subroutine minimize_multidim(method, dim, x, step, line_tol, tolgrad, toldr, maxiter, f, write_iter_info, minimum, ierr)
 
subroutine minimize_sd(dim, x, step, maxiter, f, write_iter_info, minimum, ierr)
 
subroutine minimize_fire(dim, x, step, tolgrad, maxiter, f, write_iter_info, en, ierr, mass, integrator)
Implementation of the Fast Inertial Relaxation Engine (FIRE)
 
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 minmethod_nmsimplex
 
integer minmethod_nlopt_lbfgs
 
integer minmethod_nlopt_bobyqa
 
integer minmethod_sd_native