(git:ed6f26b)
Loading...
Searching...
No Matches
bse_main.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Main routines for GW + Bethe-Salpeter for computing electronic excitations
10!> \par History
11!> 04.2024 created [Maximilian Graml]
12! **************************************************************************************************
13
15
16 USE bse_full_diag, ONLY: create_a,&
17 create_b,&
29 USE cp_fm_types, ONLY: cp_fm_release,&
35 USE input_constants, ONLY: bse_abba,&
36 bse_both,&
42 USE kinds, ONLY: dp
44 USE mp2_types, ONLY: mp2_type
46#include "./base/base_uses.f90"
47
48 IMPLICIT NONE
49
50 PRIVATE
51
52 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'bse_main'
53
54 PUBLIC :: start_bse_calculation
55
56CONTAINS
57
58! **************************************************************************************************
59!> \brief Main subroutine managing BSE calculations
60!> \param fm_mat_S_ia_bse ...
61!> \param fm_mat_S_ij_bse ...
62!> \param fm_mat_S_ab_bse ...
63!> \param fm_mat_Q_static_bse_gemm ...
64!> \param Eigenval ...
65!> \param Eigenval_scf ...
66!> \param homo ...
67!> \param virtual ...
68!> \param dimen_RI ...
69!> \param dimen_RI_red ...
70!> \param bse_lev_virt ...
71!> \param gd_array ...
72!> \param color_sub ...
73!> \param mp2_env ...
74!> \param qs_env ...
75!> \param mo_coeff ...
76!> \param unit_nr ...
77! **************************************************************************************************
78 SUBROUTINE start_bse_calculation(fm_mat_S_ia_bse, fm_mat_S_ij_bse, fm_mat_S_ab_bse, &
79 fm_mat_Q_static_bse_gemm, &
80 Eigenval, Eigenval_scf, &
81 homo, virtual, dimen_RI, dimen_RI_red, bse_lev_virt, &
82 gd_array, color_sub, mp2_env, qs_env, mo_coeff, unit_nr)
83
84 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_s_ia_bse, fm_mat_s_ij_bse, &
85 fm_mat_s_ab_bse
86 TYPE(cp_fm_type), INTENT(INOUT) :: fm_mat_q_static_bse_gemm
87 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :), &
88 INTENT(IN) :: eigenval, eigenval_scf
89 INTEGER, DIMENSION(:), INTENT(IN) :: homo, virtual
90 INTEGER, INTENT(IN) :: dimen_ri, dimen_ri_red, bse_lev_virt
91 TYPE(group_dist_d1_type), INTENT(IN) :: gd_array
92 INTEGER, INTENT(IN) :: color_sub
93 TYPE(mp2_type) :: mp2_env
94 TYPE(qs_environment_type), POINTER :: qs_env
95 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
96 INTEGER, INTENT(IN) :: unit_nr
97
98 CHARACTER(LEN=*), PARAMETER :: routinen = 'start_bse_calculation'
99
100 INTEGER :: handle, homo_red, virtual_red
101 LOGICAL :: my_do_abba, my_do_fulldiag, &
102 my_do_iterat_diag, my_do_tda
103 REAL(kind=dp) :: diag_runtime_est
104 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eigenval_reduced
105 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: b_abq_bse_local, b_bar_iaq_bse_local, &
106 b_bar_ijq_bse_local, b_iaq_bse_local
107 TYPE(cp_fm_type) :: fm_a_bse, fm_b_bse, fm_c_bse, fm_inv_sqrt_a_minus_b, fm_mat_s_ab_trunc, &
108 fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, fm_mat_s_ia_trunc, fm_mat_s_ij_trunc, &
109 fm_sqrt_a_minus_b
110 TYPE(cp_logger_type), POINTER :: logger
111 TYPE(mp_para_env_type), POINTER :: para_env
112
113 CALL timeset(routinen, handle)
114
115 para_env => fm_mat_s_ia_bse%matrix_struct%para_env
116
117 my_do_fulldiag = .false.
118 my_do_iterat_diag = .false.
119 my_do_tda = .false.
120 my_do_abba = .false.
121 !Method: Iterative or full diagonalization
122 SELECT CASE (mp2_env%bse%bse_diag_method)
123 CASE (bse_iterdiag)
124 my_do_iterat_diag = .true.
125 !MG: Basics of the Davidson solver are implemented, but not rigorously checked.
126 cpabort("Iterative BSE not yet implemented")
127 CASE (bse_fulldiag)
128 my_do_fulldiag = .true.
129 END SELECT
130 !Approximation: TDA and/or full ABBA matrix
131 SELECT CASE (mp2_env%bse%flag_tda)
132 CASE (bse_tda)
133 my_do_tda = .true.
134 CASE (bse_abba)
135 my_do_abba = .true.
136 CASE (bse_both)
137 my_do_tda = .true.
138 my_do_abba = .true.
139 END SELECT
140
141 CALL print_bse_start_flag(my_do_tda, my_do_abba, unit_nr)
142
143 ! Link BSE debug flag against debug print level
144 logger => cp_get_default_logger()
145 IF (logger%iter_info%print_level == debug_print_level) THEN
146 mp2_env%bse%bse_debug_print = .true.
147 END IF
148
149 CALL fm_mat_s_ia_bse%matrix_struct%para_env%sync()
150 ! We apply the BSE cutoffs using the DFT Eigenenergies
151 ! Reduce matrices in case of energy cutoff for occupied and unoccupied in A/B-BSE-matrices
152 CALL truncate_bse_matrices(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, &
153 fm_mat_s_ia_trunc, fm_mat_s_ij_trunc, fm_mat_s_ab_trunc, &
154 eigenval_scf(:, 1, 1), eigenval(:, 1, 1), eigenval_reduced, &
155 homo(1), virtual(1), dimen_ri, unit_nr, &
156 bse_lev_virt, &
157 homo_red, virtual_red, &
158 mp2_env)
159 ! \bar{B}^P_rs = \sum_R W_PR B^R_rs where B^R_rs = \sum_T [1/sqrt(v)]_RT (T|rs)
160 ! r,s: MO-index, P,R,T: RI-index
161 ! B: fm_mat_S_..., W: fm_mat_Q_...
162 CALL mult_b_with_w(fm_mat_s_ij_trunc, fm_mat_s_ia_trunc, fm_mat_s_bar_ia_bse, &
163 fm_mat_s_bar_ij_bse, fm_mat_q_static_bse_gemm, &
164 dimen_ri_red, homo_red, virtual_red)
165
166 IF (my_do_iterat_diag) THEN
167 CALL fill_local_3c_arrays(fm_mat_s_ab_trunc, fm_mat_s_ia_trunc, &
168 fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, &
169 b_bar_ijq_bse_local, b_abq_bse_local, b_bar_iaq_bse_local, &
170 b_iaq_bse_local, dimen_ri_red, homo_red, virtual_red, &
171 gd_array, color_sub, para_env)
172 END IF
173
174 CALL adapt_bse_input_params(homo_red, virtual_red, unit_nr, mp2_env, qs_env)
175
176 IF (my_do_fulldiag) THEN
177 ! Quick estimate of memory consumption and runtime of diagonalizations
178 CALL estimate_bse_resources(homo_red, virtual_red, unit_nr, my_do_abba, &
179 para_env, diag_runtime_est)
180 ! Matrix A constructed from GW energies and 3c-B-matrices (cf. subroutine mult_B_with_W)
181 ! A_ia,jb = (ε_a-ε_i) δ_ij δ_ab + α * v_ia,jb - W_ij,ab
182 ! ε_a, ε_i are GW singleparticle energies from Eigenval_reduced
183 ! α is a spin-dependent factor
184 ! v_ia,jb = \sum_P B^P_ia B^P_jb (unscreened Coulomb interaction)
185 ! W_ij,ab = \sum_P \bar{B}^P_ij B^P_ab (screened Coulomb interaction)
186
187 ! For unscreened W matrix, we need fm_mat_S_ij_trunc
188 IF (mp2_env%bse%screening_method == bse_screening_tdhf .OR. &
189 mp2_env%bse%screening_method == bse_screening_alpha) THEN
190 CALL create_a(fm_mat_s_ia_trunc, fm_mat_s_ij_trunc, fm_mat_s_ab_trunc, &
191 fm_a_bse, eigenval_reduced, unit_nr, &
192 homo_red, virtual_red, dimen_ri, mp2_env, &
193 para_env)
194 ELSE
195 CALL create_a(fm_mat_s_ia_trunc, fm_mat_s_bar_ij_bse, fm_mat_s_ab_trunc, &
196 fm_a_bse, eigenval_reduced, unit_nr, &
197 homo_red, virtual_red, dimen_ri, mp2_env, &
198 para_env)
199 END IF
200 IF (my_do_abba) THEN
201 ! Matrix B constructed from 3c-B-matrices (cf. subroutine mult_B_with_W)
202 ! B_ia,jb = α * v_ia,jb - W_ib,aj
203 ! α is a spin-dependent factor
204 ! v_ia,jb = \sum_P B^P_ia B^P_jb (unscreened Coulomb interaction)
205 ! W_ib,aj = \sum_P \bar{B}^P_ib B^P_aj (screened Coulomb interaction)
206
207 ! For unscreened W matrix, we need fm_mat_S_ia_trunc
208 IF (mp2_env%bse%screening_method == bse_screening_tdhf .OR. &
209 mp2_env%bse%screening_method == bse_screening_alpha) THEN
210 CALL create_b(fm_mat_s_ia_trunc, fm_mat_s_ia_trunc, fm_b_bse, &
211 homo_red, virtual_red, dimen_ri, unit_nr, mp2_env)
212 ELSE
213 CALL create_b(fm_mat_s_ia_trunc, fm_mat_s_bar_ia_bse, fm_b_bse, &
214 homo_red, virtual_red, dimen_ri, unit_nr, mp2_env)
215 END IF
216 ! Construct Matrix C=(A-B)^0.5 (A+B) (A-B)^0.5 to solve full BSE matrix as a hermitian problem
217 ! (cf. Eq. (A7) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001)).
218 ! We keep fm_sqrt_A_minus_B and fm_inv_sqrt_A_minus_B for print of singleparticle transitions
219 ! of ABBA as described in Eq. (A10) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001).
220 CALL create_hermitian_form_of_abba(fm_a_bse, fm_b_bse, fm_c_bse, &
221 fm_sqrt_a_minus_b, fm_inv_sqrt_a_minus_b, &
222 homo_red, virtual_red, unit_nr, mp2_env, diag_runtime_est)
223 END IF
224 CALL cp_fm_release(fm_b_bse)
225 IF (my_do_tda) THEN
226 ! Solving the hermitian eigenvalue equation A X^n = Ω^n X^n
227 CALL diagonalize_a(fm_a_bse, homo_red, virtual_red, homo(1), &
228 unit_nr, diag_runtime_est, mp2_env, qs_env, mo_coeff)
229 END IF
230 ! Release to avoid faulty use of changed A matrix
231 CALL cp_fm_release(fm_a_bse)
232 IF (my_do_abba) THEN
233 ! Solving eigenvalue equation C Z^n = (Ω^n)^2 Z^n .
234 ! Here, the eigenvectors Z^n relate to X^n via
235 ! Eq. (A10) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001).
236 CALL diagonalize_c(fm_c_bse, homo_red, virtual_red, homo(1), &
237 fm_sqrt_a_minus_b, fm_inv_sqrt_a_minus_b, &
238 unit_nr, diag_runtime_est, mp2_env, qs_env, mo_coeff)
239 END IF
240 ! Release to avoid faulty use of changed C matrix
241 CALL cp_fm_release(fm_c_bse)
242 END IF
243
244 CALL deallocate_matrices_bse(fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, &
245 fm_mat_s_ia_trunc, fm_mat_s_ij_trunc, fm_mat_s_ab_trunc, &
246 fm_mat_q_static_bse_gemm, mp2_env)
247 DEALLOCATE (eigenval_reduced)
248 IF (my_do_iterat_diag) THEN
249 ! Contains untested Block-Davidson algorithm
250 CALL do_subspace_iterations(b_bar_ijq_bse_local, b_abq_bse_local, b_bar_iaq_bse_local, &
251 b_iaq_bse_local, homo(1), virtual(1), mp2_env%bse%bse_spin_config, unit_nr, &
252 eigenval(:, 1, 1), para_env, mp2_env)
253 ! Deallocate local 3c-B-matrices
254 DEALLOCATE (b_bar_ijq_bse_local, b_abq_bse_local, b_bar_iaq_bse_local, b_iaq_bse_local)
255 END IF
256
257 IF (unit_nr > 0) THEN
258 WRITE (unit_nr, '(T2,A4,T7,A53)') 'BSE|', 'The BSE was successfully calculated. Have a nice day!'
259 END IF
260
261 CALL timestop(handle)
262
263 END SUBROUTINE start_bse_calculation
264
265END MODULE bse_main
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_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)
Matrix A constructed from GW energies and 3c-B-matrices (cf. subroutine mult_B_with_W) A_ia,...
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....
Iterative routines for GW + Bethe-Salpeter for computing electronic excitations.
subroutine, public do_subspace_iterations(b_bar_ijq_bse_local, b_abq_bse_local, b_bar_iaq_bse_local, b_iaq_bse_local, homo, virtual, bse_spin_config, unit_nr, eigenval, para_env, mp2_env)
...
subroutine, public fill_local_3c_arrays(fm_mat_s_ab_bse, fm_mat_s, fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, b_bar_ijq_bse_local, b_abq_bse_local, b_bar_iaq_bse_local, b_iaq_bse_local, dimen_ri, homo, virtual, gd_array, color_sub, para_env)
...
Main routines for GW + Bethe-Salpeter for computing electronic excitations.
Definition bse_main.F:14
subroutine, public start_bse_calculation(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, fm_mat_q_static_bse_gemm, eigenval, eigenval_scf, homo, virtual, dimen_ri, dimen_ri_red, bse_lev_virt, gd_array, color_sub, mp2_env, qs_env, mo_coeff, unit_nr)
Main subroutine managing BSE calculations.
Definition bse_main.F:83
Routines for printing information in context of the BSE calculation.
Definition bse_print.F:13
subroutine, public print_bse_start_flag(bse_tda, bse_abba, unit_nr)
...
Definition bse_print.F:53
Auxiliary routines for GW + Bethe-Salpeter for computing electronic excitations.
Definition bse_util.F:13
subroutine, public estimate_bse_resources(homo_red, virtual_red, unit_nr, bse_abba, para_env, diag_runtime_est)
Roughly estimates the needed runtime and memory during the BSE run.
Definition bse_util.F:936
subroutine, public deallocate_matrices_bse(fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, fm_mat_s_trunc, fm_mat_s_ij_trunc, fm_mat_s_ab_trunc, fm_mat_q_static_bse_gemm, mp2_env)
...
Definition bse_util.F:752
subroutine, public mult_b_with_w(fm_mat_s_ij_bse, fm_mat_s_ia_bse, fm_mat_s_bar_ia_bse, fm_mat_s_bar_ij_bse, fm_mat_q_static_bse_gemm, dimen_ri, homo, virtual)
Multiplies B-matrix (RI-3c-Integrals) with W (screening) to obtain \bar{B}.
Definition bse_util.F:113
subroutine, public adapt_bse_input_params(homo, virtual, unit_nr, mp2_env, qs_env)
Checks BSE input section and adapts them if necessary.
Definition bse_util.F:1520
subroutine, public truncate_bse_matrices(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, fm_mat_s_trunc, fm_mat_s_ij_trunc, fm_mat_s_ab_trunc, eigenval_scf, eigenval, eigenval_reduced, homo, virtual, dimen_ri, unit_nr, bse_lev_virt, homo_red, virt_red, mp2_env)
Determines indices within the given energy cutoffs and truncates Eigenvalues and matrices.
Definition bse_util.F:1211
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
Types to describe group distributions.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public bse_screening_tdhf
integer, parameter, public bse_iterdiag
integer, parameter, public bse_fulldiag
integer, parameter, public bse_tda
integer, parameter, public bse_both
integer, parameter, public bse_screening_alpha
integer, parameter, public bse_abba
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Interface to the message passing library MPI.
Types needed for MP2 calculations.
Definition mp2_types.F:14
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment