57#include "./base/base_uses.f90"
63 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'bse_full_diag'
90 SUBROUTINE create_a(fm_mat_S_ia_bse, fm_mat_S_bar_ij_bse, fm_mat_S_ab_bse, &
91 fm_A, Eigenval, unit_nr, &
92 homo, virtual, dimen_RI, mp2_env, &
95 TYPE(
cp_fm_type),
INTENT(IN) :: fm_mat_s_ia_bse, fm_mat_s_bar_ij_bse, &
98 REAL(kind=
dp),
DIMENSION(:) :: eigenval
99 INTEGER,
INTENT(IN) :: unit_nr, homo, virtual, dimen_ri
100 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
104 CHARACTER(LEN=*),
PARAMETER :: routinen =
'create_A'
106 INTEGER :: a_virt_row, handle, i_occ_row, &
107 i_row_global, ii, j_col_global, jj, &
108 ncol_local_a, nrow_local_a, sizeeigen
109 INTEGER,
DIMENSION(4) :: reordering
110 INTEGER,
DIMENSION(:),
POINTER :: col_indices_a, row_indices_a
111 REAL(kind=
dp) :: alpha, alpha_screening, eigen_diff
119 CALL timeset(routinen, handle)
121 NULLIFY (dft_control, tddfpt_control)
122 CALL get_qs_env(qs_env, dft_control=dft_control)
123 tddfpt_control => dft_control%tddfpt2_control
125 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
126 WRITE (unit_nr,
'(T2,A10,T13,A10)')
'BSE|DEBUG|',
'Creating A'
130 SELECT CASE (mp2_env%bse%bse_spin_config)
138 alpha_screening = mp2_env%bse%screening_factor
140 alpha_screening = 1.0_dp
152 CALL cp_fm_struct_create(fm_struct_a, context=fm_mat_s_ia_bse%matrix_struct%context, nrow_global=homo*virtual, &
153 ncol_global=homo*virtual, para_env=fm_mat_s_ia_bse%matrix_struct%para_env)
156 IF (tddfpt_control%do_bse_w_only)
THEN
157 CALL cp_fm_create(fm_a_copy, fm_struct_a, name=
"fm_A_iajb")
161 CALL cp_fm_struct_create(fm_struct_w, context=fm_mat_s_ab_bse%matrix_struct%context, nrow_global=homo**2, &
162 ncol_global=virtual**2, para_env=fm_mat_s_ab_bse%matrix_struct%para_env)
169 IF ((.NOT. tddfpt_control%do_bse) .AND. (.NOT. tddfpt_control%do_bse_w_only))
THEN
170 CALL parallel_gemm(transa=
"T", transb=
"N", m=homo*virtual, n=homo*virtual, k=dimen_ri, alpha=alpha, &
171 matrix_a=fm_mat_s_ia_bse, matrix_b=fm_mat_s_ia_bse, beta=0.0_dp, &
175 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
176 WRITE (unit_nr,
'(T2,A10,T13,A16)')
'BSE|DEBUG|',
'Allocated A_iajb'
182 CALL parallel_gemm(transa=
"T", transb=
"N", m=homo**2, n=virtual**2, k=dimen_ri, alpha=alpha_screening, &
183 matrix_a=fm_mat_s_bar_ij_bse, matrix_b=fm_mat_s_ab_bse, beta=0.0_dp, &
187 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
188 WRITE (unit_nr,
'(T2,A10,T13,A16)')
'BSE|DEBUG|',
'Allocated W_ijab'
193 nrow_local=nrow_local_a, &
194 ncol_local=ncol_local_a, &
195 row_indices=row_indices_a, &
196 col_indices=col_indices_a)
203 reordering = [1, 3, 2, 4]
205 virtual, virtual, unit_nr, reordering, mp2_env)
206 IF (tddfpt_control%do_bse_w_only)
THEN
208 virtual, virtual, unit_nr, reordering, mp2_env)
212 IF (tddfpt_control%do_bse .OR. tddfpt_control%do_bse_w_only .OR. &
213 tddfpt_control%do_bse_gw_only)
THEN
216 IF (.NOT. tddfpt_control%do_bse_gw_only)
THEN
217 ALLOCATE (ex_env%bse_w_matrix_MO(1, 1))
218 ALLOCATE (ex_env%bse_a_matrix_MO(1, 1))
219 CALL cp_fm_create(ex_env%bse_w_matrix_MO(1, 1), fm_struct_w)
220 CALL cp_fm_create(ex_env%bse_a_matrix_MO(1, 1), fm_struct_a)
221 CALL cp_fm_to_fm(fm_w, ex_env%bse_w_matrix_MO(1, 1))
222 IF (tddfpt_control%do_bse_w_only)
THEN
223 CALL cp_fm_to_fm(fm_a_copy, ex_env%bse_a_matrix_MO(1, 1))
225 CALL cp_fm_to_fm(fm_a, ex_env%bse_a_matrix_MO(1, 1))
230 IF (tddfpt_control%do_bse_w_only)
CALL cp_fm_release(fm_a_copy)
233 IF (.NOT. tddfpt_control%do_bse)
THEN
234 DO ii = 1, nrow_local_a
236 i_row_global = row_indices_a(ii)
238 DO jj = 1, ncol_local_a
240 j_col_global = col_indices_a(jj)
242 IF (i_row_global == j_col_global)
THEN
243 i_occ_row = (i_row_global - 1)/virtual + 1
244 a_virt_row = mod(i_row_global - 1, virtual) + 1
245 eigen_diff = eigenval(a_virt_row + homo) - eigenval(i_occ_row)
246 fm_a%local_data(ii, jj) = fm_a%local_data(ii, jj) + eigen_diff
253 IF (tddfpt_control%do_bse .OR. tddfpt_control%do_bse_w_only .OR. tddfpt_control%do_bse_gw_only)
THEN
254 sizeeigen =
SIZE(eigenval)
255 ALLOCATE (ex_env%gw_eigen(sizeeigen))
256 ex_env%gw_eigen(:) = eigenval(:)
264 CALL timestop(handle)
283 SUBROUTINE create_b(fm_mat_S_ia_bse, fm_mat_S_bar_ia_bse, fm_B, &
284 homo, virtual, dimen_RI, unit_nr, mp2_env)
286 TYPE(
cp_fm_type),
INTENT(IN) :: fm_mat_s_ia_bse, fm_mat_s_bar_ia_bse
288 INTEGER,
INTENT(IN) :: homo, virtual, dimen_ri, unit_nr
289 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
291 CHARACTER(LEN=*),
PARAMETER :: routinen =
'create_B'
294 INTEGER,
DIMENSION(4) :: reordering
295 REAL(kind=
dp) :: alpha, alpha_screening
299 CALL timeset(routinen, handle)
301 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
302 WRITE (unit_nr,
'(T2,A10,T13,A10)')
'BSE|DEBUG|',
'Creating B'
306 SELECT CASE (mp2_env%bse%bse_spin_config)
314 alpha_screening = mp2_env%bse%screening_factor
316 alpha_screening = 1.0_dp
319 CALL cp_fm_struct_create(fm_struct_v, context=fm_mat_s_ia_bse%matrix_struct%context, nrow_global=homo*virtual, &
320 ncol_global=homo*virtual, para_env=fm_mat_s_ia_bse%matrix_struct%para_env)
327 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
328 WRITE (unit_nr,
'(T2,A10,T13,A16)')
'BSE|DEBUG|',
'Allocated B_iajb'
331 CALL parallel_gemm(transa=
"T", transb=
"N", m=homo*virtual, n=homo*virtual, k=dimen_ri, alpha=alpha, &
332 matrix_a=fm_mat_s_ia_bse, matrix_b=fm_mat_s_ia_bse, beta=0.0_dp, &
338 CALL parallel_gemm(transa=
"T", transb=
"N", m=homo*virtual, n=homo*virtual, k=dimen_ri, alpha=alpha_screening, &
339 matrix_a=fm_mat_s_bar_ia_bse, matrix_b=fm_mat_s_ia_bse, beta=0.0_dp, &
346 reordering = [1, 4, 3, 2]
348 virtual, virtual, unit_nr, reordering, mp2_env)
353 CALL timestop(handle)
374 fm_sqrt_A_minus_B, fm_inv_sqrt_A_minus_B, &
375 homo, virtual, unit_nr, mp2_env, diag_est)
378 TYPE(
cp_fm_type),
INTENT(INOUT) :: fm_c, fm_sqrt_a_minus_b, &
379 fm_inv_sqrt_a_minus_b
380 INTEGER,
INTENT(IN) :: homo, virtual, unit_nr
381 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
382 REAL(kind=
dp),
INTENT(IN) :: diag_est
384 CHARACTER(LEN=*),
PARAMETER :: routinen =
'create_hermitian_form_of_ABBA'
386 INTEGER :: dim_mat, handle, n_dependent
387 REAL(kind=
dp),
DIMENSION(2) :: eigvals_ab_diff
388 TYPE(
cp_fm_type) :: fm_a_minus_b, fm_a_plus_b, fm_dummy, &
391 CALL timeset(routinen, handle)
393 IF (unit_nr > 0)
THEN
394 WRITE (unit_nr,
'(T2,A4,T7,A25,A39,ES6.0,A3)')
'BSE|',
'Diagonalizing aux. matrix', &
395 ' with size of A. This will take around ', diag_est,
" s."
412 CALL cp_fm_create(fm_sqrt_a_minus_b, fm_a%matrix_struct)
414 CALL cp_fm_create(fm_inv_sqrt_a_minus_b, fm_a%matrix_struct)
419 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
420 WRITE (unit_nr,
'(T2,A10,T13,A19)')
'BSE|DEBUG|',
'Created work arrays'
428 CALL cp_fm_to_fm(fm_a_minus_b, fm_inv_sqrt_a_minus_b)
436 CALL cp_fm_power(fm_inv_sqrt_a_minus_b, fm_dummy, -0.5_dp, 0.0_dp, n_dependent, eigvals=eigvals_ab_diff)
440 IF (eigvals_ab_diff(1) < 0)
THEN
441 CALL cp_abort(__location__, &
442 "Matrix (A-B) is not positive definite. "// &
443 "Hermitian diagonalization of full ABBA matrix is ill-defined.")
449 dim_mat = homo*virtual
450 CALL parallel_gemm(
"N",
"N", dim_mat, dim_mat, dim_mat, 1.0_dp, fm_inv_sqrt_a_minus_b, fm_a_minus_b, 0.0_dp, &
454 CALL parallel_gemm(
"N",
"N", dim_mat, dim_mat, dim_mat, 1.0_dp, fm_sqrt_a_minus_b, fm_a_plus_b, 0.0_dp, &
465 CALL parallel_gemm(
"N",
"N", dim_mat, dim_mat, dim_mat, 1.0_dp, fm_work_product, fm_sqrt_a_minus_b, 0.0_dp, &
469 IF (unit_nr > 0 .AND. mp2_env%bse%bse_debug_print)
THEN
470 WRITE (unit_nr,
'(T2,A10,T13,A36)')
'BSE|DEBUG|',
'Filled C=(A-B)^0.5 (A+B) (A-B)^0.5'
473 CALL timestop(handle)
493 fm_sqrt_A_minus_B, fm_inv_sqrt_A_minus_B, &
494 unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
497 INTEGER,
INTENT(IN) :: homo, virtual, homo_irred
498 TYPE(
cp_fm_type),
INTENT(INOUT) :: fm_sqrt_a_minus_b, fm_inv_sqrt_a_minus_b
499 INTEGER,
INTENT(IN) :: unit_nr
500 REAL(kind=
dp),
INTENT(IN) :: diag_est
501 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
503 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: mo_coeff
505 CHARACTER(LEN=*),
PARAMETER :: routinen =
'diagonalize_C'
507 INTEGER :: diag_info, handle
508 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: exc_ens
509 TYPE(
cp_fm_type) :: fm_eigvec_x, fm_eigvec_y, fm_eigvec_z, &
510 fm_mat_eigvec_transform_diff, &
511 fm_mat_eigvec_transform_sum
513 CALL timeset(routinen, handle)
515 IF (unit_nr > 0)
THEN
516 WRITE (unit_nr,
'(T2,A4,T7,A17,A22,ES6.0,A3)')
'BSE|',
'Diagonalizing C. ', &
517 'This will take around ', diag_est,
' s.'
524 ALLOCATE (exc_ens(homo*virtual))
528 IF (diag_info /= 0)
THEN
529 CALL cp_abort(__location__, &
530 "Diagonalization of C=(A-B)^0.5 (A+B) (A-B)^0.5 failed in BSE")
536 IF (exc_ens(1) < 0)
THEN
537 IF (unit_nr > 0)
THEN
538 CALL cp_abort(__location__, &
539 "Matrix C=(A-B)^0.5 (A+B) (A-B)^0.5 has negative eigenvalues, i.e. "// &
540 "(A+B) is not positive definite.")
543 exc_ens = sqrt(exc_ens)
554 CALL cp_fm_create(fm_mat_eigvec_transform_sum, fm_c%matrix_struct)
556 CALL parallel_gemm(transa=
"N", transb=
"N", m=homo*virtual, n=homo*virtual, k=homo*virtual, alpha=1.0_dp, &
557 matrix_a=fm_sqrt_a_minus_b, matrix_b=fm_eigvec_z, beta=0.0_dp, &
558 matrix_c=fm_mat_eigvec_transform_sum)
564 CALL cp_fm_create(fm_mat_eigvec_transform_diff, fm_c%matrix_struct)
566 CALL parallel_gemm(transa=
"N", transb=
"N", m=homo*virtual, n=homo*virtual, k=homo*virtual, alpha=1.0_dp, &
567 matrix_a=fm_inv_sqrt_a_minus_b, matrix_b=fm_eigvec_z, beta=0.0_dp, &
568 matrix_c=fm_mat_eigvec_transform_diff)
578 CALL cp_fm_to_fm(fm_mat_eigvec_transform_sum, fm_eigvec_x)
584 CALL cp_fm_to_fm(fm_mat_eigvec_transform_sum, fm_eigvec_y)
591 CALL postprocess_bse(exc_ens, fm_eigvec_x, mp2_env, qs_env, mo_coeff, &
592 homo, virtual, homo_irred, unit_nr, &
593 .false., fm_eigvec_y)
599 CALL timestop(handle)
616 unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
619 INTEGER,
INTENT(IN) :: homo, virtual, homo_irred, unit_nr
620 REAL(kind=
dp),
INTENT(IN) :: diag_est
621 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
623 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: mo_coeff
625 CHARACTER(LEN=*),
PARAMETER :: routinen =
'diagonalize_A'
627 INTEGER :: diag_info, handle
628 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: exc_ens
631 CALL timeset(routinen, handle)
634 IF (unit_nr > 0)
THEN
635 WRITE (unit_nr,
'(T2,A4,T7,A17,A22,ES6.0,A3)')
'BSE|',
'Diagonalizing A. ', &
636 'This will take around ', diag_est,
' s.'
643 ALLOCATE (exc_ens(homo*virtual))
647 IF (diag_info /= 0)
THEN
648 CALL cp_abort(__location__, &
649 "Diagonalization of A failed in TDA-BSE")
652 CALL postprocess_bse(exc_ens, fm_eigvec, mp2_env, qs_env, mo_coeff, &
653 homo, virtual, homo_irred, unit_nr, .true.)
658 CALL timestop(handle)
676 SUBROUTINE postprocess_bse(Exc_ens, fm_eigvec_X, mp2_env, qs_env, mo_coeff, &
677 homo, virtual, homo_irred, unit_nr, &
678 flag_TDA, fm_eigvec_Y)
680 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: exc_ens
682 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
684 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: mo_coeff
685 INTEGER :: homo, virtual, homo_irred, unit_nr
686 LOGICAL,
OPTIONAL :: flag_tda
687 TYPE(
cp_fm_type),
INTENT(IN),
OPTIONAL :: fm_eigvec_y
689 CHARACTER(LEN=*),
PARAMETER :: routinen =
'postprocess_bse'
691 CHARACTER(LEN=10) :: info_approximation, multiplet
692 INTEGER :: handle, i_exc, idir, n_moments_di, &
694 REAL(kind=
dp) :: alpha
695 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: oscill_str, ref_point_multipole
696 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: polarizability_residues, trans_mom_bse
698 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: fm_dipole_ab_trunc, fm_dipole_ai_trunc, &
699 fm_dipole_ij_trunc, fm_quadpole_ab_trunc, fm_quadpole_ai_trunc, fm_quadpole_ij_trunc
701 DIMENSION(:) :: exc_descr
703 CALL timeset(routinen, handle)
706 IF (mp2_env%bse%bse_spin_config == 0)
THEN
707 multiplet =
"Singlet"
710 multiplet =
"Triplet"
713 IF (.NOT.
PRESENT(flag_tda))
THEN
717 info_approximation =
" -TDA- "
719 info_approximation =
"-ABBA-"
726 ALLOCATE (fm_dipole_ij_trunc(n_moments_di))
727 ALLOCATE (fm_dipole_ab_trunc(n_moments_di))
728 ALLOCATE (fm_dipole_ai_trunc(n_moments_di))
729 ALLOCATE (ref_point_multipole(3))
731 CALL get_multipoles_mo(fm_dipole_ai_trunc, fm_dipole_ij_trunc, fm_dipole_ab_trunc, &
732 qs_env, mo_coeff, ref_point_multipole, 1, &
733 homo, virtual, fm_eigvec_x%matrix_struct%context)
735 IF (mp2_env%bse%num_print_exc_descr > 0)
THEN
737 ALLOCATE (fm_quadpole_ij_trunc(n_moments_quad))
738 ALLOCATE (fm_quadpole_ab_trunc(n_moments_quad))
739 ALLOCATE (fm_quadpole_ai_trunc(n_moments_quad))
740 CALL get_multipoles_mo(fm_quadpole_ai_trunc, fm_quadpole_ij_trunc, fm_quadpole_ab_trunc, &
741 qs_env, mo_coeff, ref_point_multipole, 2, &
742 homo, virtual, fm_eigvec_x%matrix_struct%context)
744 ALLOCATE (exc_descr(mp2_env%bse%num_print_exc_descr))
745 DO i_exc = 1, mp2_env%bse%num_print_exc_descr
747 .false., unit_nr, mp2_env)
748 IF (.NOT. flag_tda)
THEN
750 .false., unit_nr, mp2_env)
753 fm_quadpole_ij_trunc, fm_quadpole_ab_trunc, &
754 fm_quadpole_ai_trunc, &
755 i_exc, homo, virtual, &
759 fm_quadpole_ij_trunc, fm_quadpole_ab_trunc, &
760 fm_quadpole_ai_trunc, &
761 i_exc, homo, virtual)
764 IF (.NOT. flag_tda)
THEN
770 IF (mp2_env%bse%bse_spin_config == 0)
THEN
772 trans_mom_bse, oscill_str, polarizability_residues, &
773 mp2_env, homo, virtual, unit_nr, &
779 multiplet, alpha, mp2_env, unit_nr)
783 info_approximation, mp2_env, unit_nr)
787 info_approximation, mp2_env, unit_nr, fm_eigvec_y)
791 homo, virtual, homo_irred, flag_tda, &
792 info_approximation, mp2_env, unit_nr)
794 IF (mp2_env%bse%num_print_exc_descr > 0)
THEN
796 mp2_env%bse%num_print_exc_descr, mp2_env%bse%bse_debug_print, &
797 mp2_env%bse%print_directional_exc_descr, &
802 IF (mp2_env%bse%do_nto_analysis)
THEN
803 IF (unit_nr > 0)
THEN
804 WRITE (unit_nr,
'(T2,A4)')
'BSE|'
805 WRITE (unit_nr,
'(T2,A4,T7,A47)') &
806 'BSE|',
"Calculating Natural Transition Orbitals (NTOs)."
807 WRITE (unit_nr,
'(T2,A4)')
'BSE|'
810 mo_coeff, homo, virtual, &
811 info_approximation, &
813 qs_env, unit_nr, mp2_env)
816 DO idir = 1, n_moments_di
821 IF (mp2_env%bse%num_print_exc_descr > 0)
THEN
822 DO idir = 1, n_moments_quad
827 DEALLOCATE (fm_quadpole_ai_trunc, fm_quadpole_ij_trunc, fm_quadpole_ab_trunc)
828 DEALLOCATE (exc_descr)
830 DEALLOCATE (fm_dipole_ai_trunc, fm_dipole_ij_trunc, fm_dipole_ab_trunc)
831 DEALLOCATE (ref_point_multipole)
832 IF (mp2_env%bse%bse_spin_config == 0)
THEN
833 DEALLOCATE (oscill_str, trans_mom_bse, polarizability_residues)
835 IF (unit_nr > 0)
THEN
836 WRITE (unit_nr,
'(T2,A4)')
'BSE|'
837 WRITE (unit_nr,
'(T2,A4)')
'BSE|'
840 CALL timestop(handle)
842 END SUBROUTINE postprocess_bse
Routines for the full diagonalization of GW + Bethe-Salpeter for computing electronic excitations.
subroutine, public diagonalize_a(fm_a, homo, virtual, homo_irred, unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
Solving hermitian eigenvalue equation A X^n = Ω^n X^n.
subroutine, public diagonalize_c(fm_c, homo, virtual, homo_irred, fm_sqrt_a_minus_b, fm_inv_sqrt_a_minus_b, unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
Solving eigenvalue equation C Z^n = (Ω^n)^2 Z^n . Here, the eigenvectors Z^n relate to X^n via Eq....
subroutine, public create_b(fm_mat_s_ia_bse, fm_mat_s_bar_ia_bse, fm_b, homo, virtual, dimen_ri, unit_nr, mp2_env)
Matrix B constructed from 3c-B-matrices (cf. subroutine mult_B_with_W) B_ia,jb = α * v_ia,...
subroutine, public create_hermitian_form_of_abba(fm_a, fm_b, fm_c, fm_sqrt_a_minus_b, fm_inv_sqrt_a_minus_b, homo, virtual, unit_nr, mp2_env, diag_est)
Construct Matrix C=(A-B)^0.5 (A+B) (A-B)^0.5 to solve full BSE matrix as a hermitian problem (cf....
subroutine, public create_a(fm_mat_s_ia_bse, fm_mat_s_bar_ij_bse, fm_mat_s_ab_bse, fm_a, eigenval, unit_nr, homo, virtual, dimen_ri, mp2_env, para_env, qs_env)
Matrix A constructed from GW energies and 3c-B-matrices (cf. subroutine mult_B_with_W) A_ia,...
Routines for printing information in context of the BSE calculation.
subroutine, public print_output_header(homo, virtual, homo_irred, flag_tda, multiplet, alpha, mp2_env, unit_nr)
...
subroutine, public print_optical_properties(exc_ens, oscill_str, trans_mom_bse, polarizability_residues, homo, virtual, homo_irred, flag_tda, info_approximation, mp2_env, unit_nr)
...
subroutine, public print_exciton_descriptors(exc_descr, ref_point_multipole, unit_nr, num_print_exc_descr, print_checkvalue, print_directional_exc_descr, prefix_output, qs_env)
Prints exciton descriptors, cf. Mewes et al., JCTC 14, 710-725 (2018)
subroutine, public print_transition_amplitudes(fm_eigvec_x, homo, virtual, homo_irred, info_approximation, mp2_env, unit_nr, fm_eigvec_y)
...
subroutine, public print_excitation_energies(exc_ens, homo, virtual, flag_tda, multiplet, info_approximation, mp2_env, unit_nr)
...
Routines for computing excitonic properties, e.g. exciton diameter, from the BSE.
subroutine, public get_oscillator_strengths(fm_eigvec_x, exc_ens, fm_dipole_ai_trunc, trans_mom_bse, oscill_str, polarizability_residues, mp2_env, homo_red, virtual_red, unit_nr, fm_eigvec_y)
Compute and return BSE dipoles d_r^n = sqrt(2) sum_ia < ψ_i | r | ψ_a > ( X_ia^n + Y_ia^n ) and oscil...
subroutine, public get_exciton_descriptors(exc_descr, fm_x_ia, fm_multipole_ij_trunc, fm_multipole_ab_trunc, fm_multipole_ai_trunc, i_exc, homo, virtual, fm_y_ia)
...
subroutine, public calculate_ntos(fm_x, fm_y, mo_coeff, homo, virtual, info_approximation, oscill_str, qs_env, unit_nr, mp2_env)
...
Auxiliary routines for GW + Bethe-Salpeter for computing electronic excitations.
subroutine, public get_multipoles_mo(fm_multipole_ai_trunc, fm_multipole_ij_trunc, fm_multipole_ab_trunc, qs_env, mo_coeff, rpoint, n_moments, homo_red, virtual_red, context_bse)
...
subroutine, public reshuffle_eigvec(fm_eigvec, fm_eigvec_reshuffled, homo, virtual, n_exc, do_transpose, unit_nr, mp2_env)
...
subroutine, public comp_eigvec_coeff_bse(fm_work, eig_vals, beta, gamma, do_transpose)
Routine for computing the coefficients of the eigenvectors of the BSE matrix from a multiplication wi...
subroutine, public fm_general_add_bse(fm_out, fm_in, beta, nrow_secidx_in, ncol_secidx_in, nrow_secidx_out, ncol_secidx_out, unit_nr, reordering, mp2_env)
Adds and reorders full matrices with a combined index structure, e.g. adding W_ij,...
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
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_scale_and_add(alpha, matrix_a, beta, matrix_b)
calc A <- alpha*A + beta*B optimized for alpha == 1.0 (just add beta*B) and beta == 0....
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
subroutine, public cp_fm_power(matrix, work, exponent, threshold, n_dependent, verbose, eigvals)
...
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...
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
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_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
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
Types for excited states potential energies.
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Defines the basic variable types.
integer, parameter, public dp
Interface to the message passing library MPI.
Types needed for MP2 calculations.
basic linear algebra operations for full matrixes
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
Contains information on the excited states energy.
stores all the informations relevant to an mpi environment