34#include "./base/base_uses.f90"
42 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qs_ot_minimizer'
63 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
66 CHARACTER(len=*),
PARAMETER :: routinen =
'ot_mini'
68 INTEGER :: handle, ispin, nspin
69 LOGICAL :: do_ener, do_ks
72 CALL timeset(routinen, handle)
74 nspin =
SIZE(qs_ot_env)
76 do_ks = qs_ot_env(1)%settings%ks
77 do_ener = qs_ot_env(1)%settings%do_ener
79 qs_ot_env(1)%OT_METHOD_FULL =
""
82 IF (.NOT. qs_ot_env(1)%energy_only)
THEN
83 qs_ot_env(1)%gradient = 0.0_dp
86 SELECT CASE (qs_ot_env(1)%settings%ot_algorithm)
89 qs_ot_env(ispin)%matrix_sx, &
90 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin))
93 qs_ot_env(ispin)%matrix_x, qs_ot_env(ispin)%matrix_sx, &
94 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin))
96 cpabort(
"ALGORITHM NYI")
100 IF (qs_ot_env(1)%use_dx)
THEN
102 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_dx, tmp)
103 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient + tmp
104 IF (qs_ot_env(1)%settings%do_rotation)
THEN
105 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_dx, tmp)
106 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient + 0.5_dp*tmp
110 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_dx)
111 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient + tmp
115 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx, tmp)
116 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient - tmp
117 IF (qs_ot_env(1)%settings%do_rotation)
THEN
118 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx, tmp)
119 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient - 0.5_dp*tmp
123 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_gx)
124 qs_ot_env(1)%gradient = qs_ot_env(1)%gradient - tmp
130 SELECT CASE (qs_ot_env(1)%settings%OT_METHOD)
132 IF (current_point_is_fine(qs_ot_env))
THEN
133 qs_ot_env(1)%OT_METHOD_FULL =
"OT CG"
134 CALL ot_new_cg_direction(qs_ot_env)
135 qs_ot_env(1)%line_search_count = 0
137 qs_ot_env(1)%OT_METHOD_FULL =
"OT LS"
139 CALL do_line_search(qs_ot_env)
141 IF (current_point_is_fine(qs_ot_env))
THEN
142 qs_ot_env(1)%OT_METHOD_FULL =
"OT SD"
143 CALL ot_new_sd_direction(qs_ot_env)
144 qs_ot_env(1)%line_search_count = 0
146 qs_ot_env(1)%OT_METHOD_FULL =
"OT LS"
148 CALL do_line_search(qs_ot_env)
150 qs_ot_env(1)%OT_METHOD_FULL =
"OT DIIS"
151 CALL ot_diis_step(qs_ot_env)
153 qs_ot_env(1)%OT_METHOD_FULL =
"OT BROY"
154 CALL ot_broyden_step(qs_ot_env)
156 cpabort(
"OT_METHOD NYI")
159 CALL timestop(handle)
172 FUNCTION current_point_is_fine(qs_ot_env)
RESULT(res)
173 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
179 IF (.NOT. qs_ot_env(1)%energy_only)
THEN
182 IF (qs_ot_env(1)%line_search_count .EQ. 0)
THEN
187 IF (qs_ot_env(1)%line_search_might_be_done)
THEN
195 END FUNCTION current_point_is_fine
204 SUBROUTINE do_line_search(qs_ot_env)
205 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
207 SELECT CASE (qs_ot_env(1)%settings%line_search_method)
209 CALL do_line_search_gold(qs_ot_env)
211 CALL do_line_search_3pnt(qs_ot_env)
213 CALL do_line_search_2pnt(qs_ot_env)
215 CALL do_line_search_none(qs_ot_env)
219 END SUBROUTINE do_line_search
228 SUBROUTINE take_step(ds, qs_ot_env)
230 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
232 CHARACTER(len=*),
PARAMETER :: routinen =
'take_step'
234 INTEGER :: handle, ispin, nspin
235 LOGICAL :: do_ener, do_ks
237 CALL timeset(routinen, handle)
239 nspin =
SIZE(qs_ot_env)
241 do_ks = qs_ot_env(1)%settings%ks
242 do_ener = qs_ot_env(1)%settings%do_ener
246 IF (qs_ot_env(1)%use_dx)
THEN
249 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, qs_ot_env(ispin)%matrix_dx, &
250 alpha_scalar=1.0_dp, beta_scalar=ds)
251 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
252 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, qs_ot_env(ispin)%rot_mat_dx, &
253 alpha_scalar=1.0_dp, beta_scalar=ds)
259 qs_ot_env(ispin)%ener_x = qs_ot_env(ispin)%ener_x + ds*qs_ot_env(ispin)%ener_dx
265 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, qs_ot_env(ispin)%matrix_gx, &
266 alpha_scalar=1.0_dp, beta_scalar=-ds)
267 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
268 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, qs_ot_env(ispin)%rot_mat_gx, &
269 alpha_scalar=1.0_dp, beta_scalar=-ds)
275 qs_ot_env(ispin)%ener_x = qs_ot_env(ispin)%ener_x - ds*qs_ot_env(ispin)%ener_gx
279 CALL timestop(handle)
280 END SUBROUTINE take_step
287 SUBROUTINE do_line_search_gold(qs_ot_env)
289 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
291 CHARACTER(len=*),
PARAMETER :: routinen =
'do_line_search_gold'
292 REAL(kind=
dp),
PARAMETER :: gold_sec = 0.3819_dp
294 INTEGER :: count, handle
297 CALL timeset(routinen, handle)
299 qs_ot_env(1)%line_search_count = qs_ot_env(1)%line_search_count + 1
300 count = qs_ot_env(1)%line_search_count
301 qs_ot_env(1)%line_search_might_be_done = .false.
302 qs_ot_env(1)%energy_only = .true.
304 IF (count + 1 .GT.
SIZE(qs_ot_env(1)%OT_pos))
THEN
307 cpabort(
"MAX ITER EXCEEDED : FATAL")
310 IF (qs_ot_env(1)%line_search_count .EQ. 1)
THEN
311 qs_ot_env(1)%line_search_left = 1
312 qs_ot_env(1)%line_search_right = 0
313 qs_ot_env(1)%line_search_mid = 1
314 qs_ot_env(1)%ot_pos(1) = 0.0_dp
315 qs_ot_env(1)%ot_energy(1) = qs_ot_env(1)%etotal
316 qs_ot_env(1)%ot_pos(2) = qs_ot_env(1)%ds_min/gold_sec
318 qs_ot_env(1)%ot_energy(count) = qs_ot_env(1)%etotal
322 IF (qs_ot_env(1)%line_search_right .EQ. 0)
THEN
323 IF (qs_ot_env(1)%ot_energy(count - 1) .LT. qs_ot_env(1)%ot_energy(count))
THEN
324 qs_ot_env(1)%line_search_right = count
325 qs_ot_env(1)%ot_pos(count + 1) = qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid) + &
326 (qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_right) - &
327 qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid))*gold_sec
329 qs_ot_env(1)%line_search_left = qs_ot_env(1)%line_search_mid
330 qs_ot_env(1)%line_search_mid = count
331 qs_ot_env(1)%ot_pos(count + 1) = qs_ot_env(1)%ot_pos(count)/gold_sec
335 IF (qs_ot_env(1)%ot_pos(count) .LT. qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid))
THEN
336 IF (qs_ot_env(1)%ot_energy(count) .LT. qs_ot_env(1)%ot_energy(qs_ot_env(1)%line_search_mid))
THEN
337 qs_ot_env(1)%line_search_right = qs_ot_env(1)%line_search_mid
338 qs_ot_env(1)%line_search_mid = count
340 qs_ot_env(1)%line_search_left = count
343 IF (qs_ot_env(1)%ot_energy(count) .LT. qs_ot_env(1)%ot_energy(qs_ot_env(1)%line_search_mid))
THEN
344 qs_ot_env(1)%line_search_left = qs_ot_env(1)%line_search_mid
345 qs_ot_env(1)%line_search_mid = count
347 qs_ot_env(1)%line_search_right = count
351 IF ((qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_right) &
352 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid)) .GT. &
353 (qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid) &
354 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_left)))
THEN
355 qs_ot_env(1)%ot_pos(count + 1) = &
356 qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid) + &
357 gold_sec*(qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_right) &
358 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid))
360 qs_ot_env(1)%ot_pos(count + 1) = &
361 qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_left) + &
362 gold_sec*(qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid) &
363 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_left))
366 IF (((qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_right) &
367 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid)) .LT. &
368 qs_ot_env(1)%ds_min*qs_ot_env(1)%settings%gold_target) .AND. &
369 ((qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_mid) &
370 - qs_ot_env(1)%ot_pos(qs_ot_env(1)%line_search_left)) .LT. &
371 qs_ot_env(1)%ds_min*qs_ot_env(1)%settings%gold_target))
THEN
372 qs_ot_env(1)%energy_only = .false.
373 qs_ot_env(1)%line_search_might_be_done = .true.
377 ds = qs_ot_env(1)%OT_pos(count + 1) - qs_ot_env(1)%OT_pos(count)
378 qs_ot_env(1)%ds_min = qs_ot_env(1)%OT_pos(count + 1)
380 CALL take_step(ds, qs_ot_env)
382 CALL timestop(handle)
384 END SUBROUTINE do_line_search_gold
390 SUBROUTINE do_line_search_3pnt(qs_ot_env)
392 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
394 CHARACTER(len=*),
PARAMETER :: routinen =
'do_line_search_3pnt'
396 INTEGER :: count, handle
397 REAL(kind=
dp) :: denom, ds, fa, fb, fc, nom, pos, val, &
400 CALL timeset(routinen, handle)
402 qs_ot_env(1)%line_search_might_be_done = .false.
403 qs_ot_env(1)%energy_only = .true.
406 qs_ot_env(1)%line_search_count = qs_ot_env(1)%line_search_count + 1
407 count = qs_ot_env(1)%line_search_count
408 qs_ot_env(1)%ot_energy(count) = qs_ot_env(1)%etotal
411 qs_ot_env(1)%ot_pos(count) = 0.0_dp
412 qs_ot_env(1)%ot_pos(count + 1) = qs_ot_env(1)%ds_min*0.8_dp
414 IF (qs_ot_env(1)%OT_energy(count) .GT. qs_ot_env(1)%OT_energy(count - 1))
THEN
415 qs_ot_env(1)%OT_pos(count + 1) = qs_ot_env(1)%ds_min*0.5_dp
417 qs_ot_env(1)%OT_pos(count + 1) = qs_ot_env(1)%ds_min*1.4_dp
420 xa = qs_ot_env(1)%OT_pos(1)
421 xb = qs_ot_env(1)%OT_pos(2)
422 xc = qs_ot_env(1)%OT_pos(3)
423 fa = qs_ot_env(1)%OT_energy(1)
424 fb = qs_ot_env(1)%OT_energy(2)
425 fc = qs_ot_env(1)%OT_energy(3)
426 nom = (xb - xa)**2*(fb - fc) - (xb -
xc)**2*(fb - fa)
427 denom = (xb - xa)*(fb - fc) - (xb -
xc)*(fb - fa)
428 IF (abs(denom) .LE. 1.0e-18_dp*max(abs(fb - fc), abs(fb - fa)))
THEN
431 pos = xb - 0.5_dp*nom/denom
433 val = (pos - xa)*(pos - xb)*fc/((
xc - xa)*(
xc - xb)) + &
434 (pos - xb)*(pos -
xc)*fa/((xa - xb)*(xa -
xc)) + &
435 (pos -
xc)*(pos - xa)*fb/((xb -
xc)*(xb - xa))
436 IF (val .LT. fa .AND. val .LE. fb .AND. val .LE. fc)
THEN
438 qs_ot_env(1)%OT_pos(count + 1) = max(maxval(qs_ot_env(1)%OT_pos(1:3))*0.01_dp, &
439 min(pos, maxval(qs_ot_env(1)%OT_pos(1:3))*4.0_dp))
441 qs_ot_env(1)%OT_pos(count + 1) = maxval(qs_ot_env(1)%OT_pos(1:3))*2.0_dp
443 qs_ot_env(1)%energy_only = .false.
444 qs_ot_env(1)%line_search_might_be_done = .true.
448 ds = qs_ot_env(1)%OT_pos(count + 1) - qs_ot_env(1)%OT_pos(count)
449 qs_ot_env(1)%ds_min = qs_ot_env(1)%OT_pos(count + 1)
451 CALL take_step(ds, qs_ot_env)
453 CALL timestop(handle)
455 END SUBROUTINE do_line_search_3pnt
461 SUBROUTINE do_line_search_2pnt(qs_ot_env)
463 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
465 CHARACTER(len=*),
PARAMETER :: routinen =
'do_line_search_2pnt'
467 INTEGER :: count, handle
468 REAL(kind=
dp) :: a, b, c, ds, pos, val, x0, x1
470 CALL timeset(routinen, handle)
472 qs_ot_env(1)%line_search_might_be_done = .false.
473 qs_ot_env(1)%energy_only = .true.
476 qs_ot_env(1)%line_search_count = qs_ot_env(1)%line_search_count + 1
477 count = qs_ot_env(1)%line_search_count
478 qs_ot_env(1)%ot_energy(count) = qs_ot_env(1)%etotal
481 qs_ot_env(1)%ot_pos(count) = 0.0_dp
482 qs_ot_env(1)%ot_grad(count) = qs_ot_env(1)%gradient
483 qs_ot_env(1)%ot_pos(count + 1) = qs_ot_env(1)%ds_min*1.0_dp
486 c = qs_ot_env(1)%ot_energy(1)
487 b = qs_ot_env(1)%ot_grad(1)
488 x1 = qs_ot_env(1)%ot_pos(2)
489 a = (qs_ot_env(1)%ot_energy(2) - b*x1 - c)/(x1**2)
490 IF (a .LE. 0.0_dp) a = 1.0e-15_dp
492 val = a*pos**2 + b*pos + c
493 qs_ot_env(1)%energy_only = .false.
494 qs_ot_env(1)%line_search_might_be_done = .true.
495 IF (val .LT. qs_ot_env(1)%ot_energy(1) .AND. val .LE. qs_ot_env(1)%ot_energy(2))
THEN
498 qs_ot_env(1)%OT_pos(count + 1) = max(maxval(qs_ot_env(1)%OT_pos(1:2))*0.01_dp, &
499 min(pos, maxval(qs_ot_env(1)%OT_pos(1:2))*4.0_dp))
501 qs_ot_env(1)%OT_pos(count + 1) = maxval(qs_ot_env(1)%OT_pos(1:2))*2.0_dp
506 ds = qs_ot_env(1)%OT_pos(count + 1) - qs_ot_env(1)%OT_pos(count)
507 qs_ot_env(1)%ds_min = qs_ot_env(1)%OT_pos(count + 1)
509 CALL take_step(ds, qs_ot_env)
511 CALL timestop(handle)
513 END SUBROUTINE do_line_search_2pnt
519 SUBROUTINE do_line_search_none(qs_ot_env)
520 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
522 CALL take_step(qs_ot_env(1)%ds_min, qs_ot_env)
524 END SUBROUTINE do_line_search_none
535 SUBROUTINE ot_new_sd_direction(qs_ot_env)
536 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
538 CHARACTER(len=*),
PARAMETER :: routinen =
'ot_new_sd_direction'
540 INTEGER :: handle, ispin, itmp, k, n, nener, nspin
541 LOGICAL :: do_ener, do_ks
545 CALL timeset(routinen, handle)
549 nspin =
SIZE(qs_ot_env)
551 do_ks = qs_ot_env(1)%settings%ks
552 do_ener = qs_ot_env(1)%settings%do_ener
554 IF (
ASSOCIATED(qs_ot_env(1)%preconditioner))
THEN
555 IF (.NOT. qs_ot_env(1)%use_dx) cpabort(
"use dx")
556 qs_ot_env(1)%gnorm = 0.0_dp
560 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_dx)
561 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_dx, tmp)
562 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
564 IF (qs_ot_env(1)%gnorm .LT. 0.0_dp)
THEN
569 CALL dbcsr_scale(qs_ot_env(ispin)%matrix_dx, -1.0_dp)
571 IF (qs_ot_env(1)%settings%do_rotation)
THEN
574 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_dx, qs_ot_env(ispin)%rot_mat_gx)
575 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_dx, tmp)
577 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
580 CALL dbcsr_scale(qs_ot_env(ispin)%rot_mat_dx, -1.0_dp)
586 qs_ot_env(ispin)%ener_dx = qs_ot_env(ispin)%ener_gx
587 tmp = dot_product(qs_ot_env(ispin)%ener_dx, qs_ot_env(ispin)%ener_gx)
588 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
589 qs_ot_env(ispin)%ener_dx = -qs_ot_env(ispin)%ener_dx
593 qs_ot_env(1)%gnorm = 0.0_dp
596 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx, tmp)
597 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
599 IF (qs_ot_env(1)%settings%do_rotation)
THEN
601 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx, tmp)
603 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
609 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_gx)
610 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
621 CALL dbcsr_get_info(qs_ot_env(ispin)%matrix_x, nfullcols_total=itmp)
627 nener = nener +
SIZE(qs_ot_env(ispin)%ener_x)
631 IF (int(n, kind=
int_8)*int(k, kind=
int_8) + nener /= 0)
THEN
632 qs_ot_env(1)%delta = sqrt(abs(qs_ot_env(1)%gnorm)/(int(n, kind=
int_8)*int(k, kind=
int_8) + nener))
633 qs_ot_env(1)%gradient = -qs_ot_env(1)%gnorm
635 qs_ot_env(1)%delta = 0.0_dp
636 qs_ot_env(1)%gradient = 0.0_dp
639 CALL timestop(handle)
641 END SUBROUTINE ot_new_sd_direction
652 SUBROUTINE ot_new_cg_direction(qs_ot_env)
653 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
655 CHARACTER(len=*),
PARAMETER :: routinen =
'ot_new_cg_direction'
657 INTEGER :: handle, ispin, itmp, k, n, nener, nspin
658 LOGICAL :: do_ener, do_ks
659 REAL(kind=
dp) :: beta_pr, gnorm_cross, test_down, tmp
662 CALL timeset(routinen, handle)
664 nspin =
SIZE(qs_ot_env)
667 do_ks = qs_ot_env(1)%settings%ks
668 do_ener = qs_ot_env(1)%settings%do_ener
673 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx_old, tmp)
674 gnorm_cross = gnorm_cross + tmp
676 IF (qs_ot_env(1)%settings%do_rotation)
THEN
678 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx_old, tmp)
680 gnorm_cross = gnorm_cross + 0.5_dp*tmp
686 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_gx_old)
687 gnorm_cross = gnorm_cross + tmp
691 IF (
ASSOCIATED(qs_ot_env(1)%preconditioner))
THEN
695 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx_old)
697 qs_ot_env(1)%gnorm = 0.0_dp
700 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx_old, tmp)
701 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
703 IF (qs_ot_env(1)%gnorm .LT. 0.0_dp)
THEN
707 CALL dbcsr_copy(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx_old)
709 IF (qs_ot_env(1)%settings%do_rotation)
THEN
712 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_gx_old, qs_ot_env(ispin)%rot_mat_gx)
713 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx_old, tmp)
715 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
718 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx_old)
724 qs_ot_env(ispin)%ener_gx_old = qs_ot_env(ispin)%ener_gx
725 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_gx_old)
726 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
727 qs_ot_env(ispin)%ener_gx = qs_ot_env(ispin)%ener_gx_old
732 qs_ot_env(1)%gnorm = 0.0_dp
734 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx, tmp)
735 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
736 CALL dbcsr_copy(qs_ot_env(ispin)%matrix_gx_old, qs_ot_env(ispin)%matrix_gx)
738 IF (qs_ot_env(1)%settings%do_rotation)
THEN
740 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx, tmp)
742 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
743 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_gx_old, qs_ot_env(ispin)%rot_mat_gx)
749 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_gx)
750 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
751 qs_ot_env(ispin)%ener_gx_old = qs_ot_env(ispin)%ener_gx
762 CALL dbcsr_get_info(qs_ot_env(ispin)%matrix_x, nfullcols_total=itmp)
768 nener = nener +
SIZE(qs_ot_env(ispin)%ener_x)
772 IF (int(n, kind=
int_8)*int(k, kind=
int_8) + nener /= 0)
THEN
773 qs_ot_env(1)%delta = sqrt(abs(qs_ot_env(1)%gnorm)/(int(n, kind=
int_8)*int(k, kind=
int_8) + nener))
774 beta_pr = (qs_ot_env(1)%gnorm - gnorm_cross)/qs_ot_env(1)%gnorm_old
776 qs_ot_env(1)%delta = 0.0_dp
779 beta_pr = max(beta_pr, 0.0_dp)
784 CALL dbcsr_add(qs_ot_env(ispin)%matrix_dx, qs_ot_env(ispin)%matrix_gx, &
785 alpha_scalar=beta_pr, beta_scalar=-1.0_dp)
786 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_dx, tmp)
787 test_down = test_down + tmp
788 IF (qs_ot_env(1)%settings%do_rotation)
THEN
789 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_dx, qs_ot_env(ispin)%rot_mat_gx, &
790 alpha_scalar=beta_pr, beta_scalar=-1.0_dp)
791 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_dx, tmp)
792 test_down = test_down + 0.5_dp*tmp
798 qs_ot_env(ispin)%ener_dx = beta_pr*qs_ot_env(ispin)%ener_dx - qs_ot_env(ispin)%ener_gx
799 tmp = dot_product(qs_ot_env(ispin)%ener_gx, qs_ot_env(ispin)%ener_dx)
800 test_down = test_down + tmp
804 IF (test_down .GE. 0.0_dp)
THEN
808 CALL dbcsr_add(qs_ot_env(ispin)%matrix_dx, qs_ot_env(ispin)%matrix_gx, &
809 alpha_scalar=beta_pr, beta_scalar=-1.0_dp)
810 IF (qs_ot_env(1)%settings%do_rotation)
THEN
811 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_dx, &
812 qs_ot_env(ispin)%rot_mat_gx, alpha_scalar=beta_pr, beta_scalar=-1.0_dp)
818 qs_ot_env(ispin)%ener_dx = beta_pr*qs_ot_env(ispin)%ener_dx - qs_ot_env(ispin)%ener_gx
823 qs_ot_env(1)%gradient = beta_pr*qs_ot_env(1)%gradient - qs_ot_env(1)%gnorm
824 qs_ot_env(1)%gnorm_old = qs_ot_env(1)%gnorm
826 CALL timestop(handle)
828 END SUBROUTINE ot_new_cg_direction
834 SUBROUTINE ot_diis_step(qs_ot_env)
835 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
837 CHARACTER(len=*),
PARAMETER :: routinen =
'ot_diis_step'
839 INTEGER :: diis_bound, diis_m, handle, i, info, &
840 ispin, itmp, j, k, n, nener, nspin
841 LOGICAL :: do_ener, do_ks, do_ot_sd
842 REAL(kind=
dp) :: overlap, tmp, tr_xnew_gx, tr_xold_gx
845 CALL timeset(routinen, handle)
849 do_ks = qs_ot_env(1)%settings%ks
850 do_ener = qs_ot_env(1)%settings%do_ener
851 nspin =
SIZE(qs_ot_env)
853 diis_m = qs_ot_env(1)%settings%diis_m
855 IF (qs_ot_env(1)%diis_iter .LT. diis_m)
THEN
856 diis_bound = qs_ot_env(1)%diis_iter + 1
861 j = mod(qs_ot_env(1)%diis_iter, diis_m) + 1
867 CALL dbcsr_copy(qs_ot_env(ispin)%matrix_h_x(j)%matrix, qs_ot_env(ispin)%matrix_x)
868 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
869 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_h_x(j)%matrix, qs_ot_env(ispin)%rot_mat_x)
875 qs_ot_env(ispin)%ener_h_x(j, :) = qs_ot_env(ispin)%ener_x(:)
878 IF (
ASSOCIATED(qs_ot_env(1)%preconditioner))
THEN
879 qs_ot_env(1)%gnorm = 0.0_dp
883 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_h_e(j)%matrix)
884 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
886 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
888 IF (qs_ot_env(1)%gnorm .LT. 0.0_dp)
THEN
892 CALL dbcsr_scale(qs_ot_env(ispin)%matrix_h_e(j)%matrix, -qs_ot_env(1)%ds_min)
894 IF (qs_ot_env(1)%settings%do_rotation)
THEN
896 CALL dbcsr_copy(qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, qs_ot_env(ispin)%rot_mat_gx)
897 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, &
899 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
902 CALL dbcsr_scale(qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, -qs_ot_env(1)%ds_min)
908 qs_ot_env(ispin)%ener_h_e(j, :) = qs_ot_env(ispin)%ener_gx(:)
909 tmp = dot_product(qs_ot_env(ispin)%ener_h_e(j, :), qs_ot_env(ispin)%ener_gx(:))
910 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
911 qs_ot_env(ispin)%ener_h_e(j, :) = -qs_ot_env(1)%ds_min*qs_ot_env(ispin)%ener_h_e(j, :)
915 qs_ot_env(1)%gnorm = 0.0_dp
918 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx, tmp)
919 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
920 CALL dbcsr_add(qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
921 qs_ot_env(ispin)%matrix_gx, alpha_scalar=0.0_dp, beta_scalar=-qs_ot_env(1)%ds_min)
923 IF (qs_ot_env(1)%settings%do_rotation)
THEN
925 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, qs_ot_env(ispin)%rot_mat_gx, tmp)
926 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + 0.5_dp*tmp
927 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, &
928 qs_ot_env(ispin)%rot_mat_gx, alpha_scalar=0.0_dp, beta_scalar=-qs_ot_env(1)%ds_min)
934 tmp = dot_product(qs_ot_env(ispin)%ener_gx(:), qs_ot_env(ispin)%ener_gx(:))
935 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
936 qs_ot_env(ispin)%ener_h_e(j, :) = -qs_ot_env(1)%ds_min*qs_ot_env(ispin)%ener_gx(:)
946 CALL dbcsr_get_info(qs_ot_env(ispin)%matrix_x, nfullcols_total=itmp)
952 nener = nener +
SIZE(qs_ot_env(ispin)%ener_x)
956 IF (int(n, kind=
int_8)*int(k, kind=
int_8) + nener /= 0)
THEN
957 qs_ot_env(1)%delta = sqrt(abs(qs_ot_env(1)%gnorm)/(int(n, kind=
int_8)*int(k, kind=
int_8) + nener))
958 qs_ot_env(1)%gradient = -qs_ot_env(1)%gnorm
960 qs_ot_env(1)%delta = 0.0_dp
961 qs_ot_env(1)%gradient = 0.0_dp
971 qs_ot_env(1)%ls_diis(i, j) = 0.0_dp
974 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
975 qs_ot_env(ispin)%matrix_h_e(i)%matrix, &
977 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) + tmp
978 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
979 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, &
980 qs_ot_env(ispin)%rot_mat_h_e(i)%matrix, &
982 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) + 0.5_dp*tmp
988 tmp = dot_product(qs_ot_env(ispin)%ener_h_e(j, :), qs_ot_env(ispin)%ener_h_e(i, :))
989 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) + tmp
993 qs_ot_env(1)%ls_diis(i, j) = 0.0_dp
996 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, &
997 qs_ot_env(ispin)%matrix_h_e(i)%matrix, &
999 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) - qs_ot_env(1)%ds_min*tmp
1000 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
1001 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_gx, &
1002 qs_ot_env(ispin)%rot_mat_h_e(i)%matrix, &
1004 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) - qs_ot_env(1)%ds_min*0.5_dp*tmp
1010 tmp = dot_product(qs_ot_env(ispin)%ener_gx(:), qs_ot_env(ispin)%ener_h_e(i, :))
1011 qs_ot_env(1)%ls_diis(i, j) = qs_ot_env(1)%ls_diis(i, j) - qs_ot_env(1)%ds_min*tmp
1015 qs_ot_env(1)%ls_diis(j, i) = qs_ot_env(1)%ls_diis(i, j)
1016 qs_ot_env(1)%ls_diis(i, diis_bound + 1) = 1.0_dp
1017 qs_ot_env(1)%ls_diis(diis_bound + 1, i) = 1.0_dp
1018 qs_ot_env(1)%c_diis(i) = 0.0_dp
1020 qs_ot_env(1)%ls_diis(diis_bound + 1, diis_bound + 1) = 0.0_dp
1021 qs_ot_env(1)%c_diis(diis_bound + 1) = 1.0_dp
1023 qs_ot_env(1)%lss_diis = qs_ot_env(1)%ls_diis
1025 CALL dgesv(diis_bound + 1, 1, qs_ot_env(1)%lss_diis, diis_m + 1, qs_ot_env(1)%ipivot, &
1026 qs_ot_env(1)%c_diis, diis_m + 1, info)
1028 IF (info .NE. 0)
THEN
1036 CALL dbcsr_set(qs_ot_env(ispin)%matrix_x, 0.0_dp)
1037 DO i = 1, diis_bound
1038 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1039 qs_ot_env(ispin)%matrix_h_e(i)%matrix, &
1040 alpha_scalar=1.0_dp, beta_scalar=qs_ot_env(1)%c_diis(i))
1042 DO i = 1, diis_bound
1043 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1044 qs_ot_env(ispin)%matrix_h_x(i)%matrix, &
1045 alpha_scalar=1.0_dp, beta_scalar=qs_ot_env(1)%c_diis(i))
1047 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
1048 CALL dbcsr_set(qs_ot_env(ispin)%rot_mat_x, 0.0_dp)
1049 DO i = 1, diis_bound
1050 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, &
1051 qs_ot_env(ispin)%rot_mat_h_e(i)%matrix, &
1052 alpha_scalar=1.0_dp, beta_scalar=qs_ot_env(1)%c_diis(i))
1054 DO i = 1, diis_bound
1055 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, &
1056 qs_ot_env(ispin)%rot_mat_h_x(i)%matrix, &
1057 alpha_scalar=1.0_dp, beta_scalar=qs_ot_env(1)%c_diis(i))
1064 qs_ot_env(ispin)%ener_x(:) = 0.0_dp
1065 DO i = 1, diis_bound
1066 qs_ot_env(ispin)%ener_x(:) = qs_ot_env(ispin)%ener_x(:) &
1067 + qs_ot_env(1)%c_diis(i)*qs_ot_env(ispin)%ener_h_e(i, :)
1069 DO i = 1, diis_bound
1070 qs_ot_env(ispin)%ener_x(:) = qs_ot_env(ispin)%ener_x(:) &
1071 + qs_ot_env(1)%c_diis(i)*qs_ot_env(ispin)%ener_h_x(i, :)
1075 qs_ot_env(1)%diis_iter = qs_ot_env(1)%diis_iter + 1
1076 IF (qs_ot_env(1)%settings%safer_diis)
THEN
1084 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_h_x(j)%matrix, &
1085 qs_ot_env(ispin)%matrix_gx, tmp)
1086 tr_xold_gx = tr_xold_gx + tmp
1087 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_x, &
1088 qs_ot_env(ispin)%matrix_gx, tmp)
1089 tr_xnew_gx = tr_xnew_gx + tmp
1090 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
1091 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_h_x(j)%matrix, &
1092 qs_ot_env(ispin)%rot_mat_gx, tmp)
1093 tr_xold_gx = tr_xold_gx + 0.5_dp*tmp
1094 CALL dbcsr_dot(qs_ot_env(ispin)%rot_mat_x, &
1095 qs_ot_env(ispin)%rot_mat_gx, tmp)
1096 tr_xnew_gx = tr_xnew_gx + 0.5_dp*tmp
1102 tmp = dot_product(qs_ot_env(ispin)%ener_h_x(j, :), qs_ot_env(ispin)%ener_gx(:))
1103 tr_xold_gx = tr_xold_gx + tmp
1104 tmp = dot_product(qs_ot_env(ispin)%ener_x(:), qs_ot_env(ispin)%ener_gx(:))
1105 tr_xnew_gx = tr_xnew_gx + tmp
1108 overlap = (tr_xnew_gx - tr_xold_gx)
1110 IF (overlap .GT. 0.0_dp)
THEN
1117 qs_ot_env(1)%OT_METHOD_FULL =
"OT SD"
1120 CALL dbcsr_set(qs_ot_env(ispin)%matrix_x, 0.0_dp)
1121 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1122 qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
1124 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1125 qs_ot_env(ispin)%matrix_h_x(j)%matrix, &
1127 IF (qs_ot_env(ispin)%settings%do_rotation)
THEN
1128 CALL dbcsr_set(qs_ot_env(ispin)%rot_mat_x, 0.0_dp)
1129 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, &
1130 qs_ot_env(ispin)%rot_mat_h_e(j)%matrix, &
1132 CALL dbcsr_add(qs_ot_env(ispin)%rot_mat_x, &
1133 qs_ot_env(ispin)%rot_mat_h_x(j)%matrix, &
1140 qs_ot_env(ispin)%ener_x(:) = 0._dp
1141 qs_ot_env(ispin)%ener_x(:) = qs_ot_env(ispin)%ener_x(:) + qs_ot_env(ispin)%ener_h_e(j, :)
1142 qs_ot_env(ispin)%ener_x(:) = qs_ot_env(ispin)%ener_x(:) + qs_ot_env(ispin)%ener_h_x(j, :)
1147 CALL timestop(handle)
1149 END SUBROUTINE ot_diis_step
1156 SUBROUTINE ot_broyden_step(qs_ot_env)
1157 TYPE(
qs_ot_type),
DIMENSION(:),
POINTER :: qs_ot_env
1159 INTEGER :: diis_bound, diis_m, i, ispin, itmp, j, &
1161 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: circ_index
1162 LOGICAL :: adaptive_sigma, do_ener, do_ks, &
1163 enable_flip, forget_history
1164 REAL(kind=
dp) :: beta, eta,
gamma, omega, sigma, &
1165 sigma_dec, sigma_min, tmp, tmp2
1166 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: f, x
1167 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: g, s
1170 eta = qs_ot_env(1)%settings%broyden_eta
1171 omega = qs_ot_env(1)%settings%broyden_omega
1172 sigma_dec = qs_ot_env(1)%settings%broyden_sigma_decrease
1173 sigma_min = qs_ot_env(1)%settings%broyden_sigma_min
1174 forget_history = qs_ot_env(1)%settings%broyden_forget_history
1175 adaptive_sigma = qs_ot_env(1)%settings%broyden_adaptive_sigma
1176 enable_flip = qs_ot_env(1)%settings%broyden_enable_flip
1178 do_ks = qs_ot_env(1)%settings%ks
1179 do_ener = qs_ot_env(1)%settings%do_ener
1181 beta = qs_ot_env(1)%settings%broyden_beta
1182 gamma = qs_ot_env(1)%settings%broyden_gamma
1183 IF (adaptive_sigma)
THEN
1184 IF (qs_ot_env(1)%broyden_adaptive_sigma .LT. 0.0_dp)
THEN
1185 sigma = qs_ot_env(1)%settings%broyden_sigma
1187 sigma = qs_ot_env(1)%broyden_adaptive_sigma
1190 sigma = qs_ot_env(1)%settings%broyden_sigma
1194 IF (do_ener .OR. .NOT. do_ks .OR. qs_ot_env(1)%settings%do_rotation)
THEN
1195 cpabort(
"Not yet implemented")
1198 nspin =
SIZE(qs_ot_env)
1200 diis_m = qs_ot_env(1)%settings%diis_m
1202 IF (qs_ot_env(1)%diis_iter .LT. diis_m)
THEN
1203 diis_bound = qs_ot_env(1)%diis_iter + 1
1209 k = 2*diis_bound + 1
1214 ALLOCATE (circ_index(diis_bound))
1221 j = mod(qs_ot_env(1)%diis_iter, diis_m) + 1
1224 CALL dbcsr_copy(qs_ot_env(ispin)%matrix_h_x(j)%matrix, qs_ot_env(ispin)%matrix_x)
1227 IF (
ASSOCIATED(qs_ot_env(1)%preconditioner))
THEN
1228 qs_ot_env(1)%gnorm = 0.0_dp
1231 qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_h_e(j)%matrix)
1232 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
1234 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
1236 IF (qs_ot_env(1)%gnorm .LT. 0.0_dp)
THEN
1240 CALL dbcsr_scale(qs_ot_env(ispin)%matrix_h_e(j)%matrix, -1.0_dp)
1243 qs_ot_env(1)%gnorm = 0.0_dp
1245 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, qs_ot_env(ispin)%matrix_gx, tmp)
1246 qs_ot_env(1)%gnorm = qs_ot_env(1)%gnorm + tmp
1247 CALL dbcsr_add(qs_ot_env(ispin)%matrix_h_e(j)%matrix, &
1248 qs_ot_env(ispin)%matrix_gx, alpha_scalar=0.0_dp, beta_scalar=-1.0_dp)
1256 CALL dbcsr_get_info(qs_ot_env(ispin)%matrix_x, nfullcols_total=itmp)
1261 IF (int(n, kind=
int_8)*int(k, kind=
int_8) /= 0)
THEN
1262 qs_ot_env(1)%delta = sqrt(abs(qs_ot_env(1)%gnorm)/(int(n, kind=
int_8)*int(k, kind=
int_8)))
1263 qs_ot_env(1)%gradient = -qs_ot_env(1)%gnorm
1265 qs_ot_env(1)%delta = 0.0_dp
1266 qs_ot_env(1)%gradient = 0.0_dp
1269 IF (diis_bound == diis_m)
THEN
1270 DO i = 1, diis_bound
1271 circ_index(i) = mod(j + i - 1, diis_m) + 1
1274 DO i = 1, diis_bound
1282 DO i = 1, diis_bound
1284 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1285 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, tmp)
1286 s(i, i) = s(i, i) + tmp
1288 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1289 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, tmp)
1290 s(i + diis_bound, i + diis_bound) = s(i + diis_bound, i + diis_bound) + tmp
1292 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1293 qs_ot_env(ispin)%matrix_x, tmp)
1294 s(i, 2*diis_bound + 1) = s(i, 2*diis_bound + 1) + tmp
1295 s(i, 2*diis_bound + 1) = s(2*diis_bound + 1, i)
1297 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1298 qs_ot_env(ispin)%matrix_x, tmp)
1299 s(i + diis_bound, 2*diis_bound + 1) = s(i + diis_bound, 2*diis_bound + 1) + tmp
1300 s(i + diis_bound, 2*diis_bound + 1) = s(2*diis_bound + 1, diis_bound + i)
1301 DO k = (i + 1), diis_bound
1303 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1304 qs_ot_env(ispin)%matrix_h_x(circ_index(k))%matrix, &
1306 s(i, k) = s(i, k) + tmp
1309 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1310 qs_ot_env(ispin)%matrix_h_e(circ_index(k))%matrix, &
1312 s(diis_bound + i, diis_bound + k) = s(diis_bound + i, diis_bound + k) + tmp
1313 s(diis_bound + k, diis_bound + i) = s(diis_bound + i, diis_bound + k)
1315 DO k = 1, diis_bound
1317 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1318 qs_ot_env(ispin)%matrix_h_e(circ_index(k))%matrix, tmp)
1319 s(i, k + diis_bound) = s(i, k + diis_bound) + tmp
1320 s(k + diis_bound, i) = s(i, k + diis_bound)
1323 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_x, qs_ot_env(ispin)%matrix_x, tmp)
1324 s(2*diis_bound + 1, 2*diis_bound + 1) = s(2*diis_bound + 1, 2*diis_bound + 1) + tmp
1328 k = 2*diis_bound + 1
1330 s(k, :) = s(k, :)/tmp
1331 s(:, k) = s(:, k)/tmp
1333 IF (diis_bound .GT. 1)
THEN
1340 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1341 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1345 qs_ot_env(ispin)%matrix_h_x(circ_index(i - 1))%matrix, &
1346 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1350 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1351 qs_ot_env(ispin)%matrix_h_e(circ_index(i - 1))%matrix, &
1355 qs_ot_env(ispin)%matrix_h_x(circ_index(i - 1))%matrix, &
1356 qs_ot_env(ispin)%matrix_h_e(circ_index(i - 1))%matrix, &
1360 qs_ot_env(1)%c_broy(i - 1) = tmp2
1363 qs_ot_env(1)%energy_h(j) = qs_ot_env(1)%etotal
1366 i = minloc(qs_ot_env(1)%energy_h(1:diis_bound), dim=1)
1368 sigma = sigma_dec*sigma
1369 qs_ot_env(1)%OT_METHOD_FULL =
"OT BTRK"
1371 CALL dbcsr_set(qs_ot_env(ispin)%matrix_x, 0.0_dp)
1372 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1373 qs_ot_env(ispin)%matrix_h_x(i)%matrix, &
1374 alpha_scalar=1.0_dp, beta_scalar=(1.0_dp -
gamma))
1375 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1376 qs_ot_env(ispin)%matrix_h_x(circ_index(diis_bound))%matrix, &
1377 alpha_scalar=1.0_dp, beta_scalar=
gamma)
1381 DO i = 2, diis_bound
1388 f(diis_bound + i) = 1.0
1389 f(diis_bound + i - 1) = -1.0
1392 IF (enable_flip)
THEN
1393 IF (qs_ot_env(1)%c_broy(i - 1) .GT. 0)
THEN
1400 x(:) = tmp*x - matmul(g, f)
1403 tmp = dot_product(f, matmul(s, f))
1405 f(:) = matmul(s, f)/tmp
1407 g(:, :) = g + spread(x, dim=2, ncopies=
SIZE(f))*spread(f, dim=1, ncopies=
SIZE(x))
1410 f(2*diis_bound) = 1.0_dp
1411 x(:) = -beta*matmul(g, f)
1415 CALL dbcsr_set(qs_ot_env(ispin)%matrix_x, 0.0_dp)
1416 DO i = 1, diis_bound
1417 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1418 qs_ot_env(ispin)%matrix_h_e(circ_index(i))%matrix, &
1419 alpha_scalar=1.0_dp, beta_scalar=-x(i + diis_bound))
1421 DO i = 1, diis_bound
1422 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1423 qs_ot_env(ispin)%matrix_h_x(circ_index(i))%matrix, &
1424 alpha_scalar=1.0_dp, beta_scalar=x(i))
1428 IF (adaptive_sigma)
THEN
1429 tmp = new_sigma(g, s, diis_bound)
1432 sigma = min(omega*sigma, tmp)
1438 CALL dbcsr_dot(qs_ot_env(ispin)%matrix_gx, &
1439 qs_ot_env(ispin)%matrix_x, &
1447 IF (tmp .GE. 0.0_dp)
THEN
1448 qs_ot_env(1)%OT_METHOD_FULL =
"OT TURN"
1449 IF (forget_history)
THEN
1450 qs_ot_env(1)%diis_iter = 0
1452 sigma = sigma*sigma_dec
1453 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1454 qs_ot_env(ispin)%matrix_h_x(circ_index(diis_bound))%matrix, &
1455 alpha_scalar=-1.0_dp, beta_scalar=1.0_dp)
1457 CALL dbcsr_add(qs_ot_env(ispin)%matrix_x, &
1458 qs_ot_env(ispin)%matrix_h_x(circ_index(diis_bound))%matrix, &
1459 alpha_scalar=1.0_dp, beta_scalar=1.0_dp)
1465 DEALLOCATE (s, g, f, x, circ_index)
1468 qs_ot_env(1)%diis_iter = qs_ot_env(1)%diis_iter + 1
1469 qs_ot_env(1)%broyden_adaptive_sigma = max(sigma, sigma_min)
1471 END SUBROUTINE ot_broyden_step
1480 FUNCTION new_sigma(G, S, n)
RESULT(sigma)
1486 REAL(kind=
dp),
DIMENSION(:, :) :: g, s
1488 REAL(kind=
dp) :: sigma
1490 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eigv
1491 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: h
1494 CALL hess_g(g, s, h, n)
1501 sigma = sum(abs(eigv**2))/sum(abs(eigv))
1504 sigma = sum(abs(eigv))/max(1,
SIZE(eigv))
1507 sigma = (maxval(abs(eigv)) + minval(abs(eigv)))*0.5_dp
1510 DEALLOCATE (h, eigv)
1511 END FUNCTION new_sigma
1520 SUBROUTINE hess_g(G, S, H, n)
1527 REAL(kind=
dp),
DIMENSION(:, :) :: g, s, h
1531 REAL(kind=
dp) :: tmp
1532 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: v
1533 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :) :: q
1543 tmp = sqrt(dot_product(q(:, 1), matmul(s, q(:, 1))))
1544 q(:, :) = q(:, :)/tmp
1547 v(:) = matmul(g, q(:, i))
1549 h(j, i) = dot_product(q(:, j), matmul(s, v))
1550 v(:) = v - h(j, i)*q(:, j)
1553 tmp = dot_product(v, matmul(s, v))
1554 IF (tmp .LE. 0.0_dp)
THEN
1560 IF (abs(tmp) .LT. 1e-9_dp)
THEN
1565 q(:, i + 1) = v/h(i + 1, i)
1570 END SUBROUTINE hess_g
subroutine, public dbcsr_scale(matrix, alpha_scalar)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
subroutine, public dbcsr_init_random(matrix, keep_sparsity)
Fills the given matrix with random numbers.
subroutine, public dbcsr_dot(matrix_a, matrix_b, trace)
Computes the dot product of two matrices, also known as the trace of their matrix product.
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
Collection of simple mathematical functions and subroutines.
subroutine, public diamat_all(a, eigval, dac)
Diagonalize the symmetric n by n matrix a using the LAPACK library. Only the upper triangle of matrix...
computes preconditioners, and implements methods to apply them currently used in qs_ot
subroutine, public ot_mini(qs_ot_env, matrix_hc)
...
subroutine, public qs_ot_get_derivative_ref(matrix_hc, matrix_x, matrix_sx, matrix_gx, qs_ot_env)
...
subroutine, public qs_ot_get_derivative(matrix_hc, matrix_x, matrix_sx, matrix_gx, qs_ot_env)
this routines computes dE/dx=dx, with dx ortho to sc0 needs dE/dC=hc,C0,X,SX,p if preconditioned it w...
Exchange and Correlation functional calculations.
type of a logger, at the moment it contains just a print level starting at which level it should be l...