49#include "./base/base_uses.f90"
55 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'dm_ls_scf_methods'
77 CHARACTER(len=*),
PARAMETER :: routinen =
'ls_scf_init_matrix_S'
79 INTEGER :: handle, unit_nr
80 REAL(kind=
dp) :: frob_matrix, frob_matrix_base
84 CALL timeset(routinen, handle)
88 IF (logger%para_env%is_source())
THEN
95 IF (ls_scf_env%has_unit_metric)
THEN
96 CALL dbcsr_set(ls_scf_env%matrix_s, 0.0_dp)
99 CALL matrix_qs_to_ls(ls_scf_env%matrix_s, matrix_s, ls_scf_env%ls_mstruct, covariant=.true.)
102 CALL dbcsr_filter(ls_scf_env%matrix_s, ls_scf_env%eps_filter)
105 IF (ls_scf_env%has_s_preconditioner)
THEN
106 CALL dbcsr_create(ls_scf_env%matrix_bs_sqrt, template=ls_scf_env%matrix_s, &
107 matrix_type=dbcsr_type_no_symmetry)
108 CALL dbcsr_create(ls_scf_env%matrix_bs_sqrt_inv, template=ls_scf_env%matrix_s, &
109 matrix_type=dbcsr_type_no_symmetry)
111 ls_scf_env%s_preconditioner_type, ls_scf_env%ls_mstruct, &
112 ls_scf_env%matrix_bs_sqrt, ls_scf_env%matrix_bs_sqrt_inv, &
113 ls_scf_env%eps_filter, ls_scf_env%s_sqrt_order, &
114 ls_scf_env%eps_lanczos, ls_scf_env%max_iter_lanczos)
118 IF (ls_scf_env%has_s_preconditioner)
THEN
120 ls_scf_env%matrix_bs_sqrt, ls_scf_env%matrix_bs_sqrt_inv)
124 IF (ls_scf_env%use_s_sqrt)
THEN
126 CALL dbcsr_create(ls_scf_env%matrix_s_sqrt, template=ls_scf_env%matrix_s, &
127 matrix_type=dbcsr_type_no_symmetry)
128 CALL dbcsr_create(ls_scf_env%matrix_s_sqrt_inv, template=ls_scf_env%matrix_s, &
129 matrix_type=dbcsr_type_no_symmetry)
131 SELECT CASE (ls_scf_env%s_sqrt_method)
134 ls_scf_env%matrix_s, ls_scf_env%eps_filter, &
135 ls_scf_env%s_sqrt_order, &
136 ls_scf_env%eps_lanczos, ls_scf_env%max_iter_lanczos, &
140 ls_scf_env%matrix_s, ls_scf_env%eps_filter, &
141 ls_scf_env%s_sqrt_order, &
142 ls_scf_env%eps_lanczos, ls_scf_env%max_iter_lanczos, &
145 cpabort(
"Unknown sqrt method.")
148 IF (ls_scf_env%check_s_inv)
THEN
149 CALL dbcsr_create(matrix_tmp1, template=ls_scf_env%matrix_s, &
150 matrix_type=dbcsr_type_no_symmetry)
151 CALL dbcsr_create(matrix_tmp2, template=ls_scf_env%matrix_s, &
152 matrix_type=dbcsr_type_no_symmetry)
154 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, ls_scf_env%matrix_s, &
155 0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
157 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_tmp1, ls_scf_env%matrix_s_sqrt_inv, &
158 0.0_dp, matrix_tmp2, filter_eps=ls_scf_env%eps_filter)
163 IF (unit_nr > 0)
THEN
164 WRITE (unit_nr, *)
"Error for (inv(sqrt(S))*S*inv(sqrt(S))-I)", frob_matrix/frob_matrix_base
173 IF (ls_scf_env%needs_s_inv)
THEN
174 CALL dbcsr_create(ls_scf_env%matrix_s_inv, template=ls_scf_env%matrix_s, &
175 matrix_type=dbcsr_type_no_symmetry)
176 IF (.NOT. ls_scf_env%use_s_sqrt)
THEN
177 CALL invert_hotelling(ls_scf_env%matrix_s_inv, ls_scf_env%matrix_s, ls_scf_env%eps_filter)
179 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, ls_scf_env%matrix_s_sqrt_inv, &
180 0.0_dp, ls_scf_env%matrix_s_inv, filter_eps=ls_scf_env%eps_filter)
182 IF (ls_scf_env%check_s_inv)
THEN
183 CALL dbcsr_create(matrix_tmp1, template=ls_scf_env%matrix_s, &
184 matrix_type=dbcsr_type_no_symmetry)
185 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, ls_scf_env%matrix_s_inv, ls_scf_env%matrix_s, &
186 0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
190 IF (unit_nr > 0)
THEN
191 WRITE (unit_nr, *)
"Error for (inv(S)*S-I)", frob_matrix/frob_matrix_base
197 CALL timestop(handle)
217 matrix_bs_sqrt, matrix_bs_sqrt_inv, threshold, order, &
218 eps_lanczos, max_iter_lanczos)
221 INTEGER :: preconditioner_type
223 TYPE(
dbcsr_type),
INTENT(INOUT) :: matrix_bs_sqrt, matrix_bs_sqrt_inv
224 REAL(kind=
dp) :: threshold
226 REAL(kind=
dp) :: eps_lanczos
227 INTEGER :: max_iter_lanczos
229 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_matrix_preconditioner'
231 INTEGER :: handle, iblock_col, iblock_row
232 LOGICAL :: block_needed
233 REAL(
dp),
DIMENSION(:, :),
POINTER :: block_dp
237 CALL timeset(routinen, handle)
242 SELECT CASE (preconditioner_type)
254 block_needed = .false.
256 IF (iblock_row == iblock_col)
THEN
257 block_needed = .true.
261 IF (ls_mstruct%atom_to_molecule(iblock_row) == ls_mstruct%atom_to_molecule(iblock_col)) block_needed = .true.
266 IF (block_needed)
THEN
267 CALL dbcsr_put_block(matrix=matrix_bs, row=iblock_row, col=iblock_col, block=block_dp)
276 SELECT CASE (preconditioner_type)
284 CALL dbcsr_copy(matrix_bs_sqrt_inv, matrix_bs)
285 CALL dbcsr_set(matrix_bs_sqrt_inv, 0.0_dp)
289 CALL dbcsr_copy(matrix_bs_sqrt_inv, matrix_bs)
295 threshold=min(threshold, 1.0e-10_dp), order=order, &
296 eps_lanczos=eps_lanczos, max_iter_lanczos=max_iter_lanczos, &
302 CALL timestop(handle)
321 CHARACTER(LEN=*) :: direction
322 TYPE(
dbcsr_type),
INTENT(INOUT) :: matrix_bs_sqrt, matrix_bs_sqrt_inv
324 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_matrix_preconditioner'
329 CALL timeset(routinen, handle)
330 CALL dbcsr_create(matrix_tmp, template=matrix, matrix_type=dbcsr_type_no_symmetry)
332 SELECT CASE (direction)
334 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix, matrix_bs_sqrt_inv, &
336 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_bs_sqrt_inv, matrix_tmp, &
341 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_bs_sqrt, matrix_tmp, &
349 CALL timestop(handle)
375 matrix_s, matrix_s_inv, nelectron, threshold, sign_symmetric, submatrix_sign_method, &
379 REAL(kind=
dp),
INTENT(INOUT) :: mu
381 INTEGER :: sign_method, sign_order
382 TYPE(
dbcsr_type),
INTENT(INOUT) :: matrix_ks, matrix_s, matrix_s_inv
383 INTEGER,
INTENT(IN) :: nelectron
384 REAL(kind=
dp),
INTENT(IN) :: threshold
385 LOGICAL,
OPTIONAL :: sign_symmetric
386 INTEGER,
OPTIONAL :: submatrix_sign_method
387 TYPE(
dbcsr_type),
INTENT(IN),
OPTIONAL :: matrix_s_sqrt_inv
389 CHARACTER(LEN=*),
PARAMETER :: routinen =
'density_matrix_sign'
390 REAL(kind=
dp),
PARAMETER :: initial_increment = 0.01_dp
392 INTEGER :: handle, iter, unit_nr, &
393 used_submatrix_sign_method
394 LOGICAL :: do_sign_symmetric, has_mu_high, &
395 has_mu_low, internal_mu_adjust
396 REAL(kind=
dp) :: increment, mu_high, mu_low, trace
399 CALL timeset(routinen, handle)
402 IF (logger%para_env%is_source())
THEN
408 do_sign_symmetric = .false.
409 IF (
PRESENT(sign_symmetric)) do_sign_symmetric = sign_symmetric
412 IF (
PRESENT(submatrix_sign_method)) used_submatrix_sign_method = submatrix_sign_method
418 IF (internal_mu_adjust)
THEN
419 CALL density_matrix_sign_internal_mu(matrix_p, trace, mu, sign_method, &
420 matrix_ks, matrix_s, threshold, &
421 used_submatrix_sign_method, &
422 nelectron, matrix_s_sqrt_inv)
424 increment = initial_increment
427 has_mu_high = .false.
431 IF (has_mu_low .AND. has_mu_high)
THEN
432 mu = (mu_low + mu_high)/2
433 IF (abs(mu_high - mu_low) < threshold)
EXIT
437 matrix_ks, matrix_s, matrix_s_inv, threshold, &
438 do_sign_symmetric, used_submatrix_sign_method, &
440 IF (unit_nr > 0)
WRITE (unit_nr,
'(T2,A,I2,1X,F13.9,1X,F15.9)') &
441 "Density matrix: iter, mu, trace error: ", iter, mu, trace - nelectron
445 IF (abs(trace - nelectron) < 0.5_dp .OR. fixed_mu)
EXIT
447 IF (trace < nelectron)
THEN
451 increment = increment*2
456 increment = increment*2
462 CALL timestop(handle)
485 matrix_s, matrix_s_inv, threshold, sign_symmetric, submatrix_sign_method, &
489 REAL(kind=
dp),
INTENT(OUT) :: trace
490 REAL(kind=
dp),
INTENT(INOUT) :: mu
491 INTEGER :: sign_method, sign_order
492 TYPE(
dbcsr_type),
INTENT(INOUT) :: matrix_ks, matrix_s, matrix_s_inv
493 REAL(kind=
dp),
INTENT(IN) :: threshold
494 LOGICAL :: sign_symmetric
495 INTEGER :: submatrix_sign_method
496 TYPE(
dbcsr_type),
INTENT(IN),
OPTIONAL :: matrix_s_sqrt_inv
498 CHARACTER(LEN=*),
PARAMETER :: routinen =
'density_matrix_sign_fixed_mu'
500 INTEGER :: handle, unit_nr
501 REAL(kind=
dp) :: frob_matrix
503 TYPE(
dbcsr_type) :: matrix_p_ud, matrix_sign, matrix_sinv_ks, matrix_ssqrtinv_ks_ssqrtinv, &
504 matrix_ssqrtinv_ks_ssqrtinv2, matrix_tmp
506 CALL timeset(routinen, handle)
509 IF (logger%para_env%is_source())
THEN
515 CALL dbcsr_create(matrix_sign, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
517 IF (sign_symmetric)
THEN
519 IF (.NOT.
PRESENT(matrix_s_sqrt_inv)) &
520 cpabort(
"Argument matrix_s_sqrt_inv required if sign_symmetric is set")
522 CALL dbcsr_create(matrix_ssqrtinv_ks_ssqrtinv, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
523 CALL dbcsr_create(matrix_ssqrtinv_ks_ssqrtinv2, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
524 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_ks, &
525 0.0_dp, matrix_ssqrtinv_ks_ssqrtinv2, filter_eps=threshold)
526 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_ssqrtinv_ks_ssqrtinv2, matrix_s_sqrt_inv, &
527 0.0_dp, matrix_ssqrtinv_ks_ssqrtinv, filter_eps=threshold)
530 SELECT CASE (sign_method)
534 CALL matrix_sign_proot(matrix_sign, matrix_ssqrtinv_ks_ssqrtinv, threshold, sign_order)
536 CALL matrix_sign_submatrix(matrix_sign, matrix_ssqrtinv_ks_ssqrtinv, threshold, sign_order, submatrix_sign_method)
538 cpabort(
"Unkown sign method.")
545 CALL dbcsr_create(matrix_sinv_ks, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
547 0.0_dp, matrix_sinv_ks, filter_eps=threshold)
551 SELECT CASE (sign_method)
557 CALL matrix_sign_submatrix(matrix_sign, matrix_sinv_ks, threshold, sign_order, submatrix_sign_method)
559 cpabort(
"Unkown sign method.")
565 CALL dbcsr_create(matrix_p_ud, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
575 CALL dbcsr_create(matrix_tmp, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
577 0.0_dp, matrix_tmp, filter_eps=threshold)
578 CALL dbcsr_add(matrix_tmp, matrix_p_ud, 1.0_dp, -1.0_dp)
580 IF (unit_nr > 0 .AND. frob_matrix > 0.001_dp) &
581 WRITE (unit_nr,
'(T2,A,F20.12)')
"Deviation from idempotency: ", frob_matrix
583 IF (sign_symmetric)
THEN
584 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_p_ud, &
585 0.0_dp, matrix_tmp, filter_eps=threshold)
586 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_tmp, matrix_s_sqrt_inv, &
587 0.0_dp, matrix_p, filter_eps=threshold)
591 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_p_ud, matrix_s_inv, &
592 0.0_dp, matrix_p, filter_eps=threshold)
597 CALL timestop(handle)
617 SUBROUTINE density_matrix_sign_internal_mu(matrix_p, trace, mu, sign_method, matrix_ks, &
618 matrix_s, threshold, submatrix_sign_method, &
619 nelectron, matrix_s_sqrt_inv)
622 REAL(kind=
dp),
INTENT(OUT) :: trace
623 REAL(kind=
dp),
INTENT(INOUT) :: mu
624 INTEGER :: sign_method
625 TYPE(
dbcsr_type),
INTENT(INOUT) :: matrix_ks, matrix_s
626 REAL(kind=
dp),
INTENT(IN) :: threshold
627 INTEGER :: submatrix_sign_method
628 INTEGER,
INTENT(IN) :: nelectron
629 TYPE(
dbcsr_type),
INTENT(IN) :: matrix_s_sqrt_inv
631 CHARACTER(LEN=*),
PARAMETER :: routinen =
'density_matrix_sign_internal_mu'
633 INTEGER :: handle, unit_nr
634 REAL(kind=
dp) :: frob_matrix
636 TYPE(
dbcsr_type) :: matrix_p_ud, matrix_sign, &
637 matrix_ssqrtinv_ks_ssqrtinv, &
638 matrix_ssqrtinv_ks_ssqrtinv2, &
641 CALL timeset(routinen, handle)
644 IF (logger%para_env%is_source())
THEN
650 CALL dbcsr_create(matrix_sign, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
652 CALL dbcsr_create(matrix_ssqrtinv_ks_ssqrtinv, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
653 CALL dbcsr_create(matrix_ssqrtinv_ks_ssqrtinv2, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
654 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_ks, &
655 0.0_dp, matrix_ssqrtinv_ks_ssqrtinv2, filter_eps=threshold)
656 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_ssqrtinv_ks_ssqrtinv2, matrix_s_sqrt_inv, &
657 0.0_dp, matrix_ssqrtinv_ks_ssqrtinv, filter_eps=threshold)
660 SELECT CASE (sign_method)
662 SELECT CASE (submatrix_sign_method)
665 submatrix_sign_method)
667 cpabort(
"density_matrix_sign_internal_mu called with invalid submatrix sign method")
670 cpabort(
"density_matrix_sign_internal_mu called with invalid sign method.")
676 CALL dbcsr_create(matrix_p_ud, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
686 CALL dbcsr_create(matrix_tmp, template=matrix_s, matrix_type=dbcsr_type_no_symmetry)
688 0.0_dp, matrix_tmp, filter_eps=threshold)
689 CALL dbcsr_add(matrix_tmp, matrix_p_ud, 1.0_dp, -1.0_dp)
691 IF (unit_nr > 0 .AND. frob_matrix > 0.001_dp) &
692 WRITE (unit_nr,
'(T2,A,F20.12)')
"Deviation from idempotency: ", frob_matrix
694 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_p_ud, &
695 0.0_dp, matrix_tmp, filter_eps=threshold)
696 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_tmp, matrix_s_sqrt_inv, &
697 0.0_dp, matrix_p, filter_eps=threshold)
701 CALL timestop(handle)
703 END SUBROUTINE density_matrix_sign_internal_mu
726 nelectron, threshold, e_homo, e_lumo, e_mu, &
727 dynamic_threshold, matrix_ks_deviation, &
728 max_iter_lanczos, eps_lanczos, converged, iounit)
731 TYPE(
dbcsr_type),
INTENT(IN) :: matrix_ks, matrix_s_sqrt_inv
732 INTEGER,
INTENT(IN) :: nelectron
733 REAL(kind=
dp),
INTENT(IN) :: threshold
734 REAL(kind=
dp),
INTENT(INOUT) :: e_homo, e_lumo, e_mu
735 LOGICAL,
INTENT(IN),
OPTIONAL :: dynamic_threshold
736 TYPE(
dbcsr_type),
INTENT(INOUT),
OPTIONAL :: matrix_ks_deviation
737 INTEGER,
INTENT(IN) :: max_iter_lanczos
738 REAL(kind=
dp),
INTENT(IN) :: eps_lanczos
739 LOGICAL,
INTENT(OUT),
OPTIONAL :: converged
740 INTEGER,
INTENT(IN),
OPTIONAL :: iounit
742 CHARACTER(LEN=*),
PARAMETER :: routinen =
'density_matrix_trs4'
743 INTEGER,
PARAMETER :: max_iter = 100
744 REAL(kind=
dp),
PARAMETER :: gamma_max = 6.0_dp, gamma_min = 0.0_dp
746 INTEGER :: branch, estimated_steps, handle, i, j, &
748 INTEGER(kind=int_8) :: flop1, flop2
749 LOGICAL :: arnoldi_converged, do_dyn_threshold
750 REAL(kind=
dp) :: current_threshold, delta_n, eps_max, eps_min, est_threshold, frob_id, &
751 frob_x, gam, homo, lumo, max_eig, max_threshold, maxdev, maxev, min_eig, minev, mmin, mu, &
752 mu_a, mu_b, mu_c, mu_fa, mu_fc, occ_matrix, scaled_homo_bound, scaled_lumo_bound, t1, t2, &
753 trace_fx, trace_gx, xi
754 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: gamma_values
756 TYPE(
dbcsr_type) :: matrix_k0, matrix_x, matrix_x_nosym, &
757 matrix_xidsq, matrix_xsq, tmp_gx
759 IF (nelectron == 0)
THEN
764 CALL timeset(routinen, handle)
766 IF (
PRESENT(iounit))
THEN
770 IF (logger%para_env%is_source())
THEN
777 do_dyn_threshold = .false.
778 IF (
PRESENT(dynamic_threshold)) do_dyn_threshold = dynamic_threshold
780 IF (
PRESENT(converged)) converged = .false.
783 CALL dbcsr_create(matrix_x, template=matrix_ks, matrix_type=
"S")
786 CALL dbcsr_create(matrix_x_nosym, template=matrix_ks, matrix_type=dbcsr_type_no_symmetry)
788 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_ks, &
789 0.0_dp, matrix_x_nosym, filter_eps=threshold)
790 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_x_nosym, matrix_s_sqrt_inv, &
791 0.0_dp, matrix_x, filter_eps=threshold)
794 CALL dbcsr_create(matrix_k0, template=matrix_ks, matrix_type=dbcsr_type_no_symmetry)
798 IF (do_dyn_threshold)
THEN
799 cpassert(
PRESENT(matrix_ks_deviation))
800 CALL dbcsr_add(matrix_ks_deviation, matrix_x_nosym, -1.0_dp, 1.0_dp)
801 CALL arnoldi_extremal(matrix_ks_deviation, maxev, minev, max_iter=max_iter_lanczos, threshold=eps_lanczos, &
802 converged=arnoldi_converged)
803 maxdev = max(abs(maxev), abs(minev))
804 IF (unit_nr > 0)
THEN
805 WRITE (unit_nr,
'(T6,A,1X,L12)')
"Lanczos converged: ", arnoldi_converged
806 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
"change in mixed matrix: ", maxdev
807 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
"HOMO upper bound: ", e_homo + maxdev
808 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
"LUMO lower bound: ", e_lumo - maxdev
809 WRITE (unit_nr,
'(T6,A,1X,L12)')
"Predicts a gap ? ", ((e_lumo - maxdev) - (e_homo + maxdev)) > 0
812 CALL dbcsr_copy(matrix_ks_deviation, matrix_x_nosym)
817 CALL arnoldi_extremal(matrix_x_nosym, max_eig, min_eig, max_iter=max_iter_lanczos, threshold=eps_lanczos, &
818 converged=arnoldi_converged)
819 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,1X,2F12.5,1X,A,1X,L1)')
"Est. extremal eigenvalues", &
820 min_eig, max_eig,
" converged: ", arnoldi_converged
825 IF (eps_max == eps_min)
THEN
829 CALL dbcsr_scale(matrix_x, -1.0_dp/(eps_max - eps_min))
832 current_threshold = threshold
833 IF (do_dyn_threshold)
THEN
835 scaled_homo_bound = (eps_max - (e_homo + maxdev))/(eps_max - eps_min)
836 scaled_lumo_bound = (eps_max - (e_lumo - maxdev))/(eps_max - eps_min)
839 CALL dbcsr_create(matrix_xsq, template=matrix_ks, matrix_type=
"S")
841 CALL dbcsr_create(matrix_xidsq, template=matrix_ks, matrix_type=
"S")
843 CALL dbcsr_create(tmp_gx, template=matrix_ks, matrix_type=
"S")
845 ALLOCATE (gamma_values(max_iter))
853 0.0_dp, matrix_xsq, &
854 filter_eps=current_threshold, flop=flop1)
858 CALL dbcsr_add(matrix_xidsq, matrix_xsq, -1.0_dp, 1.0_dp)
865 CALL dbcsr_add(matrix_xidsq, matrix_xsq, -2.0_dp, 1.0_dp)
870 CALL dbcsr_add(tmp_gx, matrix_xsq, 4.0_dp, -3.0_dp)
874 CALL dbcsr_dot(matrix_xsq, matrix_xidsq, trace_gx)
875 CALL dbcsr_dot(matrix_xsq, tmp_gx, trace_fx)
880 delta_n = nelectron - trace_fx
882 IF (((frob_id*frob_id) < (threshold*frob_x*frob_x)) .AND. (abs(delta_n) < 0.5_dp))
THEN
884 ELSE IF (abs(delta_n) < 1e-14_dp)
THEN
888 gam = delta_n/max(trace_gx, abs(delta_n)/100)
890 gamma_values(i) = gam
892 IF (unit_nr > 0 .AND. .false.)
THEN
893 WRITE (unit_nr, *)
"trace_fx", trace_fx,
"trace_gx", trace_gx,
"gam", gam, &
894 "frob_id", frob_id,
"conv", abs(frob_id/frob_x)
897 IF (do_dyn_threshold)
THEN
899 xi = (scaled_homo_bound - scaled_lumo_bound)
900 IF (xi > 0.0_dp)
THEN
901 mmin = 0.5*(scaled_homo_bound + scaled_lumo_bound)
902 max_threshold = abs(1 - 2*mmin)*xi
904 scaled_homo_bound = evaluate_trs4_polynomial(scaled_homo_bound, gamma_values(i:), 1)
905 scaled_lumo_bound = evaluate_trs4_polynomial(scaled_lumo_bound, gamma_values(i:), 1)
906 estimated_steps = estimate_steps(scaled_homo_bound, scaled_lumo_bound, threshold)
908 est_threshold = (threshold/(estimated_steps + i + 1))*xi/(1 + threshold/(estimated_steps + i + 1))
909 est_threshold = min(max_threshold, est_threshold)
910 IF (i > 1) est_threshold = max(est_threshold, 0.1_dp*current_threshold)
911 current_threshold = est_threshold
913 current_threshold = threshold
917 IF (gam > gamma_max)
THEN
919 CALL dbcsr_add(matrix_x, matrix_xsq, 2.0_dp, -1.0_dp)
922 ELSE IF (gam < gamma_min)
THEN
928 CALL dbcsr_add(tmp_gx, matrix_xidsq, 1.0_dp, gam)
931 flop=flop2, filter_eps=current_threshold)
937 IF (unit_nr > 0)
THEN
939 '(T6,A,I3,1X,F10.8,E12.3,F12.3,F13.3,E12.3)')
"TRS4 it ", &
940 i, occ_matrix, abs(trace_gx), t2 - t1, &
941 (flop1 + flop2)/(1.0e6_dp*max(t2 - t1, 0.001_dp)), current_threshold
946 cpabort(
"trace_gx is an abnormal value (NaN/Inf).")
951 IF ((frob_id*frob_id) < (threshold*frob_x*frob_x) .AND. branch == 3 .AND. (abs(delta_n) < 0.5_dp))
THEN
952 IF (
PRESENT(converged)) converged = .true.
959 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,I3,1X,F10.8,E12.3)')
'Final TRS4 iteration ', i, occ_matrix, abs(trace_gx)
967 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_x, matrix_s_sqrt_inv, &
968 0.0_dp, matrix_x_nosym, filter_eps=threshold)
969 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_x_nosym, &
970 0.0_dp, matrix_p, filter_eps=threshold)
974 mu_a = 0.0_dp; mu_b = 1.0_dp;
975 mu_fa = evaluate_trs4_polynomial(mu_a, gamma_values, i - 1) - 0.5_dp
977 mu_c = 0.5*(mu_a + mu_b)
978 mu_fc = evaluate_trs4_polynomial(mu_c, gamma_values, i - 1) - 0.5_dp
979 IF (abs(mu_fc) < 1.0e-6_dp .OR. (mu_b - mu_a)/2 < 1.0e-6_dp)
EXIT
981 IF (mu_fc*mu_fa > 0)
THEN
988 mu = (eps_min - eps_max)*mu_c + eps_max
989 DEALLOCATE (gamma_values)
990 IF (unit_nr > 0)
THEN
991 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
'Chemical potential (mu): ', mu
995 IF (do_dyn_threshold)
THEN
998 threshold, max_iter_lanczos, eps_lanczos, homo, lumo, unit_nr)
1006 CALL timestop(handle)
1028 nelectron, threshold, e_homo, e_lumo, &
1029 non_monotonic, eps_lanczos, max_iter_lanczos, iounit)
1032 TYPE(
dbcsr_type),
INTENT(IN) :: matrix_ks, matrix_s_sqrt_inv
1033 INTEGER,
INTENT(IN) :: nelectron
1034 REAL(kind=
dp),
INTENT(IN) :: threshold
1035 REAL(kind=
dp),
INTENT(INOUT) :: e_homo, e_lumo
1036 LOGICAL,
INTENT(IN),
OPTIONAL :: non_monotonic
1037 REAL(kind=
dp),
INTENT(IN) :: eps_lanczos
1038 INTEGER,
INTENT(IN) :: max_iter_lanczos
1039 INTEGER,
INTENT(IN),
OPTIONAL :: iounit
1041 CHARACTER(LEN=*),
PARAMETER :: routinen =
'density_matrix_tc2'
1042 INTEGER,
PARAMETER :: max_iter = 100
1044 INTEGER :: handle, i, j, k, unit_nr
1045 INTEGER(kind=int_8) :: flop1, flop2
1046 LOGICAL :: converged, do_non_monotonic
1047 REAL(kind=
dp) :: beta, betab, eps_max, eps_min, gama, &
1048 max_eig, min_eig, occ_matrix, t1, t2, &
1050 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: alpha, lambda, nu, poly, wu, x, y
1052 TYPE(
dbcsr_type) :: matrix_tmp, matrix_x, matrix_xsq
1054 CALL timeset(routinen, handle)
1056 IF (
PRESENT(iounit))
THEN
1060 IF (logger%para_env%is_source())
THEN
1067 do_non_monotonic = .false.
1068 IF (
PRESENT(non_monotonic)) do_non_monotonic = non_monotonic
1071 CALL dbcsr_create(matrix_x, template=matrix_ks, matrix_type=dbcsr_type_no_symmetry)
1072 CALL dbcsr_create(matrix_xsq, template=matrix_ks, matrix_type=dbcsr_type_no_symmetry)
1074 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_ks, &
1075 0.0_dp, matrix_xsq, filter_eps=threshold)
1076 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_xsq, matrix_s_sqrt_inv, &
1077 0.0_dp, matrix_x, filter_eps=threshold)
1079 IF (unit_nr > 0)
THEN
1080 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
"HOMO upper bound: ", e_homo
1081 WRITE (unit_nr,
'(T6,A,1X,F12.5)')
"LUMO lower bound: ", e_lumo
1082 WRITE (unit_nr,
'(T6,A,1X,L12)')
"Predicts a gap ? ", ((e_lumo) - (e_homo)) > 0
1086 CALL arnoldi_extremal(matrix_x, max_eig, min_eig, max_iter=max_iter_lanczos, threshold=eps_lanczos, &
1087 converged=converged)
1088 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,1X,2F12.5,1X,A,1X,L1)')
"Est. extremal eigenvalues", &
1089 min_eig, max_eig,
" converged: ", converged
1101 CALL dbcsr_create(matrix_tmp, template=matrix_ks, matrix_type=dbcsr_type_no_symmetry)
1103 ALLOCATE (poly(max_iter))
1104 ALLOCATE (nu(max_iter))
1105 ALLOCATE (wu(max_iter))
1106 ALLOCATE (alpha(max_iter))
1109 ALLOCATE (lambda(4))
1112 beta = (eps_max - abs(e_lumo))/(eps_max - eps_min)
1113 betab = (eps_max + abs(e_homo))/(eps_max - eps_min)
1115 IF ((beta - betab) < 0.005_dp)
THEN
1116 beta = beta - 0.002_dp
1117 betab = betab + 0.002_dp
1120 IF (.NOT. do_non_monotonic)
THEN
1125 IF (e_homo == 0.0_dp)
THEN
1131 trace_fx = nelectron
1136 flop1 = 0; flop2 = 0
1138 IF (abs(trace_fx - nelectron) <= abs(trace_gx - nelectron))
THEN
1141 alpha(i) = 2.0_dp/(2.0_dp - beta)
1146 0.0_dp, matrix_xsq, &
1147 filter_eps=threshold, flop=flop1)
1154 beta = (1.0_dp - alpha(i)) + alpha(i)*beta
1156 betab = (1.0_dp - alpha(i)) + alpha(i)*betab
1161 alpha(i) = 2.0_dp/(1.0_dp + betab)
1165 0.0_dp, matrix_xsq, &
1166 filter_eps=threshold, flop=flop1)
1171 CALL dbcsr_add(matrix_x, matrix_xsq, 2.0_dp, -1.0_dp)
1173 beta = alpha(i)*beta
1174 beta = 2.0_dp*beta - beta*beta
1175 betab = alpha(i)*betab
1176 betab = 2.0_dp*betab - betab*betab
1181 IF (unit_nr > 0)
THEN
1183 '(T6,A,I3,1X,F10.8,E12.3,F12.3,F13.3,E12.3)')
"TC2 it ", &
1184 i, occ_matrix, t2 - t1, &
1185 (flop1 + flop2)/(1.0e6_dp*(t2 - t1)), threshold
1193 CALL dbcsr_add(matrix_xsq, matrix_tmp, -1.0_dp, 1.0_dp)
1199 CALL dbcsr_add(matrix_xsq, matrix_tmp, 1.0_dp, 1.0_dp)
1202 IF (abs(nu(i)) < (threshold))
EXIT
1206 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,I3,1X,1F10.8,1X,1F10.8)')
'Final TC2 iteration ', i, occ_matrix, abs(nu(i))
1209 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_x, matrix_s_sqrt_inv, &
1210 0.0_dp, matrix_tmp, filter_eps=threshold)
1211 CALL dbcsr_multiply(
"N",
"N", 1.0_dp, matrix_s_sqrt_inv, matrix_tmp, &
1212 0.0_dp, matrix_p, filter_eps=threshold)
1222 gama = 6.0_dp - 4.0_dp*(sqrt(2.0_dp))
1223 gama = gama - gama*gama
1224 DO WHILE (nu(i) < gama)
1226 IF (wu(i) < 1.0e-14_dp)
THEN
1230 IF ((1.0_dp - 4.0_dp*nu(i)*nu(i)/wu(i)) < 0.0_dp)
THEN
1234 y(1) = 0.5_dp*(1.0_dp - sqrt(1.0_dp - 4.0_dp*nu(i)*nu(i)/wu(i)))
1235 y(2) = 0.5_dp*(1.0_dp - sqrt(1.0_dp - 4.0_dp*nu(i)))
1236 y(3) = 0.5_dp*(1.0_dp + sqrt(1.0_dp - 4.0_dp*nu(i)))
1237 y(4) = 0.5_dp*(1.0_dp + sqrt(1.0_dp - 4.0_dp*nu(i)*nu(i)/wu(i)))
1238 y(:) = min(1.0_dp, max(0.0_dp, y(:)))
1240 IF (poly(j) == 1.0_dp)
THEN
1243 y(k) = (y(k) - 1.0_dp + alpha(j))/alpha(j)
1247 y(k) = 1.0_dp - sqrt(1.0_dp - y(k))
1248 y(k) = y(k)/alpha(j)
1252 x(1) = min(x(1), y(1))
1253 x(2) = min(x(2), y(2))
1254 x(3) = max(x(3), y(3))
1255 x(4) = max(x(4), y(4))
1260 lambda(k) = eps_max - (eps_max - eps_min)*x(k)
1265 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,3E12.4)')
"outer homo/lumo/gap", e_homo, e_lumo, (e_lumo - e_homo)
1266 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,3E12.4)')
"inner homo/lumo/gap", lambda(3), lambda(2), (lambda(2) - lambda(3))
1277 CALL timestop(handle)
1298 SUBROUTINE compute_homo_lumo(matrix_k, matrix_p, eps_min, eps_max, threshold, max_iter_lanczos, eps_lanczos, homo, lumo, unit_nr)
1300 REAL(kind=
dp) :: eps_min, eps_max, threshold
1301 INTEGER,
INTENT(IN) :: max_iter_lanczos
1302 REAL(kind=
dp),
INTENT(IN) :: eps_lanczos
1303 REAL(kind=
dp) :: homo, lumo
1306 LOGICAL :: converged
1307 REAL(kind=
dp) :: max_eig, min_eig, shift1, shift2
1312 CALL dbcsr_create(tmp1, template=matrix_k, matrix_type=dbcsr_type_no_symmetry)
1314 CALL dbcsr_create(tmp2, template=matrix_k, matrix_type=dbcsr_type_no_symmetry)
1316 CALL dbcsr_create(tmp3, template=matrix_k, matrix_type=dbcsr_type_no_symmetry)
1325 0.0_dp, tmp1, filter_eps=threshold)
1327 threshold=eps_lanczos, max_iter=max_iter_lanczos)
1328 homo = max_eig - shift1
1329 IF (unit_nr > 0)
THEN
1330 WRITE (unit_nr,
'(T6,A,1X,L12)')
"Lanczos converged: ", converged
1340 0.0_dp, tmp1, filter_eps=threshold)
1342 threshold=eps_lanczos, max_iter=max_iter_lanczos)
1343 lumo = -max_eig + shift2
1345 IF (unit_nr > 0)
THEN
1346 WRITE (unit_nr,
'(T6,A,1X,L12)')
"Lanczos converged: ", converged
1347 WRITE (unit_nr,
'(T6,A,1X,3F12.5)')
'HOMO/LUMO/gap', homo, lumo, lumo - homo
1362 FUNCTION evaluate_trs4_polynomial(x, gamma_values, i)
RESULT(xr)
1364 REAL(kind=
dp),
DIMENSION(:) :: gamma_values
1368 REAL(kind=
dp),
PARAMETER :: gam_max = 6.0_dp, gam_min = 0.0_dp
1374 IF (gamma_values(k) > gam_max)
THEN
1376 ELSE IF (gamma_values(k) < gam_min)
THEN
1379 xr = (xr*xr)*(4*xr - 3*xr*xr) + gamma_values(k)*xr*xr*((1 - xr)**2)
1382 END FUNCTION evaluate_trs4_polynomial
1391 FUNCTION estimate_steps(homo, lumo, threshold)
RESULT(steps)
1392 REAL(kind=
dp) :: homo, lumo, threshold
1396 REAL(kind=
dp) :: h, l, m
1402 IF (abs(l) < threshold .AND. abs(1 - h) < threshold)
EXIT
1404 IF (m > 0.5_dp)
THEN
1413 END FUNCTION estimate_steps
arnoldi iteration using dbcsr
subroutine, public arnoldi_extremal(matrix_a, max_ev, min_ev, converged, threshold, max_iter)
simple wrapper to estimate extremal eigenvalues with arnoldi, using the old lanczos interface this hi...
subroutine, public dbcsr_scale(matrix, alpha_scalar)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_desymmetrize(matrix_a, matrix_b)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, first_row, last_row, first_column, last_column, first_k, last_k, retain_sparsity, filter_eps, flop)
...
subroutine, public dbcsr_filter(matrix, eps)
...
real(kind=dp) function, public dbcsr_get_occupation(matrix)
...
subroutine, public dbcsr_finalize(matrix)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_set(matrix, alpha)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_put_block(matrix, row, col, block, summation)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
subroutine, public dbcsr_add_on_diag(matrix, alpha)
Adds the given scalar to the diagonal of the matrix. Reserves any missing diagonal blocks.
subroutine, public dbcsr_trace(matrix, trace)
Computes the trace of the given matrix, also known as the sum of its diagonal elements.
real(dp) function, public dbcsr_frobenius_norm(matrix)
Compute the frobenius norm of a dbcsr matrix.
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
lower level routines for linear scaling SCF
subroutine, public density_matrix_trs4(matrix_p, matrix_ks, matrix_s_sqrt_inv, nelectron, threshold, e_homo, e_lumo, e_mu, dynamic_threshold, matrix_ks_deviation, max_iter_lanczos, eps_lanczos, converged, iounit)
compute the density matrix using a trace-resetting algorithm
subroutine, public ls_scf_init_matrix_s(matrix_s, ls_scf_env)
initialize S matrix related properties (sqrt, inverse...) Might be factored-out since this seems comm...
subroutine, public density_matrix_sign_fixed_mu(matrix_p, trace, mu, sign_method, sign_order, matrix_ks, matrix_s, matrix_s_inv, threshold, sign_symmetric, submatrix_sign_method, matrix_s_sqrt_inv)
for a fixed mu, compute the corresponding density matrix and its trace
subroutine, public apply_matrix_preconditioner(matrix, direction, matrix_bs_sqrt, matrix_bs_sqrt_inv)
apply a preconditioner either forward (precondition) inv(sqrt(bs)) * A * inv(sqrt(bs)) backward (rest...
subroutine, public density_matrix_sign(matrix_p, mu, fixed_mu, sign_method, sign_order, matrix_ks, matrix_s, matrix_s_inv, nelectron, threshold, sign_symmetric, submatrix_sign_method, matrix_s_sqrt_inv)
compute the density matrix with a trace that is close to nelectron. take a mu as input,...
subroutine, public compute_matrix_preconditioner(matrix_s, preconditioner_type, ls_mstruct, matrix_bs_sqrt, matrix_bs_sqrt_inv, threshold, order, eps_lanczos, max_iter_lanczos)
compute for a block positive definite matrix s (bs) the sqrt(bs) and inv(sqrt(bs))
subroutine, public density_matrix_tc2(matrix_p, matrix_ks, matrix_s_sqrt_inv, nelectron, threshold, e_homo, e_lumo, non_monotonic, eps_lanczos, max_iter_lanczos, iounit)
compute the density matrix using a non monotonic trace conserving algorithm based on SIAM DOI....
subroutine, public compute_homo_lumo(matrix_k, matrix_p, eps_min, eps_max, threshold, max_iter_lanczos, eps_lanczos, homo, lumo, unit_nr)
compute the homo and lumo given a KS matrix and a density matrix in the orthonormalized basis and the...
Routines for a linear scaling quickstep SCF run based on the density matrix, with a focus on the inte...
subroutine, public matrix_qs_to_ls(matrix_ls, matrix_qs, ls_mstruct, covariant)
first link to QS, copy a QS matrix to LS matrix used to isolate QS style matrices from LS style will ...
Types needed for a linear scaling quickstep SCF run based on the density matrix.
Routines useful for iterative matrix calculations.
subroutine, public matrix_sign_newton_schulz(matrix_sign, matrix, threshold, sign_order, iounit)
compute the sign a matrix using Newton-Schulz iterations
subroutine, public matrix_sign_submatrix(matrix_sign, matrix, threshold, sign_order, submatrix_sign_method)
Submatrix method.
subroutine, public matrix_sign_submatrix_mu_adjust(matrix_sign, matrix, mu, nelectron, threshold, variant)
Submatrix method with internal adjustment of chemical potential.
subroutine, public invert_hotelling(matrix_inverse, matrix, threshold, use_inv_as_guess, norm_convergence, filter_eps, accelerator_order, max_iter_lanczos, eps_lanczos, silent)
invert a symmetric positive definite matrix by Hotelling's method explicit symmetrization makes this ...
subroutine, public matrix_sign_proot(matrix_sign, matrix, threshold, sign_order)
compute the sign a matrix using the general algorithm for the p-th root of Richters et al....
subroutine, public matrix_sqrt_newton_schulz(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged, iounit)
compute the sqrt of a matrix via the sign function and the corresponding Newton-Schulz iterations the...
subroutine, public matrix_sqrt_proot(matrix_sqrt, matrix_sqrt_inv, matrix, threshold, order, eps_lanczos, max_iter_lanczos, symmetrize, converged)
compute the sqrt of a matrix via the general algorithm for the p-th root of Richters et al....
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Collection of simple mathematical functions and subroutines.
logical function, public abnormal_value(a)
determines if a value is not normal (e.g. for Inf and Nan) based on IO to work also under optimizatio...
type of a logger, at the moment it contains just a print level starting at which level it should be l...