(git:ed6f26b)
Loading...
Searching...
No Matches
bfgs_optimizer.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines for Geometry optimization using BFGS algorithm
10!> \par History
11!> Module modified by Pierre-André Cazade [pcazade] 01.2020 - University of Limerick.
12!> Modifications enable Space Group Symmetry.
13! **************************************************************************************************
15
19 USE cell_types, ONLY: cell_type, &
20 pbc
26 USE cp_files, ONLY: close_file, &
34 USE cp_fm_types, ONLY: &
39 cp_fm_type, &
45 USE cp_output_handling, ONLY: cp_iterate, &
46 cp_p_file, &
51 USE cp_subsys_types, ONLY: cp_subsys_get, &
53 USE force_env_types, ONLY: force_env_get, &
56 USE gopt_f_methods, ONLY: gopt_f_ii, &
57 gopt_f_io, &
62 USE gopt_f_types, ONLY: gopt_f_type
70 USE kinds, ONLY: default_path_length, &
71 dp
72 USE machine, ONLY: m_flush, &
76 print_spgr, &
80#include "../base/base_uses.f90"
81
82 IMPLICIT NONE
83 PRIVATE
84
85
86! **************************************************************************************************
87!> \brief evaluate the potential energy and its gradients using an array
88!> with same dimension as the particle_set
89!> \param gopt_env the geometry optimization environment
90!> \param x the position where the function should be evaluated
91!> \param f the function value
92!> \param gradient the value of its gradient
93!> \par History
94!> none
95!> \author Teodoro Laino [tlaino] - University of Zurich - 01.2008
96! **************************************************************************************************
97INTERFACE
98
99 SUBROUTINE cp_eval_at(gopt_env, x, f, gradient, master, &
100 final_evaluation, para_env)
101
103 USE gopt_f_types, ONLY: gopt_f_type
104 USE kinds, ONLY: dp
105
106 TYPE(gopt_f_type), POINTER :: gopt_env
107 REAL(KIND=dp), DIMENSION(:), POINTER :: x
108 REAL(KIND=dp), INTENT(out), OPTIONAL :: f
109 REAL(KIND=dp), DIMENSION(:), OPTIONAL, &
110 POINTER :: gradient
111 INTEGER, INTENT(IN) :: master
112 LOGICAL, INTENT(IN), OPTIONAL :: final_evaluation
113 TYPE(mp_para_env_type), POINTER :: para_env
114
115 END SUBROUTINE cp_eval_at
116
117END INTERFACE
118
119 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'bfgs_optimizer'
120 LOGICAL, PARAMETER :: debug_this_module = .true.
121
122 PUBLIC :: geoopt_bfgs
123
124CONTAINS
125
126! **************************************************************************************************
127!> \brief Main driver for BFGS geometry optimizations
128!> \param force_env ...
129!> \param gopt_param ...
130!> \param globenv ...
131!> \param geo_section ...
132!> \param gopt_env ...
133!> \param x0 ...
134!> \par History
135!> 01.2020 modified to perform Space Group Symmetry [pcazade]
136! **************************************************************************************************
137 RECURSIVE SUBROUTINE geoopt_bfgs(force_env, gopt_param, globenv, geo_section, gopt_env, x0)
138
139 TYPE(force_env_type), POINTER :: force_env
140 TYPE(gopt_param_type), POINTER :: gopt_param
141 TYPE(global_environment_type), POINTER :: globenv
142 TYPE(section_vals_type), POINTER :: geo_section
143 TYPE(gopt_f_type), POINTER :: gopt_env
144 REAL(kind=dp), DIMENSION(:), POINTER :: x0
145
146 CHARACTER(len=*), PARAMETER :: routinen = 'geoopt_bfgs'
147 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
148
149 CHARACTER(LEN=5) :: wildcard
150 CHARACTER(LEN=default_path_length) :: hes_filename
151 INTEGER :: handle, hesunit_read, indf, info, &
152 iter_nr, its, maxiter, ndf, nfree, &
153 output_unit
154 LOGICAL :: conv, hesrest, ionode, shell_present, &
155 should_stop, use_mod_hes, use_rfo
156 REAL(kind=dp) :: ediff, emin, eold, etot, pred, rad, rat, &
157 step, t_diff, t_now, t_old
158 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: dg, dr, dx, eigval, gold, work, xold
159 REAL(kind=dp), DIMENSION(:), POINTER :: g
160 TYPE(atomic_kind_list_type), POINTER :: atomic_kinds
161 TYPE(cp_blacs_env_type), POINTER :: blacs_env
162 TYPE(cp_fm_struct_type), POINTER :: fm_struct_hes
163 TYPE(cp_fm_type) :: eigvec_mat, hess_mat, hess_tmp
164 TYPE(cp_logger_type), POINTER :: logger
165 TYPE(mp_para_env_type), POINTER :: para_env
166 TYPE(cp_subsys_type), POINTER :: subsys
167 TYPE(section_vals_type), POINTER :: print_key, root_section
168 TYPE(spgr_type), POINTER :: spgr
169
170 NULLIFY (logger, g, blacs_env, spgr)
171 logger => cp_get_default_logger()
172 para_env => force_env%para_env
173 root_section => force_env%root_section
174 spgr => gopt_env%spgr
175 t_old = m_walltime()
176
177 CALL timeset(routinen, handle)
178 CALL section_vals_val_get(geo_section, "BFGS%TRUST_RADIUS", r_val=rad)
179 print_key => section_vals_get_subs_vals(geo_section, "BFGS%RESTART")
180 ionode = para_env%is_source()
181 maxiter = gopt_param%max_iter
182 conv = .false.
183 rat = 0.0_dp
184 wildcard = " BFGS"
185 hes_filename = ""
186
187 ! Stop if not yet implemented
188 SELECT CASE (gopt_env%type_id)
190 cpabort("BFGS method not yet working with DIMER")
191 END SELECT
192
193 CALL section_vals_val_get(geo_section, "BFGS%USE_RAT_FUN_OPT", l_val=use_rfo)
194 CALL section_vals_val_get(geo_section, "BFGS%USE_MODEL_HESSIAN", l_val=use_mod_hes)
195 CALL section_vals_val_get(geo_section, "BFGS%RESTART_HESSIAN", l_val=hesrest)
196 output_unit = cp_print_key_unit_nr(logger, geo_section, "PRINT%PROGRAM_RUN_INFO", &
197 extension=".geoLog")
198 IF (output_unit > 0) THEN
199 IF (use_rfo) THEN
200 WRITE (unit=output_unit, fmt="(/,T2,A,T78,A3)") &
201 "BFGS| Use rational function optimization for step estimation: ", "YES"
202 ELSE
203 WRITE (unit=output_unit, fmt="(/,T2,A,T78,A3)") &
204 "BFGS| Use rational function optimization for step estimation: ", " NO"
205 END IF
206 IF (use_mod_hes) THEN
207 WRITE (unit=output_unit, fmt="(T2,A,T78,A3)") &
208 "BFGS| Use model Hessian for initial guess: ", "YES"
209 ELSE
210 WRITE (unit=output_unit, fmt="(T2,A,T78,A3)") &
211 "BFGS| Use model Hessian for initial guess: ", " NO"
212 END IF
213 IF (hesrest) THEN
214 WRITE (unit=output_unit, fmt="(T2,A,T78,A3)") &
215 "BFGS| Restart Hessian: ", "YES"
216 ELSE
217 WRITE (unit=output_unit, fmt="(T2,A,T78,A3)") &
218 "BFGS| Restart Hessian: ", " NO"
219 END IF
220 WRITE (unit=output_unit, fmt="(T2,A,T61,F20.3)") &
221 "BFGS| Trust radius: ", rad
222 END IF
223
224 ndf = SIZE(x0)
225 nfree = gopt_env%nfree
226 IF (ndf > 3000) &
227 CALL cp_warn(__location__, &
228 "The dimension of the Hessian matrix ("// &
229 trim(adjustl(cp_to_string(ndf)))//") is greater than 3000. "// &
230 "The diagonalisation of the full Hessian matrix needed for BFGS "// &
231 "is computationally expensive. You should consider to use the linear "// &
232 "scaling variant L-BFGS instead.")
233
234 ! Initialize hessian (hes = unitary matrix or model hessian )
235 CALL cp_blacs_env_create(blacs_env, para_env, globenv%blacs_grid_layout, &
236 globenv%blacs_repeatable)
237 CALL cp_fm_struct_create(fm_struct_hes, para_env=para_env, context=blacs_env, &
238 nrow_global=ndf, ncol_global=ndf)
239 CALL cp_fm_create(hess_mat, fm_struct_hes, name="hess_mat")
240 CALL cp_fm_create(hess_tmp, fm_struct_hes, name="hess_tmp")
241 CALL cp_fm_create(eigvec_mat, fm_struct_hes, name="eigvec_mat")
242 ALLOCATE (eigval(ndf))
243 eigval(:) = 0.0_dp
244
245 CALL force_env_get(force_env=force_env, subsys=subsys)
246 CALL cp_subsys_get(subsys, atomic_kinds=atomic_kinds)
247 CALL get_atomic_kind_set(atomic_kind_set=atomic_kinds%els, shell_present=shell_present)
248 IF (use_mod_hes) THEN
249 IF (shell_present) THEN
250 CALL cp_warn(__location__, &
251 "No model Hessian is available for core-shell models. "// &
252 "A unit matrix is used as the initial Hessian.")
253 use_mod_hes = .false.
254 END IF
255 IF (gopt_env%type_id == default_cell_method_id) THEN
256 CALL cp_warn(__location__, &
257 "No model Hessian is available for cell optimizations. "// &
258 "A unit matrix is used as the initial Hessian.")
259 use_mod_hes = .false.
260 END IF
261 END IF
262
263 IF (use_mod_hes) THEN
264 CALL cp_fm_set_all(hess_mat, alpha=zero, beta=0.00_dp)
265 CALL construct_initial_hess(gopt_env%force_env, hess_mat)
266 CALL cp_fm_to_fm(hess_mat, hess_tmp)
267 CALL choose_eigv_solver(hess_tmp, eigvec_mat, eigval, info=info)
268 ! In rare cases the diagonalization of hess_mat fails (bug in scalapack?)
269 IF (info /= 0) THEN
270 CALL cp_fm_set_all(hess_mat, alpha=zero, beta=one)
271 IF (output_unit > 0) WRITE (output_unit, *) &
272 "BFGS: Matrix diagonalization failed, using unity as model Hessian."
273 ELSE
274 DO its = 1, SIZE(eigval)
275 IF (eigval(its) .LT. 0.1_dp) eigval(its) = 0.1_dp
276 END DO
277 CALL cp_fm_to_fm(eigvec_mat, hess_tmp)
278 CALL cp_fm_column_scale(eigvec_mat, eigval)
279 CALL parallel_gemm("N", "T", ndf, ndf, ndf, one, hess_tmp, eigvec_mat, zero, hess_mat)
280 END IF
281 ELSE
282 CALL cp_fm_set_all(hess_mat, alpha=zero, beta=one)
283 END IF
284
285 ALLOCATE (xold(ndf))
286 xold(:) = x0(:)
287
288 ALLOCATE (g(ndf))
289 g(:) = 0.0_dp
290
291 ALLOCATE (gold(ndf))
292 gold(:) = 0.0_dp
293
294 ALLOCATE (dx(ndf))
295 dx(:) = 0.0_dp
296
297 ALLOCATE (dg(ndf))
298 dg(:) = 0.0_dp
299
300 ALLOCATE (work(ndf))
301 work(:) = 0.0_dp
302
303 ALLOCATE (dr(ndf))
304 dr(:) = 0.0_dp
305
306 ! find space_group
307 CALL section_vals_val_get(geo_section, "KEEP_SPACE_GROUP", l_val=spgr%keep_space_group)
308 IF (spgr%keep_space_group) THEN
309 CALL identify_space_group(subsys, geo_section, gopt_env, output_unit)
310 CALL spgr_apply_rotations_coord(spgr, x0)
311 CALL print_spgr(spgr)
312 END IF
313
314 ! Geometry optimization starts now
315 CALL cp_iterate(logger%iter_info, increment=0, iter_nr_out=iter_nr)
316 CALL print_geo_opt_header(gopt_env, output_unit, wildcard)
317
318 ! Calculate Energy & Gradients
319 CALL cp_eval_at(gopt_env, x0, etot, g, gopt_env%force_env%para_env%mepos, &
320 .false., gopt_env%force_env%para_env)
321
322 ! Symmetrize coordinates and forces
323 IF (spgr%keep_space_group) THEN
324 CALL spgr_apply_rotations_coord(spgr, x0)
325 CALL spgr_apply_rotations_force(spgr, g)
326 END IF
327
328 ! Print info at time 0
329 emin = etot
330 t_now = m_walltime()
331 t_diff = t_now - t_old
332 t_old = t_now
333 CALL gopt_f_io_init(gopt_env, output_unit, etot, wildcard=wildcard, its=iter_nr, used_time=t_diff)
334 DO its = iter_nr + 1, maxiter
335 CALL cp_iterate(logger%iter_info, last=(its == maxiter))
336 CALL section_vals_val_set(geo_section, "STEP_START_VAL", i_val=its)
337 CALL gopt_f_ii(its, output_unit)
338
339 ! Hessian update/restarting
340 IF (((its - iter_nr) == 1) .AND. hesrest) THEN
341 IF (ionode) THEN
342 CALL section_vals_val_get(geo_section, "BFGS%RESTART_FILE_NAME", c_val=hes_filename)
343 IF (len_trim(hes_filename) == 0) THEN
344 ! Set default Hessian restart file name if no file name is defined
345 hes_filename = trim(logger%iter_info%project_name)//"-BFGS.Hessian"
346 END IF
347 IF (output_unit > 0) THEN
348 WRITE (unit=output_unit, fmt="(/,T2,A)") &
349 "BFGS| Checking for Hessian restart file <"//trim(adjustl(hes_filename))//">"
350 END IF
351 CALL open_file(file_name=trim(hes_filename), file_status="OLD", &
352 file_form="UNFORMATTED", file_action="READ", unit_number=hesunit_read)
353 IF (output_unit > 0) THEN
354 WRITE (unit=output_unit, fmt="(T2,A)") &
355 "BFGS| Hessian restart file read"
356 END IF
357 END IF
358 CALL cp_fm_read_unformatted(hess_mat, hesunit_read)
359 IF (ionode) CALL close_file(unit_number=hesunit_read)
360 ELSE
361 IF ((its - iter_nr) > 1) THEN
362 ! Symmetrize old coordinates and old forces
363 IF (spgr%keep_space_group) THEN
364 CALL spgr_apply_rotations_coord(spgr, xold)
365 CALL spgr_apply_rotations_force(spgr, gold)
366 END IF
367
368 DO indf = 1, ndf
369 dx(indf) = x0(indf) - xold(indf)
370 dg(indf) = g(indf) - gold(indf)
371 END DO
372
373 CALL bfgs(ndf, dx, dg, hess_mat, work, para_env)
374
375 ! Symmetrize coordinates and forces change
376 IF (spgr%keep_space_group) THEN
377 CALL spgr_apply_rotations_force(spgr, dx)
378 CALL spgr_apply_rotations_force(spgr, dg)
379 END IF
380
381 !Possibly dump the Hessian file
382 IF (btest(cp_print_key_should_output(logger%iter_info, print_key), cp_p_file)) THEN
383 CALL write_bfgs_hessian(geo_section, hess_mat, logger)
384 END IF
385 END IF
386 END IF
387
388 ! Symmetrize coordinates and forces
389 IF (spgr%keep_space_group) THEN
390 CALL spgr_apply_rotations_coord(spgr, x0)
391 CALL spgr_apply_rotations_force(spgr, g)
392 END IF
393
394 ! Setting the present positions & gradients as old
395 xold(:) = x0
396 gold(:) = g
397
398 ! Copying hessian hes to (ndf x ndf) matrix hes_mat for diagonalization
399 CALL cp_fm_to_fm(hess_mat, hess_tmp)
400
401 CALL choose_eigv_solver(hess_tmp, eigvec_mat, eigval, info=info)
402
403 ! In rare cases the diagonalization of hess_mat fails (bug in scalapack?)
404 IF (info /= 0) THEN
405 IF (output_unit > 0) WRITE (output_unit, *) &
406 "BFGS: Matrix diagonalization failed, resetting Hessian to unity."
407 CALL cp_fm_set_all(hess_mat, alpha=zero, beta=one)
408 CALL cp_fm_to_fm(hess_mat, hess_tmp)
409 CALL choose_eigv_solver(hess_tmp, eigvec_mat, eigval)
410 END IF
411
412 IF (use_rfo) THEN
413 CALL set_hes_eig(ndf, eigval, work)
414 dx(:) = eigval
415 CALL rat_fun_opt(ndf, dg, eigval, work, eigvec_mat, g, para_env)
416 END IF
417 CALL geoopt_get_step(ndf, eigval, eigvec_mat, hess_tmp, dr, g, para_env, use_rfo)
418
419 ! Symmetrize dr
420 IF (spgr%keep_space_group) THEN
421 CALL spgr_apply_rotations_force(spgr, dr)
422 END IF
423
424 CALL trust_radius(ndf, step, rad, rat, dr, output_unit)
425
426 ! Update the atomic positions
427 x0 = x0 + dr
428
429 ! Symmetrize coordinates
430 IF (spgr%keep_space_group) THEN
431 CALL spgr_apply_rotations_coord(spgr, x0)
432 END IF
433
434 CALL energy_predict(ndf, work, hess_mat, dr, g, conv, pred, para_env)
435 eold = etot
436
437 ! Energy & Gradients at new step
438 CALL cp_eval_at(gopt_env, x0, etot, g, gopt_env%force_env%para_env%mepos, &
439 .false., gopt_env%force_env%para_env)
440
441 ediff = etot - eold
442
443 ! Symmetrize forces
444 IF (spgr%keep_space_group) THEN
445 CALL spgr_apply_rotations_force(spgr, g)
446 END IF
447
448 ! check for an external exit command
449 CALL external_control(should_stop, "GEO", globenv=globenv)
450 IF (should_stop) EXIT
451
452 ! Some IO and Convergence check
453 t_now = m_walltime()
454 t_diff = t_now - t_old
455 t_old = t_now
456 CALL gopt_f_io(gopt_env, force_env, root_section, its, etot, output_unit, &
457 eold, emin, wildcard, gopt_param, ndf, dr, g, conv, pred, rat, &
458 step, rad, used_time=t_diff)
459
460 IF (conv .OR. (its == maxiter)) EXIT
461 IF (etot < emin) emin = etot
462 IF (use_rfo) CALL update_trust_rad(rat, rad, step, ediff)
463 END DO
464
465 IF (its == maxiter .AND. (.NOT. conv)) THEN
466 CALL print_geo_opt_nc(gopt_env, output_unit)
467 END IF
468
469 ! Write final information, if converged
470 CALL cp_iterate(logger%iter_info, last=.true., increment=0)
471 CALL write_bfgs_hessian(geo_section, hess_mat, logger)
472 CALL gopt_f_io_finalize(gopt_env, force_env, x0, conv, its, root_section, &
473 gopt_env%force_env%para_env, gopt_env%force_env%para_env%mepos, output_unit)
474
475 CALL cp_fm_struct_release(fm_struct_hes)
476 CALL cp_fm_release(hess_mat)
477 CALL cp_fm_release(eigvec_mat)
478 CALL cp_fm_release(hess_tmp)
479
480 CALL cp_blacs_env_release(blacs_env)
481 DEALLOCATE (xold)
482 DEALLOCATE (g)
483 DEALLOCATE (gold)
484 DEALLOCATE (dx)
485 DEALLOCATE (dg)
486 DEALLOCATE (eigval)
487 DEALLOCATE (work)
488 DEALLOCATE (dr)
489
490 CALL cp_print_key_finished_output(output_unit, logger, geo_section, &
491 "PRINT%PROGRAM_RUN_INFO")
492 CALL timestop(handle)
493
494 END SUBROUTINE geoopt_bfgs
495
496! **************************************************************************************************
497!> \brief ...
498!> \param ndf ...
499!> \param dg ...
500!> \param eigval ...
501!> \param work ...
502!> \param eigvec_mat ...
503!> \param g ...
504!> \param para_env ...
505! **************************************************************************************************
506 SUBROUTINE rat_fun_opt(ndf, dg, eigval, work, eigvec_mat, g, para_env)
507
508 INTEGER, INTENT(IN) :: ndf
509 REAL(kind=dp), INTENT(INOUT) :: dg(ndf), eigval(ndf), work(ndf)
510 TYPE(cp_fm_type), INTENT(IN) :: eigvec_mat
511 REAL(kind=dp), INTENT(INOUT) :: g(ndf)
512 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
513
514 CHARACTER(LEN=*), PARAMETER :: routinen = 'rat_fun_opt'
515 REAL(kind=dp), PARAMETER :: one = 1.0_dp
516
517 INTEGER :: handle, i, indf, iref, iter, j, k, l, &
518 maxit, ncol_local, nrow_local
519 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
520 LOGICAL :: bisec, fail, set
521 REAL(kind=dp) :: fun, fun1, fun2, fun3, fung, lam1, lam2, &
522 ln, lp, ssize, step, stol
523 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), POINTER :: local_data
524
525 CALL timeset(routinen, handle)
526
527 stol = 1.0e-8_dp
528 ssize = 0.2_dp
529 maxit = 999
530 fail = .false.
531 bisec = .false.
532
533 dg = 0._dp
534
535 CALL cp_fm_get_info(eigvec_mat, row_indices=row_indices, col_indices=col_indices, &
536 local_data=local_data, nrow_local=nrow_local, ncol_local=ncol_local)
537
538 DO i = 1, nrow_local
539 j = row_indices(i)
540 DO k = 1, ncol_local
541 l = col_indices(k)
542 dg(l) = dg(l) + local_data(i, k)*g(j)
543 END DO
544 END DO
545 CALL para_env%sum(dg)
546
547 set = .false.
548
549 DO
550
551! calculating Lambda
552
553 lp = 0.0_dp
554 iref = 1
555 ln = 0.0_dp
556 IF (eigval(iref) < 0.0_dp) ln = eigval(iref) - 0.01_dp
557
558 iter = 0
559 DO
560 iter = iter + 1
561 fun = 0.0_dp
562 fung = 0.0_dp
563 DO indf = 1, ndf
564 fun = fun + dg(indf)**2/(ln - eigval(indf))
565 fung = fung - dg(indf)**2/(ln - eigval(indf)**2)
566 END DO
567 fun = fun - ln
568 fung = fung - one
569 step = fun/fung
570 ln = ln - step
571 IF (abs(step) < stol) GOTO 200
572 IF (iter >= maxit) EXIT
573 END DO
574100 CONTINUE
575 bisec = .true.
576 iter = 0
577 maxit = 9999
578 lam1 = 0.0_dp
579 IF (eigval(iref) < 0.0_dp) lam1 = eigval(iref) - 0.01_dp
580 fun1 = 0.0_dp
581 DO indf = 1, ndf
582 fun1 = fun1 + dg(indf)**2/(lam1 - eigval(indf))
583 END DO
584 fun1 = fun1 - lam1
585 step = abs(lam1)/1000.0_dp
586 IF (step < ssize) step = ssize
587 DO
588 iter = iter + 1
589 IF (iter > maxit) THEN
590 ln = 0.0_dp
591 lp = 0.0_dp
592 fail = .true.
593 GOTO 300
594 END IF
595 fun2 = 0.0_dp
596 lam2 = lam1 - iter*step
597 DO indf = 1, ndf
598 fun2 = fun2 + eigval(indf)**2/(lam2 - eigval(indf))
599 END DO
600 fun2 = fun2 - lam2
601 IF (fun2*fun1 < 0.0_dp) THEN
602 iter = 0
603 DO
604 iter = iter + 1
605 IF (iter > maxit) THEN
606 ln = 0.0_dp
607 lp = 0.0_dp
608 fail = .true.
609 GOTO 300
610 END IF
611 step = (lam1 + lam2)/2
612 fun3 = 0.0_dp
613 DO indf = 1, ndf
614 fun3 = fun3 + dg(indf)**2/(step - eigval(indf))
615 END DO
616 fun3 = fun3 - step
617
618 IF (abs(step - lam2) < stol) THEN
619 ln = step
620 GOTO 200
621 END IF
622
623 IF (fun3*fun1 < stol) THEN
624 lam2 = step
625 ELSE
626 lam1 = step
627 END IF
628 END DO
629 END IF
630 END DO
631
632200 CONTINUE
633 IF ((ln > eigval(iref)) .OR. ((ln > 0.0_dp) .AND. &
634 (eigval(iref) > 0.0_dp))) THEN
635
636 IF (.NOT. bisec) GOTO 100
637 ln = 0.0_dp
638 lp = 0.0_dp
639 fail = .true.
640 END IF
641
642300 CONTINUE
643
644 IF (fail .AND. .NOT. set) THEN
645 set = .true.
646 DO indf = 1, ndf
647 eigval(indf) = eigval(indf)*work(indf)
648 END DO
649 cycle
650 END IF
651
652 IF (.NOT. set) THEN
653 work(1:ndf) = one
654 END IF
655
656 DO indf = 1, ndf
657 eigval(indf) = eigval(indf) - ln
658 END DO
659 EXIT
660 END DO
661
662 CALL timestop(handle)
663
664 END SUBROUTINE rat_fun_opt
665
666! **************************************************************************************************
667!> \brief ...
668!> \param ndf ...
669!> \param dx ...
670!> \param dg ...
671!> \param hess_mat ...
672!> \param work ...
673!> \param para_env ...
674! **************************************************************************************************
675 SUBROUTINE bfgs(ndf, dx, dg, hess_mat, work, para_env)
676 INTEGER, INTENT(IN) :: ndf
677 REAL(kind=dp), INTENT(INOUT) :: dx(ndf), dg(ndf)
678 TYPE(cp_fm_type), INTENT(IN) :: hess_mat
679 REAL(kind=dp), INTENT(INOUT) :: work(ndf)
680 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
681
682 CHARACTER(LEN=*), PARAMETER :: routinen = 'bfgs'
683 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
684
685 INTEGER :: handle, i, j, k, l, ncol_local, &
686 nrow_local
687 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
688 REAL(kind=dp) :: ddot, dxw, gdx
689 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), POINTER :: local_hes
690
691 CALL timeset(routinen, handle)
692
693 CALL cp_fm_get_info(hess_mat, row_indices=row_indices, col_indices=col_indices, &
694 local_data=local_hes, nrow_local=nrow_local, ncol_local=ncol_local)
695
696 work = zero
697 DO i = 1, nrow_local
698 j = row_indices(i)
699 DO k = 1, ncol_local
700 l = col_indices(k)
701 work(j) = work(j) + local_hes(i, k)*dx(l)
702 END DO
703 END DO
704
705 CALL para_env%sum(work)
706
707 gdx = ddot(ndf, dg, 1, dx, 1)
708 gdx = one/gdx
709 dxw = ddot(ndf, dx, 1, work, 1)
710 dxw = one/dxw
711
712 DO i = 1, nrow_local
713 j = row_indices(i)
714 DO k = 1, ncol_local
715 l = col_indices(k)
716 local_hes(i, k) = local_hes(i, k) + gdx*dg(j)*dg(l) - &
717 dxw*work(j)*work(l)
718 END DO
719 END DO
720
721 CALL timestop(handle)
722
723 END SUBROUTINE bfgs
724
725! **************************************************************************************************
726!> \brief ...
727!> \param ndf ...
728!> \param eigval ...
729!> \param work ...
730! **************************************************************************************************
731 SUBROUTINE set_hes_eig(ndf, eigval, work)
732 INTEGER, INTENT(IN) :: ndf
733 REAL(kind=dp), INTENT(INOUT) :: eigval(ndf), work(ndf)
734
735 CHARACTER(LEN=*), PARAMETER :: routinen = 'set_hes_eig'
736 REAL(kind=dp), PARAMETER :: max_neg = -0.5_dp, max_pos = 5.0_dp, &
737 min_eig = 0.005_dp, one = 1.0_dp
738
739 INTEGER :: handle, indf
740 LOGICAL :: neg
741
742 CALL timeset(routinen, handle)
743
744 DO indf = 1, ndf
745 IF (eigval(indf) < 0.0_dp) neg = .true.
746 IF (eigval(indf) > 1000.0_dp) eigval(indf) = 1000.0_dp
747 END DO
748 DO indf = 1, ndf
749 IF (eigval(indf) < 0.0_dp) THEN
750 IF (eigval(indf) < max_neg) THEN
751 eigval(indf) = max_neg
752 ELSE IF (eigval(indf) > -min_eig) THEN
753 eigval(indf) = -min_eig
754 END IF
755 ELSE IF (eigval(indf) < 1000.0_dp) THEN
756 IF (eigval(indf) < min_eig) THEN
757 eigval(indf) = min_eig
758 ELSE IF (eigval(indf) > max_pos) THEN
759 eigval(indf) = max_pos
760 END IF
761 END IF
762 END DO
763
764 DO indf = 1, ndf
765 IF (eigval(indf) < 0.0_dp) THEN
766 work(indf) = -one
767 ELSE
768 work(indf) = one
769 END IF
770 END DO
771
772 CALL timestop(handle)
773
774 END SUBROUTINE set_hes_eig
775
776! **************************************************************************************************
777!> \brief ...
778!> \param ndf ...
779!> \param eigval ...
780!> \param eigvec_mat ...
781!> \param hess_tmp ...
782!> \param dr ...
783!> \param g ...
784!> \param para_env ...
785!> \param use_rfo ...
786! **************************************************************************************************
787 SUBROUTINE geoopt_get_step(ndf, eigval, eigvec_mat, hess_tmp, dr, g, para_env, use_rfo)
788
789 INTEGER, INTENT(IN) :: ndf
790 REAL(kind=dp), INTENT(INOUT) :: eigval(ndf)
791 TYPE(cp_fm_type), INTENT(IN) :: eigvec_mat, hess_tmp
792 REAL(kind=dp), INTENT(INOUT) :: dr(ndf), g(ndf)
793 TYPE(mp_para_env_type), OPTIONAL, POINTER :: para_env
794 LOGICAL :: use_rfo
795
796 REAL(kind=dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
797
798 INTEGER :: i, indf, j, k, l, ncol_local, nrow_local
799 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
800 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), POINTER :: local_data
801 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
802 TYPE(cp_fm_type) :: tmp
803
804 CALL cp_fm_to_fm(eigvec_mat, hess_tmp)
805 IF (use_rfo) THEN
806 DO indf = 1, ndf
807 eigval(indf) = one/eigval(indf)
808 END DO
809 ELSE
810 DO indf = 1, ndf
811 eigval(indf) = one/max(0.0001_dp, eigval(indf))
812 END DO
813 END IF
814
815 CALL cp_fm_column_scale(hess_tmp, eigval)
816 CALL cp_fm_get_info(eigvec_mat, matrix_struct=matrix_struct)
817 CALL cp_fm_create(tmp, matrix_struct, name="tmp")
818
819 CALL parallel_gemm("N", "T", ndf, ndf, ndf, one, hess_tmp, eigvec_mat, zero, tmp)
820
821 CALL cp_fm_transpose(tmp, hess_tmp)
822 CALL cp_fm_release(tmp)
823
824! ** New step **
825
826 CALL cp_fm_get_info(hess_tmp, row_indices=row_indices, col_indices=col_indices, &
827 local_data=local_data, nrow_local=nrow_local, ncol_local=ncol_local)
828
829 dr = 0.0_dp
830 DO i = 1, nrow_local
831 j = row_indices(i)
832 DO k = 1, ncol_local
833 l = col_indices(k)
834 dr(j) = dr(j) - local_data(i, k)*g(l)
835 END DO
836 END DO
837
838 CALL para_env%sum(dr)
839
840 END SUBROUTINE geoopt_get_step
841
842! **************************************************************************************************
843!> \brief ...
844!> \param ndf ...
845!> \param step ...
846!> \param rad ...
847!> \param rat ...
848!> \param dr ...
849!> \param output_unit ...
850! **************************************************************************************************
851 SUBROUTINE trust_radius(ndf, step, rad, rat, dr, output_unit)
852 INTEGER, INTENT(IN) :: ndf
853 REAL(kind=dp), INTENT(INOUT) :: step, rad, rat, dr(ndf)
854 INTEGER, INTENT(IN) :: output_unit
855
856 CHARACTER(LEN=*), PARAMETER :: routinen = 'trust_radius'
857 REAL(kind=dp), PARAMETER :: one = 1.0_dp
858
859 INTEGER :: handle
860 REAL(kind=dp) :: scal
861
862 CALL timeset(routinen, handle)
863
864 step = maxval(abs(dr))
865 scal = max(one, rad/step)
866
867 IF (step > rad) THEN
868 rat = rad/step
869 CALL dscal(ndf, rat, dr, 1)
870 step = rad
871 IF (output_unit > 0) THEN
872 WRITE (unit=output_unit, fmt="(/,T2,A,F8.5)") &
873 " Step is scaled; Scaling factor = ", rat
874 CALL m_flush(output_unit)
875 END IF
876 END IF
877 CALL timestop(handle)
878
879 END SUBROUTINE trust_radius
880
881! **************************************************************************************************
882!> \brief ...
883!> \param ndf ...
884!> \param work ...
885!> \param hess_mat ...
886!> \param dr ...
887!> \param g ...
888!> \param conv ...
889!> \param pred ...
890!> \param para_env ...
891! **************************************************************************************************
892 SUBROUTINE energy_predict(ndf, work, hess_mat, dr, g, conv, pred, para_env)
893
894 INTEGER, INTENT(IN) :: ndf
895 REAL(kind=dp), INTENT(INOUT) :: work(ndf)
896 TYPE(cp_fm_type), INTENT(IN) :: hess_mat
897 REAL(kind=dp), INTENT(INOUT) :: dr(ndf), g(ndf)
898 LOGICAL, INTENT(INOUT) :: conv
899 REAL(kind=dp), INTENT(INOUT) :: pred
900 TYPE(mp_para_env_type), POINTER :: para_env
901
902 CHARACTER(LEN=*), PARAMETER :: routinen = 'energy_predict'
903 REAL(kind=dp), PARAMETER :: zero = 0.0_dp
904
905 INTEGER :: handle, i, j, k, l, ncol_local, &
906 nrow_local
907 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
908 REAL(kind=dp) :: ddot, ener1, ener2
909 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), POINTER :: local_data
910
911 CALL timeset(routinen, handle)
912
913 ener1 = ddot(ndf, g, 1, dr, 1)
914
915 CALL cp_fm_get_info(hess_mat, row_indices=row_indices, col_indices=col_indices, &
916 local_data=local_data, nrow_local=nrow_local, ncol_local=ncol_local)
917
918 work = zero
919 DO i = 1, nrow_local
920 j = row_indices(i)
921 DO k = 1, ncol_local
922 l = col_indices(k)
923 work(j) = work(j) + local_data(i, k)*dr(l)
924 END DO
925 END DO
926
927 CALL para_env%sum(work)
928 ener2 = ddot(ndf, dr, 1, work, 1)
929 pred = ener1 + 0.5_dp*ener2
930 conv = .false.
931 CALL timestop(handle)
932
933 END SUBROUTINE energy_predict
934
935! **************************************************************************************************
936!> \brief ...
937!> \param rat ...
938!> \param rad ...
939!> \param step ...
940!> \param ediff ...
941! **************************************************************************************************
942 SUBROUTINE update_trust_rad(rat, rad, step, ediff)
943
944 REAL(kind=dp), INTENT(INOUT) :: rat, rad, step, ediff
945
946 CHARACTER(LEN=*), PARAMETER :: routinen = 'update_trust_rad'
947 REAL(kind=dp), PARAMETER :: max_trust = 1.0_dp, min_trust = 0.1_dp
948
949 INTEGER :: handle
950
951 CALL timeset(routinen, handle)
952
953 IF (rat > 4.0_dp) THEN
954 IF (ediff < 0.0_dp) THEN
955 rad = step*0.5_dp
956 ELSE
957 rad = step*0.25_dp
958 END IF
959 ELSE IF (rat > 2.0_dp) THEN
960 IF (ediff < 0.0_dp) THEN
961 rad = step*0.75_dp
962 ELSE
963 rad = step*0.5_dp
964 END IF
965 ELSE IF (rat > 4.0_dp/3.0_dp) THEN
966 IF (ediff < 0.0_dp) THEN
967 rad = step
968 ELSE
969 rad = step*0.75_dp
970 END IF
971 ELSE IF (rat > 10.0_dp/9.0_dp) THEN
972 IF (ediff < 0.0_dp) THEN
973 rad = step*1.25_dp
974 ELSE
975 rad = step
976 END IF
977 ELSE IF (rat > 0.9_dp) THEN
978 IF (ediff < 0.0_dp) THEN
979 rad = step*1.5_dp
980 ELSE
981 rad = step*1.25_dp
982 END IF
983 ELSE IF (rat > 0.75_dp) THEN
984 IF (ediff < 0.0_dp) THEN
985 rad = step*1.25_dp
986 ELSE
987 rad = step
988 END IF
989 ELSE IF (rat > 0.5_dp) THEN
990 IF (ediff < 0.0_dp) THEN
991 rad = step
992 ELSE
993 rad = step*0.75_dp
994 END IF
995 ELSE IF (rat > 0.25_dp) THEN
996 IF (ediff < 0.0_dp) THEN
997 rad = step*0.75_dp
998 ELSE
999 rad = step*0.5_dp
1000 END IF
1001 ELSE IF (ediff < 0.0_dp) THEN
1002 rad = step*0.5_dp
1003 ELSE
1004 rad = step*0.25_dp
1005 END IF
1006
1007 rad = max(rad, min_trust)
1008 rad = min(rad, max_trust)
1009 CALL timestop(handle)
1010
1011 END SUBROUTINE update_trust_rad
1012
1013! **************************************************************************************************
1014
1015! **************************************************************************************************
1016!> \brief ...
1017!> \param geo_section ...
1018!> \param hess_mat ...
1019!> \param logger ...
1020! **************************************************************************************************
1021 SUBROUTINE write_bfgs_hessian(geo_section, hess_mat, logger)
1022
1023 TYPE(section_vals_type), POINTER :: geo_section
1024 TYPE(cp_fm_type), INTENT(IN) :: hess_mat
1025 TYPE(cp_logger_type), POINTER :: logger
1026
1027 CHARACTER(LEN=*), PARAMETER :: routinen = 'write_bfgs_hessian'
1028
1029 INTEGER :: handle, hesunit
1030
1031 CALL timeset(routinen, handle)
1032
1033 hesunit = cp_print_key_unit_nr(logger, geo_section, "BFGS%RESTART", &
1034 extension=".Hessian", file_form="UNFORMATTED", file_action="WRITE", &
1035 file_position="REWIND")
1036
1037 CALL cp_fm_write_unformatted(hess_mat, hesunit)
1038
1039 CALL cp_print_key_finished_output(hesunit, logger, geo_section, "BFGS%RESTART")
1040
1041 CALL timestop(handle)
1042
1043 END SUBROUTINE write_bfgs_hessian
1044! **************************************************************************************************
1045
1046! **************************************************************************************************
1047!> \brief ...
1048!> \param force_env ...
1049!> \param hess_mat ...
1050! **************************************************************************************************
1051 SUBROUTINE construct_initial_hess(force_env, hess_mat)
1052
1053 TYPE(force_env_type), POINTER :: force_env
1054 TYPE(cp_fm_type), INTENT(IN) :: hess_mat
1055
1056 INTEGER :: i, iat_col, iat_row, iglobal, iind, j, &
1057 jat_row, jglobal, jind, k, natom, &
1058 ncol_local, nrow_local, z
1059 INTEGER, ALLOCATABLE, DIMENSION(:) :: at_row
1060 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
1061 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: d_ij, rho_ij
1062 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: r_ij
1063 REAL(kind=dp), DIMENSION(3, 3) :: alpha, r0
1064 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), POINTER :: fixed, local_data
1065 TYPE(cell_type), POINTER :: cell
1066 TYPE(cp_subsys_type), POINTER :: subsys
1067 TYPE(particle_list_type), POINTER :: particles
1068
1069 CALL force_env_get(force_env=force_env, subsys=subsys, cell=cell)
1070 CALL cp_subsys_get(subsys, &
1071 particles=particles)
1072
1073 alpha(1, :) = (/1._dp, 0.3949_dp, 0.3949_dp/)
1074 alpha(2, :) = (/0.3494_dp, 0.2800_dp, 0.2800_dp/)
1075 alpha(3, :) = (/0.3494_dp, 0.2800_dp, 0.1800_dp/)
1076
1077 r0(1, :) = (/1.35_dp, 2.10_dp, 2.53_dp/)
1078 r0(2, :) = (/2.10_dp, 2.87_dp, 3.40_dp/)
1079 r0(3, :) = (/2.53_dp, 3.40_dp, 3.40_dp/)
1080
1081 CALL cp_fm_get_info(hess_mat, row_indices=row_indices, col_indices=col_indices, &
1082 local_data=local_data, nrow_local=nrow_local, ncol_local=ncol_local)
1083 natom = particles%n_els
1084 ALLOCATE (at_row(natom))
1085 ALLOCATE (rho_ij(natom, natom))
1086 ALLOCATE (d_ij(natom, natom))
1087 ALLOCATE (r_ij(natom, natom, 3))
1088 ALLOCATE (fixed(3, natom))
1089 fixed = 1.0_dp
1090 CALL fix_atom_control(force_env, fixed)
1091 DO i = 1, 3
1092 CALL hess_mat%matrix_struct%para_env%min(fixed(i, :))
1093 END DO
1094 rho_ij = 0
1095 !XXXX insert proper rows !XXX
1096 at_row = 3
1097 DO i = 1, natom
1098 CALL get_atomic_kind(atomic_kind=particles%els(i)%atomic_kind, z=z)
1099 IF (z .LE. 10) at_row(i) = 2
1100 IF (z .LE. 2) at_row(i) = 1
1101 END DO
1102 DO i = 2, natom
1103 iat_row = at_row(i)
1104 DO j = 1, i - 1
1105 jat_row = at_row(j)
1106 !pbc for a distance vector
1107 r_ij(j, i, :) = pbc(particles%els(i)%r, particles%els(j)%r, cell)
1108 r_ij(i, j, :) = -r_ij(j, i, :)
1109 d_ij(j, i) = sqrt(dot_product(r_ij(j, i, :), r_ij(j, i, :)))
1110 d_ij(i, j) = d_ij(j, i)
1111 rho_ij(j, i) = exp(alpha(jat_row, iat_row)*(r0(jat_row, iat_row)**2 - d_ij(j, i)**2))
1112 rho_ij(i, j) = rho_ij(j, i)
1113 END DO
1114 END DO
1115 DO i = 1, ncol_local
1116 iglobal = col_indices(i)
1117 iind = mod(iglobal - 1, 3) + 1
1118 iat_col = (iglobal + 2)/3
1119 IF (iat_col .GT. natom) cycle
1120 DO j = 1, nrow_local
1121 jglobal = row_indices(j)
1122 jind = mod(jglobal - 1, 3) + 1
1123 iat_row = (jglobal + 2)/3
1124 IF (iat_row .GT. natom) cycle
1125 IF (iat_row .NE. iat_col) THEN
1126 IF (d_ij(iat_row, iat_col) .LT. 6.0_dp) &
1127 local_data(j, i) = local_data(j, i) + &
1128 angle_second_deriv(r_ij, d_ij, rho_ij, iind, jind, iat_col, iat_row, natom)
1129 ELSE
1130 local_data(j, i) = local_data(j, i) + &
1131 angle_second_deriv(r_ij, d_ij, rho_ij, iind, jind, iat_col, iat_row, natom)
1132 END IF
1133 IF (iat_col .NE. iat_row) THEN
1134 IF (d_ij(iat_row, iat_col) .LT. 6.0_dp) &
1135 local_data(j, i) = local_data(j, i) - &
1136 dist_second_deriv(r_ij(iat_col, iat_row, :), &
1137 iind, jind, d_ij(iat_row, iat_col), rho_ij(iat_row, iat_col))
1138 ELSE
1139 DO k = 1, natom
1140 IF (k == iat_col) cycle
1141 IF (d_ij(iat_row, k) .LT. 6.0_dp) &
1142 local_data(j, i) = local_data(j, i) + &
1143 dist_second_deriv(r_ij(iat_col, k, :), &
1144 iind, jind, d_ij(iat_row, k), rho_ij(iat_row, k))
1145 END DO
1146 END IF
1147 IF (fixed(jind, iat_row) .LT. 0.5_dp .OR. fixed(iind, iat_col) .LT. 0.5_dp) THEN
1148 local_data(j, i) = 0.0_dp
1149 IF (jind == iind .AND. iat_row == iat_col) local_data(j, i) = 1.0_dp
1150 END IF
1151 END DO
1152 END DO
1153 DEALLOCATE (fixed)
1154 DEALLOCATE (rho_ij)
1155 DEALLOCATE (d_ij)
1156 DEALLOCATE (r_ij)
1157 DEALLOCATE (at_row)
1158
1159 END SUBROUTINE construct_initial_hess
1160
1161! **************************************************************************************************
1162!> \brief ...
1163!> \param r1 ...
1164!> \param i ...
1165!> \param j ...
1166!> \param d ...
1167!> \param rho ...
1168!> \return ...
1169! **************************************************************************************************
1170 FUNCTION dist_second_deriv(r1, i, j, d, rho) RESULT(deriv)
1171 REAL(kind=dp), DIMENSION(3) :: r1
1172 INTEGER :: i, j
1173 REAL(kind=dp) :: d, rho, deriv
1174
1175 deriv = 0.45_dp*rho*(r1(i)*r1(j))/d**2
1176 END FUNCTION
1177
1178! **************************************************************************************************
1179!> \brief ...
1180!> \param r_ij ...
1181!> \param d_ij ...
1182!> \param rho_ij ...
1183!> \param idir ...
1184!> \param jdir ...
1185!> \param iat_der ...
1186!> \param jat_der ...
1187!> \param natom ...
1188!> \return ...
1189! **************************************************************************************************
1190 FUNCTION angle_second_deriv(r_ij, d_ij, rho_ij, idir, jdir, iat_der, jat_der, natom) RESULT(deriv)
1191 REAL(kind=dp), DIMENSION(:, :, :) :: r_ij
1192 REAL(kind=dp), DIMENSION(:, :) :: d_ij, rho_ij
1193 INTEGER :: idir, jdir, iat_der, jat_der, natom
1194 REAL(kind=dp) :: deriv
1195
1196 INTEGER :: i, iat, idr, j, jat, jdr
1197 REAL(kind=dp) :: d12, d23, d31, d_mat(3, 2), denom1, &
1198 denom2, denom3, ka1, ka2, ka3, rho12, &
1199 rho23, rho31, rsst1, rsst2, rsst3
1200 REAL(kind=dp), DIMENSION(3) :: r12, r23, r31
1201
1202 deriv = 0._dp
1203 IF (iat_der == jat_der) THEN
1204 DO i = 1, natom - 1
1205 IF (rho_ij(iat_der, i) .LT. 0.00001) cycle
1206 DO j = i + 1, natom
1207 IF (rho_ij(iat_der, j) .LT. 0.00001) cycle
1208 IF (i == iat_der .OR. j == iat_der) cycle
1209 IF (iat_der .LT. i .OR. iat_der .GT. j) THEN
1210 r12 = r_ij(iat_der, i, :); r23 = r_ij(i, j, :); r31 = r_ij(j, iat_der, :)
1211 d12 = d_ij(iat_der, i); d23 = d_ij(i, j); d31 = d_ij(j, iat_der)
1212 rho12 = rho_ij(iat_der, i); rho23 = rho_ij(i, j); rho31 = rho_ij(j, iat_der)
1213 ELSE
1214 r12 = r_ij(iat_der, j, :); r23 = r_ij(j, i, :); r31 = r_ij(i, iat_der, :)
1215 d12 = d_ij(iat_der, j); d23 = d_ij(j, i); d31 = d_ij(i, iat_der)
1216 rho12 = rho_ij(iat_der, j); rho23 = rho_ij(j, i); rho31 = rho_ij(i, iat_der)
1217 END IF
1218 ka1 = 0.15_dp*rho12*rho23; ka2 = 0.15_dp*rho23*rho31; ka3 = 0.15_dp*rho31*rho12
1219 rsst1 = dot_product(r12, r23); rsst2 = dot_product(r23, r31); rsst3 = dot_product(r31, r12)
1220 denom1 = 1.0_dp - rsst1**2/(d12**2*d23**2); denom2 = 1.0_dp - rsst2**2/(d23**2*d31**2)
1221 denom3 = 1.0_dp - rsst3**2/(d31**2*d12**2)
1222 denom1 = sign(1.0_dp, denom1)*max(abs(denom1), 0.01_dp)
1223 denom2 = sign(1.0_dp, denom2)*max(abs(denom2), 0.01_dp)
1224 denom3 = sign(1.0_dp, denom3)*max(abs(denom3), 0.01_dp)
1225 d_mat(1, 1) = r23(idir)/(d12*d23) - rsst1*r12(idir)/(d12**3*d23)
1226 d_mat(1, 2) = r23(jdir)/(d12*d23) - rsst1*r12(jdir)/(d12**3*d23)
1227 d_mat(2, 1) = -r23(idir)/(d23*d31) + rsst2*r31(idir)/(d23*d31**3)
1228 d_mat(2, 2) = -r23(jdir)/(d23*d31) + rsst2*r31(jdir)/(d23*d31**3)
1229 d_mat(3, 1) = (r31(idir) - r12(idir))/(d31*d12) + rsst3*r31(idir)/(d31**3*d12) - &
1230 rsst3*r12(idir)/(d31*d12**3)
1231 d_mat(3, 2) = (r31(jdir) - r12(jdir))/(d31*d12) + rsst3*r31(jdir)/(d31**3*d12) - &
1232 rsst3*r12(jdir)/(d31*d12**3)
1233 IF (abs(denom1) .LE. 0.011_dp) d_mat(1, 1) = 0.0_dp
1234 IF (abs(denom2) .LE. 0.011_dp) d_mat(2, 1) = 0.0_dp
1235 IF (abs(denom3) .LE. 0.011_dp) d_mat(3, 1) = 0.0_dp
1236 deriv = deriv + ka1*d_mat(1, 1)*d_mat(1, 2)/denom1 + &
1237 ka2*d_mat(2, 1)*d_mat(2, 2)/denom2 + &
1238 ka3*d_mat(3, 1)*d_mat(3, 2)/denom3
1239
1240 END DO
1241 END DO
1242 ELSE
1243 DO i = 1, natom
1244 IF (i == iat_der .OR. i == jat_der) cycle
1245 IF (jat_der .LT. iat_der) THEN
1246 iat = jat_der; jat = iat_der; idr = jdir; jdr = idir
1247 ELSE
1248 iat = iat_der; jat = jat_der; idr = idir; jdr = jdir
1249 END IF
1250 IF (jat .LT. i .OR. iat .GT. i) THEN
1251 r12 = r_ij(iat, jat, :); r23 = r_ij(jat, i, :); r31 = r_ij(i, iat, :)
1252 d12 = d_ij(iat, jat); d23 = d_ij(jat, i); d31 = d_ij(i, iat)
1253 rho12 = rho_ij(iat, jat); rho23 = rho_ij(jat, i); rho31 = rho_ij(i, iat)
1254 ELSE
1255 r12 = r_ij(iat, i, :); r23 = r_ij(i, jat, :); r31 = r_ij(jat, iat, :)
1256 d12 = d_ij(iat, i); d23 = d_ij(i, jat); d31 = d_ij(jat, iat)
1257 rho12 = rho_ij(iat, i); rho23 = rho_ij(i, jat); rho31 = rho_ij(jat, iat)
1258 END IF
1259 ka1 = 0.15_dp*rho12*rho23; ka2 = 0.15_dp*rho23*rho31; ka3 = 0.15_dp*rho31*rho12
1260 rsst1 = dot_product(r12, r23); rsst2 = dot_product(r23, r31); rsst3 = dot_product(r31, r12)
1261 denom1 = 1.0_dp - rsst1**2/(d12**2*d23**2); denom2 = 1.0_dp - rsst2**2/(d23**2*d31**2)
1262 denom3 = 1.0_dp - rsst3**2/(d31**2*d12**2)
1263 denom1 = sign(1.0_dp, denom1)*max(abs(denom1), 0.01_dp)
1264 denom2 = sign(1.0_dp, denom2)*max(abs(denom2), 0.01_dp)
1265 denom3 = sign(1.0_dp, denom3)*max(abs(denom3), 0.01_dp)
1266 d_mat(1, 1) = r23(idr)/(d12*d23) - rsst1*r12(idr)/(d12**3*d23)
1267 d_mat(2, 1) = -r23(idr)/(d23*d31) + rsst2*r31(idr)/(d23*d31**3)
1268 d_mat(3, 1) = (r31(idr) - r12(idr))/(d31*d12) + rsst3*r31(idr)/(d31**3*d12) - &
1269 rsst3*r12(idr)/(d31*d12**3)
1270 IF (jat .LT. i .OR. iat .GT. i) THEN
1271 d_mat(1, 2) = (r12(jdr) - r23(jdr))/(d12*d23) + rsst1*r12(jdr)/(d12**3*d23) - &
1272 rsst1*r23(jdr)/(d12*d23**3)
1273 d_mat(2, 2) = r31(jdr)/(d23*d31) - rsst2*r23(jdr)/(d23**3*d31)
1274 d_mat(3, 2) = -r31(jdr)/(d31*d12) + rsst3*r12(jdr)/(d31*d12**3)
1275 ELSE
1276 d_mat(1, 2) = -r12(jdr)/(d12*d23) + rsst1*r23(jdr)/(d12*d23**3)
1277 d_mat(2, 2) = (r23(jdr) - r31(jdr))/(d23*d31) + rsst2*r23(jdr)/(d23**3*d31) - &
1278 rsst2*r31(jdr)/(d23*d31**3)
1279 d_mat(3, 2) = r12(jdr)/(d31*d12) - rsst3*r31(jdr)/(d31**3*d12)
1280 END IF
1281 IF (abs(denom1) .LE. 0.011_dp) d_mat(1, 1) = 0.0_dp
1282 IF (abs(denom2) .LE. 0.011_dp) d_mat(2, 1) = 0.0_dp
1283 IF (abs(denom3) .LE. 0.011_dp) d_mat(3, 1) = 0.0_dp
1284
1285 deriv = deriv + ka1*d_mat(1, 1)*d_mat(1, 2)/denom1 + &
1286 ka2*d_mat(2, 1)*d_mat(2, 2)/denom2 + &
1287 ka3*d_mat(3, 1)*d_mat(3, 2)/denom3
1288 END DO
1289 END IF
1290 deriv = 0.25_dp*deriv
1291
1292 END FUNCTION angle_second_deriv
1293
1294END MODULE bfgs_optimizer
subroutine cp_eval_at(gopt_env, x, f, gradient, master, final_evaluation, para_env)
evaluete the potential energy and its gradients using an array with same dimension as the particle_se...
represent a simple array based list of the given type
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Routines for Geometry optimization using BFGS algorithm.
recursive subroutine, public geoopt_bfgs(force_env, gopt_param, globenv, geo_section, gopt_env, x0)
Main driver for BFGS geometry optimizations.
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public fix_atom_control(force_env, w)
allows for fix atom constraints
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Routines to handle the external control of CP2K.
subroutine, public external_control(should_stop, flag, globenv, target_time, start_time, force_check)
External manipulations during a run : when the <PROJECT_NAME>.EXIT_$runtype command is sent the progr...
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
Definition cp_files.F:308
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Definition cp_files.F:119
basic linear algebra operations for full matrices
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
subroutine, public cp_fm_transpose(matrix, matrixt)
transposes a matrix matrixt = matrix ^ T
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
Definition cp_fm_diag.F:229
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_write_unformatted(fm, unit)
...
subroutine, public cp_fm_read_unformatted(fm, unit)
...
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public cp_p_file
subroutine, public cp_iterate(iteration_info, last, iter_nr, increment, iter_nr_out)
adds one to the actual iteration
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
types that represent a subsys, i.e. a part of the system
subroutine, public cp_subsys_get(subsys, ref_count, atomic_kinds, atomic_kind_set, particles, particle_set, local_particles, molecules, molecule_set, molecule_kinds, molecule_kind_set, local_molecules, para_env, colvar_p, shell_particles, core_particles, gci, multipoles, natom, nparticle, ncore, nshell, nkind, atprop, virial, results, cell)
returns information about various attributes of the given subsys
Interface for the force calculations.
recursive subroutine, public force_env_get(force_env, in_use, fist_env, qs_env, meta_env, fp_env, subsys, para_env, potential_energy, additional_potential, kinetic_energy, harmonic_shell, kinetic_shell, cell, sub_force_env, qmmm_env, qmmmx_env, eip_env, pwdft_env, globenv, input, force_env_section, method_name_id, root_section, mixed_env, nnp_env, embed_env, ipi_env)
returns various attributes about the force environment
Define type storing the global information of a run. Keep the amount of stored data small....
contains a functional that calculates the energy and its derivatives for the geometry optimizer
subroutine, public print_geo_opt_header(gopt_env, output_unit, label)
...
subroutine, public gopt_f_io_init(gopt_env, output_unit, opt_energy, wildcard, its, used_time)
Handles the Output during an optimization run.
recursive subroutine, public gopt_f_io_finalize(gopt_env, force_env, x0, conv, its, root_section, para_env, master, output_unit)
Handles the Output at the end of an optimization run.
subroutine, public gopt_f_io(gopt_env, force_env, root_section, its, opt_energy, output_unit, eold, emin, wildcard, gopt_param, ndf, dx, xi, conv, pred, rat, step, rad, used_time)
Handles the Output during an optimization run.
subroutine, public print_geo_opt_nc(gopt_env, output_unit)
...
subroutine, public gopt_f_ii(its, output_unit)
Prints iteration step of the optimization procedure on screen.
contains a functional that calculates the energy and its derivatives for the geometry optimizer
contains typo and related routines to handle parameters controlling the GEO_OPT module
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public default_cell_method_id
integer, parameter, public default_ts_method_id
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_path_length
Definition kinds.F:58
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:130
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:147
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
represent a simple array based list of the given type
Space Group Symmetry Type Module (version 1.0, Ferbruary 12, 2021)
Space Group Symmetry Module (version 1.0, January 16, 2020)
subroutine, public print_spgr(spgr)
routine prints Space Group Information.
subroutine, public spgr_apply_rotations_coord(spgr, coord)
routine applies the rotation matrices to the coordinates.
subroutine, public identify_space_group(subsys, geo_section, gopt_env, iunit)
routine indentifies the space group and finds rotation matrices.
subroutine, public spgr_apply_rotations_force(spgr, force)
routine applies the rotation matrices to the forces.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
represents a system: atoms, molecules, their pos,vel,...
wrapper to abstract the force evaluation of the various methods
contains the initially parsed file and the initial parallel environment
calculates the potential energy of a system, and its derivatives
stores all the informations relevant to an mpi environment