(git:03cdb6f)
Loading...
Searching...
No Matches
bse_full_diag.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines for the full diagonalization of GW + Bethe-Salpeter for computing
10!> electronic excitations
11!> \par History
12!> 10.2023 created [Maximilian Graml]
13! **************************************************************************************************
15
40 USE cp_fm_types, ONLY: cp_fm_create,&
51 USE kinds, ONLY: dp
53 USE mp2_types, ONLY: mp2_type
57#include "./base/base_uses.f90"
58
59 IMPLICIT NONE
60
61 PRIVATE
62
63 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'bse_full_diag'
64
67
68CONTAINS
69
70! **************************************************************************************************
71!> \brief Matrix A constructed from GW energies and 3c-B-matrices (cf. subroutine mult_B_with_W)
72!> A_ia,jb = (ε_a-ε_i) δ_ij δ_ab + α * v_ia,jb - W_ij,ab
73!> ε_a, ε_i are GW singleparticle energies from Eigenval_reduced
74!> α is a spin-dependent factor
75!> v_ia,jb = \sum_P B^P_ia B^P_jb (unscreened Coulomb interaction)
76!> W_ij,ab = \sum_P \bar{B}^P_ij B^P_ab (screened Coulomb interaction)
77!> \param fm_mat_S_ia_bse ...
78!> \param fm_mat_S_bar_ij_bse ...
79!> \param fm_mat_S_ab_bse ...
80!> \param fm_A ...
81!> \param Eigenval ...
82!> \param unit_nr ...
83!> \param homo ...
84!> \param virtual ...
85!> \param dimen_RI ...
86!> \param mp2_env ...
87!> \param para_env ...
88!> \param qs_env ...
89! **************************************************************************************************
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, &
93 para_env, qs_env)
94
95 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_s_ia_bse, fm_mat_s_bar_ij_bse, &
96 fm_mat_s_ab_bse
97 TYPE(cp_fm_type), INTENT(INOUT) :: fm_a
98 REAL(kind=dp), DIMENSION(:) :: eigenval
99 INTEGER, INTENT(IN) :: unit_nr, homo, virtual, dimen_ri
100 TYPE(mp2_type), INTENT(INOUT) :: mp2_env
101 TYPE(mp_para_env_type), INTENT(INOUT) :: para_env
102 TYPE(qs_environment_type), POINTER :: qs_env
103
104 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_A'
105
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
112 TYPE(cp_blacs_env_type), POINTER :: blacs_env
113 TYPE(cp_fm_struct_type), POINTER :: fm_struct_a, fm_struct_w
114 TYPE(cp_fm_type) :: fm_a_copy, fm_w
115 TYPE(dft_control_type), POINTER :: dft_control
116 TYPE(excited_energy_type), POINTER :: ex_env
117 TYPE(tddfpt2_control_type), POINTER :: tddfpt_control
118
119 CALL timeset(routinen, handle)
120
121 NULLIFY (dft_control, tddfpt_control)
122 CALL get_qs_env(qs_env, dft_control=dft_control)
123 tddfpt_control => dft_control%tddfpt2_control
124
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'
127 END IF
128
129 !Determines factor of exchange term, depending on requested spin configuration (cf. input_constants.F)
130 SELECT CASE (mp2_env%bse%bse_spin_config)
131 CASE (bse_singlet)
132 alpha = 2.0_dp
133 CASE (bse_triplet)
134 alpha = 0.0_dp
135 END SELECT
136
137 IF (mp2_env%bse%screening_method == bse_screening_alpha) THEN
138 alpha_screening = mp2_env%bse%screening_factor
139 ELSE
140 alpha_screening = 1.0_dp
141 END IF
142
143 ! create the blacs env for ij matrices (NOT fm_mat_S_ia_bse%matrix_struct related parallel_gemms!)
144 NULLIFY (blacs_env)
145 CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=para_env)
146
147 ! We have to use the same blacs_env for A as for the matrices fm_mat_S_ia_bse from RPA
148 ! Logic: A_ia,jb = (ε_a-ε_i) δ_ij δ_ab + α * v_ia,jb - W_ij,ab
149 ! We create v_ia,jb and W_ij,ab, then we communicate entries from local W_ij,ab
150 ! to the full matrix v_ia,jb. By adding these and the energy diffenences: v_ia,jb -> A_ia,jb
151 ! We use the A matrix already from the start instead of v
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)
154 CALL cp_fm_create(fm_a, fm_struct_a, name="fm_A_iajb")
155 CALL cp_fm_set_all(fm_a, 0.0_dp)
156 IF (tddfpt_control%do_bse_w_only) THEN
157 CALL cp_fm_create(fm_a_copy, fm_struct_a, name="fm_A_iajb")
158 CALL cp_fm_set_all(fm_a_copy, 0.0_dp)
159 END IF
160
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)
163 CALL cp_fm_create(fm_w, fm_struct_w, name="fm_W_ijab")
164 CALL cp_fm_set_all(fm_w, 0.0_dp)
165
166 ! Create A matrix from GW Energies, v_ia,jb and W_ij,ab (different blacs_env!)
167 ! v_ia,jb, which is directly initialized in A (with a factor of alpha)
168 ! v_ia,jb = \sum_P B^P_ia B^P_jb
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, &
172 matrix_c=fm_a)
173 END IF
174
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'
177 END IF
178
179 ! If infinite screening is applied, fm_W is simply 0 - Otherwise it needs to be computed from 3c integrals
180 IF (mp2_env%bse%screening_method /= bse_screening_rpa) THEN
181 !W_ij,ab = \sum_P \bar{B}^P_ij B^P_ab
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, &
184 matrix_c=fm_w)
185 END IF
186
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'
189 END IF
190
191 ! We start by moving data from local parts of W_ij,ab to the full matrix A_ia,jb using buffers
192 CALL cp_fm_get_info(matrix=fm_a, &
193 nrow_local=nrow_local_a, &
194 ncol_local=ncol_local_a, &
195 row_indices=row_indices_a, &
196 col_indices=col_indices_a)
197 ! Writing -1.0_dp * W_ij,ab to A_ia,jb, i.e. beta = -1.0_dp,
198 ! W_ij,ab: nrow_secidx_in = homo, ncol_secidx_in = virtual
199 ! A_ia,jb: nrow_secidx_out = virtual, ncol_secidx_out = virtual
200
201 ! If infinite screening is applied, fm_W is simply 0 - Otherwise it needs to be computed from 3c integrals
202 IF (mp2_env%bse%screening_method /= bse_screening_rpa) THEN
203 reordering = [1, 3, 2, 4]
204 CALL fm_general_add_bse(fm_a, fm_w, -1.0_dp, homo, virtual, &
205 virtual, virtual, unit_nr, reordering, mp2_env)
206 IF (tddfpt_control%do_bse_w_only) THEN
207 CALL fm_general_add_bse(fm_a_copy, fm_w, -1.0_dp, homo, virtual, &
208 virtual, virtual, unit_nr, reordering, mp2_env)
209 END IF
210 END IF
211 !full matrix W is not needed anymore, release it to save memory
212 IF (tddfpt_control%do_bse .OR. tddfpt_control%do_bse_w_only .OR. &
213 tddfpt_control%do_bse_gw_only) THEN
214 NULLIFY (ex_env)
215 CALL get_qs_env(qs_env, exstate_env=ex_env)
216 IF (.NOT. tddfpt_control%do_bse_gw_only) THEN
217 ALLOCATE (ex_env%bse_w_matrix_MO(1, 1)) ! for now only closed-shell
218 ALLOCATE (ex_env%bse_a_matrix_MO(1, 1)) ! for now only closed-shell
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))
224 ELSE
225 CALL cp_fm_to_fm(fm_a, ex_env%bse_a_matrix_MO(1, 1))
226 END IF
227 END IF
228 END IF
229 CALL cp_fm_release(fm_w)
230 IF (tddfpt_control%do_bse_w_only) CALL cp_fm_release(fm_a_copy)
231
232 !Now add the energy differences (ε_a-ε_i) on the diagonal (i.e. δ_ij δ_ab) of A_ia,jb
233 IF (.NOT. tddfpt_control%do_bse) THEN
234 DO ii = 1, nrow_local_a
235
236 i_row_global = row_indices_a(ii)
237
238 DO jj = 1, ncol_local_a
239
240 j_col_global = col_indices_a(jj)
241
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
247
248 END IF
249 END DO
250 END DO
251 END IF
252
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)) ! for now only closed-shell
256 ex_env%gw_eigen(:) = eigenval(:)
257 END IF
258
259 CALL cp_fm_struct_release(fm_struct_a)
260 CALL cp_fm_struct_release(fm_struct_w)
261
262 CALL cp_blacs_env_release(blacs_env)
263
264 CALL timestop(handle)
265
266 END SUBROUTINE create_a
267
268! **************************************************************************************************
269!> \brief Matrix B constructed from 3c-B-matrices (cf. subroutine mult_B_with_W)
270!> B_ia,jb = α * v_ia,jb - W_ib,aj
271!> α is a spin-dependent factor
272!> v_ia,jb = \sum_P B^P_ia B^P_jb (unscreened Coulomb interaction)
273!> W_ib,aj = \sum_P \bar{B}^P_ib B^P_aj (screened Coulomb interaction)
274!> \param fm_mat_S_ia_bse ...
275!> \param fm_mat_S_bar_ia_bse ...
276!> \param fm_B ...
277!> \param homo ...
278!> \param virtual ...
279!> \param dimen_RI ...
280!> \param unit_nr ...
281!> \param mp2_env ...
282! **************************************************************************************************
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)
285
286 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_s_ia_bse, fm_mat_s_bar_ia_bse
287 TYPE(cp_fm_type), INTENT(INOUT) :: fm_b
288 INTEGER, INTENT(IN) :: homo, virtual, dimen_ri, unit_nr
289 TYPE(mp2_type), INTENT(INOUT) :: mp2_env
290
291 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_B'
292
293 INTEGER :: handle
294 INTEGER, DIMENSION(4) :: reordering
295 REAL(kind=dp) :: alpha, alpha_screening
296 TYPE(cp_fm_struct_type), POINTER :: fm_struct_v
297 TYPE(cp_fm_type) :: fm_w
298
299 CALL timeset(routinen, handle)
300
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'
303 END IF
304
305 ! Determines factor of exchange term, depending on requested spin configuration (cf. input_constants.F)
306 SELECT CASE (mp2_env%bse%bse_spin_config)
307 CASE (bse_singlet)
308 alpha = 2.0_dp
309 CASE (bse_triplet)
310 alpha = 0.0_dp
311 END SELECT
312
313 IF (mp2_env%bse%screening_method == bse_screening_alpha) THEN
314 alpha_screening = mp2_env%bse%screening_factor
315 ELSE
316 alpha_screening = 1.0_dp
317 END IF
318
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)
321 CALL cp_fm_create(fm_b, fm_struct_v, name="fm_B_iajb")
322 CALL cp_fm_set_all(fm_b, 0.0_dp)
323
324 CALL cp_fm_create(fm_w, fm_struct_v, name="fm_W_ibaj")
325 CALL cp_fm_set_all(fm_w, 0.0_dp)
326
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'
329 END IF
330 ! v_ia,jb = \sum_P B^P_ia B^P_jb
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, &
333 matrix_c=fm_b)
334
335 ! If infinite screening is applied, fm_W is simply 0 - Otherwise it needs to be computed from 3c integrals
336 IF (mp2_env%bse%screening_method /= bse_screening_rpa) THEN
337 ! W_ib,aj = \sum_P \bar{B}^P_ib B^P_aj
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, &
340 matrix_c=fm_w)
341
342 ! from W_ib,ja to A_ia,jb (formally: W_ib,aj, but our internal indexorder is different)
343 ! Writing -1.0_dp * W_ib,ja to A_ia,jb, i.e. beta = -1.0_dp,
344 ! W_ib,ja: nrow_secidx_in = virtual, ncol_secidx_in = virtual
345 ! A_ia,jb: nrow_secidx_out = virtual, ncol_secidx_out = virtual
346 reordering = [1, 4, 3, 2]
347 CALL fm_general_add_bse(fm_b, fm_w, -1.0_dp, virtual, virtual, &
348 virtual, virtual, unit_nr, reordering, mp2_env)
349 END IF
350
351 CALL cp_fm_release(fm_w)
352 CALL cp_fm_struct_release(fm_struct_v)
353 CALL timestop(handle)
354
355 END SUBROUTINE create_b
356
357 ! **************************************************************************************************
358!> \brief Construct Matrix C=(A-B)^0.5 (A+B) (A-B)^0.5 to solve full BSE matrix as a hermitian problem
359!> (cf. Eq. (A7) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001)).
360!> We keep fm_sqrt_A_minus_B and fm_inv_sqrt_A_minus_B for print of singleparticle transitions
361!> of ABBA as described in Eq. (A10) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001).
362!> \param fm_A ...
363!> \param fm_B ...
364!> \param fm_C ...
365!> \param fm_sqrt_A_minus_B ...
366!> \param fm_inv_sqrt_A_minus_B ...
367!> \param homo ...
368!> \param virtual ...
369!> \param unit_nr ...
370!> \param mp2_env ...
371!> \param diag_est ...
372! **************************************************************************************************
373 SUBROUTINE create_hermitian_form_of_abba(fm_A, fm_B, fm_C, &
374 fm_sqrt_A_minus_B, fm_inv_sqrt_A_minus_B, &
375 homo, virtual, unit_nr, mp2_env, diag_est)
376
377 TYPE(cp_fm_type), INTENT(IN) :: fm_a, fm_b
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
383
384 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_hermitian_form_of_ABBA'
385
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, &
389 fm_work_product
390
391 CALL timeset(routinen, handle)
392
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."
396 END IF
397
398 ! Create work matrices, which will hold A+B and A-B and their powers
399 ! C is created afterwards to save memory
400 ! Final result: C = (A-B)^0.5 (A+B) (A-B)^0.5 EQ.I
401 ! \_______/ \___/ \______/
402 ! fm_sqrt_A_minus_B fm_A_plus_B fm_sqrt_A_minus_B
403 ! (EQ.Ia) (EQ.Ib) (EQ.Ia)
404 ! Intermediate work matrices:
405 ! fm_inv_sqrt_A_minus_B: (A-B)^-0.5 EQ.II
406 ! fm_A_minus_B: (A-B) EQ.III
407 ! fm_work_product: (A-B)^0.5 (A+B) from (EQ.Ia) and (EQ.Ib) EQ.IV
408 CALL cp_fm_create(fm_a_plus_b, fm_a%matrix_struct)
409 CALL cp_fm_to_fm(fm_a, fm_a_plus_b)
410 CALL cp_fm_create(fm_a_minus_b, fm_a%matrix_struct)
411 CALL cp_fm_to_fm(fm_a, fm_a_minus_b)
412 CALL cp_fm_create(fm_sqrt_a_minus_b, fm_a%matrix_struct)
413 CALL cp_fm_set_all(fm_sqrt_a_minus_b, 0.0_dp)
414 CALL cp_fm_create(fm_inv_sqrt_a_minus_b, fm_a%matrix_struct)
415 CALL cp_fm_set_all(fm_inv_sqrt_a_minus_b, 0.0_dp)
416
417 CALL cp_fm_create(fm_work_product, fm_a%matrix_struct)
418
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'
421 END IF
422
423 ! Add/Substract B (cf. EQs. Ib and III)
424 CALL cp_fm_scale_and_add(1.0_dp, fm_a_plus_b, 1.0_dp, fm_b)
425 CALL cp_fm_scale_and_add(1.0_dp, fm_a_minus_b, -1.0_dp, fm_b)
426
427 ! cp_fm_power will overwrite matrix, therefore we create copies
428 CALL cp_fm_to_fm(fm_a_minus_b, fm_inv_sqrt_a_minus_b)
429
430 ! In order to avoid a second diagonalization (cp_fm_power), we create (A-B)^0.5 (EQ.Ia)
431 ! from (A-B)^-0.5 (EQ.II) by multiplication with (A-B) (EQ.III) afterwards.
432
433 ! Raise A-B to -0.5_dp, no quenching of eigenvectors, hence threshold=0.0_dp
434 CALL cp_fm_create(fm_dummy, fm_a%matrix_struct)
435 ! Create (A-B)^-0.5 (cf. EQ.II)
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)
437 CALL cp_fm_release(fm_dummy)
438 ! Raise an error in case the the matrix A-B is not positive definite (i.e. negative eigenvalues)
439 ! In this case, the procedure for hermitian form of ABBA is not applicable
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.")
444 END IF
445
446 ! We keep fm_inv_sqrt_A_minus_B for print of singleparticle transitions of ABBA
447 ! We further create (A-B)^0.5 for the singleparticle transitions of ABBA
448 ! Create (A-B)^0.5= (A-B)^-0.5 * (A-B) (EQ.Ia)
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, &
451 fm_sqrt_a_minus_b)
452
453 ! Compute and store LHS of C, i.e. (A-B)^0.5 (A+B) (EQ.IV)
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, &
455 fm_work_product)
456
457 ! Release to save memory
458 CALL cp_fm_release(fm_a_plus_b)
459 CALL cp_fm_release(fm_a_minus_b)
460
461 ! Now create full
462 CALL cp_fm_create(fm_c, fm_a%matrix_struct)
463 CALL cp_fm_set_all(fm_c, 0.0_dp)
464 ! Compute C=(A-B)^0.5 (A+B) (A-B)^0.5 (EQ.I)
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, &
466 fm_c)
467 CALL cp_fm_release(fm_work_product)
468
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'
471 END IF
472
473 CALL timestop(handle)
474 END SUBROUTINE create_hermitian_form_of_abba
475
476! **************************************************************************************************
477!> \brief Solving eigenvalue equation C Z^n = (Ω^n)^2 Z^n .
478!> Here, the eigenvectors Z^n relate to X^n via
479!> Eq. (A10) in F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001).
480!> \param fm_C ...
481!> \param homo ...
482!> \param virtual ...
483!> \param homo_irred ...
484!> \param fm_sqrt_A_minus_B ...
485!> \param fm_inv_sqrt_A_minus_B ...
486!> \param unit_nr ...
487!> \param diag_est ...
488!> \param mp2_env ...
489!> \param qs_env ...
490!> \param mo_coeff ...
491! **************************************************************************************************
492 SUBROUTINE diagonalize_c(fm_C, homo, virtual, homo_irred, &
493 fm_sqrt_A_minus_B, fm_inv_sqrt_A_minus_B, &
494 unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
495
496 TYPE(cp_fm_type), INTENT(INOUT) :: fm_c
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
502 TYPE(qs_environment_type), POINTER :: qs_env
503 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
504
505 CHARACTER(LEN=*), PARAMETER :: routinen = 'diagonalize_C'
506
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
512
513 CALL timeset(routinen, handle)
514
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.'
518 END IF
519
520 !We have now the full matrix C=(A-B)^0.5 (A+B) (A-B)^0.5
521 !Now: Diagonalize it
522 CALL cp_fm_create(fm_eigvec_z, fm_c%matrix_struct)
523
524 ALLOCATE (exc_ens(homo*virtual))
525
526 CALL choose_eigv_solver(fm_c, fm_eigvec_z, exc_ens, diag_info)
527
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")
531 END IF
532
533 ! C could have negative eigenvalues, since we do not explicitly check A+B
534 ! for positive definiteness (would make another O(N^6) Diagon. necessary)
535 ! Instead, we include a check here
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.")
541 END IF
542 END IF
543 exc_ens = sqrt(exc_ens)
544
545 ! Prepare eigenvector for interpretation of singleparticle transitions
546 ! Compare: F. Furche J. Chem. Phys., Vol. 114, No. 14, (2001)
547 ! We aim for the upper part of the vector (X,Y) for a direct comparison with the TDA result
548
549 ! Following Furche, we basically use Eqs. (A10): First, we multiply
550 ! the (A-B)^+-0.5 with eigenvectors and then the eigenvalues
551 ! One has to be careful about the index structure, since the eigenvector matrix is not symmetric anymore!
552
553 ! First, Eq. I from (A10) from Furche: (X+Y)_n = (Ω_n)^-0.5 (A-B)^0.5 T_n
554 CALL cp_fm_create(fm_mat_eigvec_transform_sum, fm_c%matrix_struct)
555 CALL cp_fm_set_all(fm_mat_eigvec_transform_sum, 0.0_dp)
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)
559 CALL cp_fm_release(fm_sqrt_a_minus_b)
560 ! This normalizes the eigenvectors
561 CALL comp_eigvec_coeff_bse(fm_mat_eigvec_transform_sum, exc_ens, -0.5_dp, gamma=2.0_dp, do_transpose=.true.)
562
563 ! Second, Eq. II from (A10) from Furche: (X-Y)_n = (Ω_n)^0.5 (A-B)^-0.5 T_n
564 CALL cp_fm_create(fm_mat_eigvec_transform_diff, fm_c%matrix_struct)
565 CALL cp_fm_set_all(fm_mat_eigvec_transform_diff, 0.0_dp)
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)
569 CALL cp_fm_release(fm_inv_sqrt_a_minus_b)
570 CALL cp_fm_release(fm_eigvec_z)
571
572 ! This normalizes the eigenvectors
573 CALL comp_eigvec_coeff_bse(fm_mat_eigvec_transform_diff, exc_ens, 0.5_dp, gamma=2.0_dp, do_transpose=.true.)
574
575 ! Now, we add the two equations to obtain X_n
576 ! Add overwrites the first argument, therefore we copy it beforehand
577 CALL cp_fm_create(fm_eigvec_x, fm_c%matrix_struct)
578 CALL cp_fm_to_fm(fm_mat_eigvec_transform_sum, fm_eigvec_x)
579 CALL cp_fm_scale_and_add(1.0_dp, fm_eigvec_x, 1.0_dp, fm_mat_eigvec_transform_diff)
580
581 ! Now, we subtract the two equations to obtain Y_n
582 ! Add overwrites the first argument, therefore we copy it beforehand
583 CALL cp_fm_create(fm_eigvec_y, fm_c%matrix_struct)
584 CALL cp_fm_to_fm(fm_mat_eigvec_transform_sum, fm_eigvec_y)
585 CALL cp_fm_scale_and_add(1.0_dp, fm_eigvec_y, -1.0_dp, fm_mat_eigvec_transform_diff)
586
587 !Cleanup
588 CALL cp_fm_release(fm_mat_eigvec_transform_diff)
589 CALL cp_fm_release(fm_mat_eigvec_transform_sum)
590
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)
594
595 DEALLOCATE (exc_ens)
596 CALL cp_fm_release(fm_eigvec_x)
597 CALL cp_fm_release(fm_eigvec_y)
598
599 CALL timestop(handle)
600
601 END SUBROUTINE diagonalize_c
602
603! **************************************************************************************************
604!> \brief Solving hermitian eigenvalue equation A X^n = Ω^n X^n
605!> \param fm_A ...
606!> \param homo ...
607!> \param virtual ...
608!> \param homo_irred ...
609!> \param unit_nr ...
610!> \param diag_est ...
611!> \param mp2_env ...
612!> \param qs_env ...
613!> \param mo_coeff ...
614! **************************************************************************************************
615 SUBROUTINE diagonalize_a(fm_A, homo, virtual, homo_irred, &
616 unit_nr, diag_est, mp2_env, qs_env, mo_coeff)
617
618 TYPE(cp_fm_type), INTENT(INOUT) :: fm_a
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
622 TYPE(qs_environment_type), POINTER :: qs_env
623 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
624
625 CHARACTER(LEN=*), PARAMETER :: routinen = 'diagonalize_A'
626
627 INTEGER :: diag_info, handle
628 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: exc_ens
629 TYPE(cp_fm_type) :: fm_eigvec
630
631 CALL timeset(routinen, handle)
632
633 !Continue with formatting of subroutine create_A
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.'
637 END IF
638
639 !We have now the full matrix A_iajb, distributed over all ranks
640 !Now: Diagonalize it
641 CALL cp_fm_create(fm_eigvec, fm_a%matrix_struct)
642
643 ALLOCATE (exc_ens(homo*virtual))
644
645 CALL choose_eigv_solver(fm_a, fm_eigvec, exc_ens, diag_info)
646
647 IF (diag_info /= 0) THEN
648 CALL cp_abort(__location__, &
649 "Diagonalization of A failed in TDA-BSE")
650 END IF
651
652 CALL postprocess_bse(exc_ens, fm_eigvec, mp2_env, qs_env, mo_coeff, &
653 homo, virtual, homo_irred, unit_nr, .true.)
654
655 CALL cp_fm_release(fm_eigvec)
656 DEALLOCATE (exc_ens)
657
658 CALL timestop(handle)
659
660 END SUBROUTINE diagonalize_a
661
662! **************************************************************************************************
663!> \brief Prints the success message (incl. energies) for full diag of BSE (TDA/full ABBA via flag)
664!> \param Exc_ens ...
665!> \param fm_eigvec_X ...
666!> \param mp2_env ...
667!> \param qs_env ...
668!> \param mo_coeff ...
669!> \param homo ...
670!> \param virtual ...
671!> \param homo_irred ...
672!> \param unit_nr ...
673!> \param flag_TDA ...
674!> \param fm_eigvec_Y ...
675! **************************************************************************************************
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)
679
680 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: exc_ens
681 TYPE(cp_fm_type), INTENT(IN) :: fm_eigvec_x
682 TYPE(mp2_type), INTENT(INOUT) :: mp2_env
683 TYPE(qs_environment_type), POINTER :: qs_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
688
689 CHARACTER(LEN=*), PARAMETER :: routinen = 'postprocess_bse'
690
691 CHARACTER(LEN=10) :: info_approximation, multiplet
692 INTEGER :: handle, i_exc, idir, n_moments_di, &
693 n_moments_quad
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
697 TYPE(cp_fm_type) :: fm_x_ia, fm_y_ia
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
700 TYPE(exciton_descr_type), ALLOCATABLE, &
701 DIMENSION(:) :: exc_descr
702
703 CALL timeset(routinen, handle)
704
705 !Prepare variables for printing
706 IF (mp2_env%bse%bse_spin_config == 0) THEN
707 multiplet = "Singlet"
708 alpha = 2.0_dp
709 ELSE
710 multiplet = "Triplet"
711 alpha = 0.0_dp
712 END IF
713 IF (.NOT. PRESENT(flag_tda)) THEN
714 flag_tda = .false.
715 END IF
716 IF (flag_tda) THEN
717 info_approximation = " -TDA- "
718 ELSE
719 info_approximation = "-ABBA-"
720 END IF
721
722 n_moments_di = 3
723 n_moments_quad = 9
724 ! Compute BSE dipoles and oscillator strengths - Keep in memory for later usage
725 ! Need dipoles also for spatial expectation values, which are well-defined also for triplets
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))
730 ! Obtain dipoles in MO basis
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)
734 ! Compute exciton descriptors from these multipoles
735 IF (mp2_env%bse%num_print_exc_descr > 0) THEN
736 ! Obtain quadrupoles in MO basis
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)
743 ! Iterate over excitation index outside of routine to make it compatible with tddft module
744 ALLOCATE (exc_descr(mp2_env%bse%num_print_exc_descr))
745 DO i_exc = 1, mp2_env%bse%num_print_exc_descr
746 CALL reshuffle_eigvec(fm_eigvec_x, fm_x_ia, homo, virtual, i_exc, &
747 .false., unit_nr, mp2_env)
748 IF (.NOT. flag_tda) THEN
749 CALL reshuffle_eigvec(fm_eigvec_y, fm_y_ia, homo, virtual, i_exc, &
750 .false., unit_nr, mp2_env)
751
752 CALL get_exciton_descriptors(exc_descr, fm_x_ia, &
753 fm_quadpole_ij_trunc, fm_quadpole_ab_trunc, &
754 fm_quadpole_ai_trunc, &
755 i_exc, homo, virtual, &
756 fm_y_ia)
757 ELSE
758 CALL get_exciton_descriptors(exc_descr, fm_x_ia, &
759 fm_quadpole_ij_trunc, fm_quadpole_ab_trunc, &
760 fm_quadpole_ai_trunc, &
761 i_exc, homo, virtual)
762 END IF
763 CALL cp_fm_release(fm_x_ia)
764 IF (.NOT. flag_tda) THEN
765 CALL cp_fm_release(fm_y_ia)
766 END IF
767 END DO
768 END IF
769
770 IF (mp2_env%bse%bse_spin_config == 0) THEN
771 CALL get_oscillator_strengths(fm_eigvec_x, exc_ens, fm_dipole_ai_trunc, &
772 trans_mom_bse, oscill_str, polarizability_residues, &
773 mp2_env, homo, virtual, unit_nr, &
774 fm_eigvec_y)
775 END IF
776
777 ! Prints basic definitions used in BSE calculation
778 CALL print_output_header(homo, virtual, homo_irred, flag_tda, &
779 multiplet, alpha, mp2_env, unit_nr)
780
781 ! Prints excitation energies up to user-specified number
782 CALL print_excitation_energies(exc_ens, homo, virtual, flag_tda, multiplet, &
783 info_approximation, mp2_env, unit_nr)
784
785 ! Print single particle transition amplitudes, i.e. components of eigenvectors X and Y
786 CALL print_transition_amplitudes(fm_eigvec_x, homo, virtual, homo_irred, &
787 info_approximation, mp2_env, unit_nr, fm_eigvec_y)
788
789 ! Prints optical properties, if state is a singlet
790 CALL print_optical_properties(exc_ens, oscill_str, trans_mom_bse, polarizability_residues, &
791 homo, virtual, homo_irred, flag_tda, &
792 info_approximation, mp2_env, unit_nr)
793 ! Print exciton descriptors if keyword is invoked
794 IF (mp2_env%bse%num_print_exc_descr > 0) THEN
795 CALL print_exciton_descriptors(exc_descr, ref_point_multipole, unit_nr, &
796 mp2_env%bse%num_print_exc_descr, mp2_env%bse%bse_debug_print, &
797 mp2_env%bse%print_directional_exc_descr, &
798 'BSE|', qs_env)
799 END IF
800
801 ! Compute and print excitation wavefunctions
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|'
808 END IF
809 CALL calculate_ntos(fm_eigvec_x, fm_eigvec_y, &
810 mo_coeff, homo, virtual, &
811 info_approximation, &
812 oscill_str, &
813 qs_env, unit_nr, mp2_env)
814 END IF
815
816 DO idir = 1, n_moments_di
817 CALL cp_fm_release(fm_dipole_ai_trunc(idir))
818 CALL cp_fm_release(fm_dipole_ij_trunc(idir))
819 CALL cp_fm_release(fm_dipole_ab_trunc(idir))
820 END DO
821 IF (mp2_env%bse%num_print_exc_descr > 0) THEN
822 DO idir = 1, n_moments_quad
823 CALL cp_fm_release(fm_quadpole_ai_trunc(idir))
824 CALL cp_fm_release(fm_quadpole_ij_trunc(idir))
825 CALL cp_fm_release(fm_quadpole_ab_trunc(idir))
826 END DO
827 DEALLOCATE (fm_quadpole_ai_trunc, fm_quadpole_ij_trunc, fm_quadpole_ab_trunc)
828 DEALLOCATE (exc_descr)
829 END IF
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)
834 END IF
835 IF (unit_nr > 0) THEN
836 WRITE (unit_nr, '(T2,A4)') 'BSE|'
837 WRITE (unit_nr, '(T2,A4)') 'BSE|'
838 END IF
839
840 CALL timestop(handle)
841
842 END SUBROUTINE postprocess_bse
843
844END MODULE bse_full_diag
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.
Definition bse_print.F:13
subroutine, public print_output_header(homo, virtual, homo_irred, flag_tda, multiplet, alpha, mp2_env, unit_nr)
...
Definition bse_print.F:98
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)
...
Definition bse_print.F:339
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)
Definition bse_print.F:507
subroutine, public print_transition_amplitudes(fm_eigvec_x, homo, virtual, homo_irred, info_approximation, mp2_env, unit_nr, fm_eigvec_y)
...
Definition bse_print.F:264
subroutine, public print_excitation_energies(exc_ens, homo, virtual, flag_tda, multiplet, info_approximation, mp2_env, unit_nr)
...
Definition bse_print.F:215
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.
Definition bse_util.F:13
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)
...
Definition bse_util.F:1715
subroutine, public reshuffle_eigvec(fm_eigvec, fm_eigvec_reshuffled, homo, virtual, n_exc, do_transpose, unit_nr, mp2_env)
...
Definition bse_util.F:1351
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...
Definition bse_util.F:787
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,...
Definition bse_util.F:199
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...
Definition cp_fm_diag.F:17
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...
Definition cp_fm_diag.F:239
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_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...
Definition gamma.F:15
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public bse_singlet
integer, parameter, public bse_triplet
integer, parameter, public bse_screening_alpha
integer, parameter, public bse_screening_rpa
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
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
represent a full matrix
Contains information on the excited states energy.
stores all the informations relevant to an mpi environment