(git:bb35279)
Loading...
Searching...
No Matches
rpa_gw_kpoints_util.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 Routines treating GW and RPA calculations with kpoints
10!> \par History
11!> since 2018 continuous development [J. Wilhelm]
12! **************************************************************************************************
14 USE cell_types, ONLY: cell_type,&
15 get_cell,&
16 pbc
23 USE cp_cfm_diag, ONLY: cp_cfm_geeig,&
26 USE cp_cfm_types, ONLY: cp_cfm_create,&
34 USE cp_dbcsr_api, ONLY: &
38 dbcsr_release, dbcsr_set, dbcsr_transposed, dbcsr_type, dbcsr_type_no_symmetry
50 USE hfx_types, ONLY: hfx_release
51 USE input_constants, ONLY: cholesky_off,&
55 USE kinds, ONLY: dp
59 USE kpoint_types, ONLY: get_kpoint_info,&
62 USE machine, ONLY: m_walltime
63 USE mathconstants, ONLY: gaussi,&
64 twopi,&
65 z_one,&
66 z_zero
67 USE mathlib, ONLY: invmat
74 USE qs_mo_types, ONLY: get_mo_set
80#include "./base/base_uses.f90"
81
82 IMPLICIT NONE
83
84 PRIVATE
85
86 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rpa_gw_kpoints_util'
87
91
92CONTAINS
93
94! **************************************************************************************************
95!> \brief ...
96!> \param dimen_RI ...
97!> \param num_integ_points ...
98!> \param jquad ...
99!> \param nkp ...
100!> \param count_ev_sc_GW ...
101!> \param para_env ...
102!> \param Erpa ...
103!> \param tau_tj ...
104!> \param tj ...
105!> \param wj ...
106!> \param weights_cos_tf_w_to_t ...
107!> \param wkp_W ...
108!> \param do_gw_im_time ...
109!> \param do_ri_Sigma_x ...
110!> \param do_kpoints_from_Gamma ...
111!> \param cfm_mat_Q ...
112!> \param ikp_local ...
113!> \param mat_P_omega ...
114!> \param mat_P_omega_kp ...
115!> \param qs_env ...
116!> \param eps_filter_im_time ...
117!> \param unit_nr ...
118!> \param kpoints ...
119!> \param fm_mat_Minv_L_kpoints ...
120!> \param fm_matrix_L_kpoints ...
121!> \param fm_mat_W ...
122!> \param fm_mat_RI_global_work ...
123!> \param mat_MinvVMinv ...
124!> \param fm_matrix_Minv ...
125!> \param fm_matrix_Minv_Vtrunc_Minv ...
126! **************************************************************************************************
127 SUBROUTINE invert_eps_compute_w_and_erpa_kp(dimen_RI, num_integ_points, jquad, nkp, count_ev_sc_GW, para_env, &
128 Erpa, tau_tj, tj, wj, weights_cos_tf_w_to_t, wkp_W, do_gw_im_time, &
129 do_ri_Sigma_x, do_kpoints_from_Gamma, &
130 cfm_mat_Q, ikp_local, mat_P_omega, mat_P_omega_kp, &
131 qs_env, eps_filter_im_time, unit_nr, kpoints, fm_mat_Minv_L_kpoints, &
132 fm_matrix_L_kpoints, fm_mat_W, &
133 fm_mat_RI_global_work, mat_MinvVMinv, fm_matrix_Minv, &
134 fm_matrix_Minv_Vtrunc_Minv)
135
136 INTEGER, INTENT(IN) :: dimen_ri, num_integ_points, jquad, nkp, &
137 count_ev_sc_gw
138 TYPE(mp_para_env_type), POINTER :: para_env
139 REAL(kind=dp), INTENT(INOUT) :: erpa
140 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: tau_tj, tj, wj
141 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
142 INTENT(IN) :: weights_cos_tf_w_to_t
143 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: wkp_w
144 LOGICAL, INTENT(IN) :: do_gw_im_time, do_ri_sigma_x, &
145 do_kpoints_from_gamma
146 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
147 INTEGER, INTENT(IN) :: ikp_local
148 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
149 TYPE(qs_environment_type), POINTER :: qs_env
150 REAL(kind=dp), INTENT(IN) :: eps_filter_im_time
151 INTEGER, INTENT(IN) :: unit_nr
152 TYPE(kpoint_type), POINTER :: kpoints
153 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_mat_minv_l_kpoints, &
154 fm_matrix_l_kpoints
155 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_w
156 TYPE(cp_fm_type) :: fm_mat_ri_global_work
157 TYPE(dbcsr_p_type), INTENT(IN) :: mat_minvvminv
158 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_minv, &
159 fm_matrix_minv_vtrunc_minv
160
161 CHARACTER(LEN=*), PARAMETER :: routinen = 'invert_eps_compute_W_and_Erpa_kp'
162
163 INTEGER :: handle, ikp
164 LOGICAL :: do_this_ikp
165 REAL(kind=dp) :: t1, t2
166 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: trace_qomega
167
168 CALL timeset(routinen, handle)
169
170 t1 = m_walltime()
171
172 IF (do_kpoints_from_gamma) THEN
173 CALL get_mat_cell_t_from_mat_gamma(mat_p_omega(jquad, :), qs_env, kpoints, jquad, unit_nr)
174 END IF
175
176 CALL transform_p_from_real_space_to_kpoints(mat_p_omega, mat_p_omega_kp, &
177 kpoints, eps_filter_im_time, jquad)
178
179 ALLOCATE (trace_qomega(dimen_ri))
180
181 IF (unit_nr > 0) WRITE (unit_nr, '(/T3,A,1X,I3)') &
182 'GW_INFO| Computing chi and W frequency point', jquad
183
184 DO ikp = 1, nkp
185
186 ! parallization, we either have all kpoints on all processors or a single kpoint per group
187 do_this_ikp = (ikp_local == -1) .OR. (ikp_local == 0 .AND. ikp == 1) .OR. (ikp_local == ikp)
188 IF (.NOT. do_this_ikp) cycle
189
190 ! 1. remove all spurious negative eigenvalues from P(iw,k), multiplication Q(iw,k) = K^H(k)P(iw,k)K(k)
191 CALL compute_q_kp_rpa(cfm_mat_q, &
192 mat_p_omega_kp, &
193 fm_mat_minv_l_kpoints(ikp, 1), &
194 fm_mat_minv_l_kpoints(ikp, 2), &
195 fm_mat_ri_global_work, &
196 dimen_ri, ikp, nkp, ikp_local, para_env, &
197 qs_env%mp2_env%ri_rpa_im_time%make_chi_pos_definite)
198
199 ! 2. Cholesky decomposition of Id + Q(iw,k)
200 CALL cholesky_decomp_q(cfm_mat_q, para_env, trace_qomega, dimen_ri)
201
202 ! 3. Computing E_c^RPA = E_c^RPA + a_w/N_k*sum_k ln[det(1+Q(iw,k))-Tr(Q(iw,k))]
203 CALL frequency_and_kpoint_integration(erpa, cfm_mat_q, para_env, trace_qomega, &
204 dimen_ri, wj(jquad), kpoints%wkp(ikp))
205
206 IF (do_gw_im_time) THEN
207
208 ! compute S^-1*V*S^-1 for exchange part of the self-energy in real space as W in real space
209 IF (do_ri_sigma_x .AND. jquad == 1 .AND. count_ev_sc_gw == 1 .AND. do_kpoints_from_gamma) THEN
210
211 CALL dbcsr_set(mat_minvvminv%matrix, 0.0_dp)
212 CALL copy_fm_to_dbcsr(fm_matrix_minv_vtrunc_minv(1, 1), mat_minvvminv%matrix, keep_sparsity=.false.)
213
214 END IF
215 IF (do_kpoints_from_gamma) THEN
216 CALL compute_wc_real_space_tau_gw(fm_mat_w, cfm_mat_q, &
217 fm_matrix_l_kpoints(ikp, 1), &
218 fm_matrix_l_kpoints(ikp, 2), &
219 dimen_ri, num_integ_points, jquad, &
220 ikp, tj, tau_tj, weights_cos_tf_w_to_t, &
221 ikp_local, para_env, kpoints, qs_env, wkp_w)
222 END IF
223
224 END IF
225 END DO
226
227 ! after the transform of (eps(iw)-1)^-1 from iw to it is done, multiply with V^1/2 to obtain W(it)
228 IF (do_gw_im_time .AND. do_kpoints_from_gamma .AND. jquad == num_integ_points) THEN
229 CALL wc_to_minv_wc_minv(fm_mat_w, fm_matrix_minv, para_env, dimen_ri, num_integ_points)
230 CALL deallocate_kp_matrices(fm_matrix_l_kpoints, fm_mat_minv_l_kpoints)
231 END IF
232
233 DEALLOCATE (trace_qomega)
234
235 t2 = m_walltime()
236
237 IF (unit_nr > 0) WRITE (unit_nr, '(T6,A,T56,F25.1)') 'Execution time (s):', t2 - t1
238
239 CALL timestop(handle)
240
242
243! **************************************************************************************************
244!> \brief ...
245!> \param fm_matrix_L_kpoints ...
246!> \param fm_mat_Minv_L_kpoints ...
247! **************************************************************************************************
248 SUBROUTINE deallocate_kp_matrices(fm_matrix_L_kpoints, fm_mat_Minv_L_kpoints)
249
250 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_l_kpoints, &
251 fm_mat_minv_l_kpoints
252
253 CHARACTER(LEN=*), PARAMETER :: routinen = 'deallocate_kp_matrices'
254
255 INTEGER :: handle
256
257 CALL timeset(routinen, handle)
258
259 CALL cp_fm_release(fm_mat_minv_l_kpoints)
260 CALL cp_fm_release(fm_matrix_l_kpoints)
261
262 CALL timestop(handle)
263
264 END SUBROUTINE deallocate_kp_matrices
265
266! **************************************************************************************************
267!> \brief ...
268!> \param matrix ...
269!> \param threshold ...
270!> \param exponent ...
271!> \param min_eigval ...
272! **************************************************************************************************
273 SUBROUTINE cp_cfm_power(matrix, threshold, exponent, min_eigval)
274 TYPE(cp_cfm_type), INTENT(INOUT) :: matrix
275 REAL(kind=dp) :: threshold, exponent
276 REAL(kind=dp), OPTIONAL :: min_eigval
277
278 CHARACTER(LEN=*), PARAMETER :: routinen = 'cp_cfm_power'
279
280 COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues_exponent
281 INTEGER :: handle, i, ncol_global, nrow_global
282 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues
283 TYPE(cp_cfm_type) :: cfm_work
284
285 CALL timeset(routinen, handle)
286
287 CALL cp_cfm_create(cfm_work, matrix%matrix_struct)
288 CALL cp_cfm_set_all(cfm_work, z_zero)
289
290 ! Test that matrix is square
291 CALL cp_cfm_get_info(matrix, nrow_global=nrow_global, ncol_global=ncol_global)
292 cpassert(nrow_global == ncol_global)
293 ALLOCATE (eigenvalues(nrow_global))
294 eigenvalues(:) = 0.0_dp
295 ALLOCATE (eigenvalues_exponent(nrow_global))
296 eigenvalues_exponent(:) = z_zero
297
298 ! Diagonalize matrix: get eigenvectors and eigenvalues
299 CALL cp_cfm_heevd(matrix, cfm_work, eigenvalues)
300
301 DO i = 1, nrow_global
302 IF (eigenvalues(i) > threshold) THEN
303 eigenvalues_exponent(i) = cmplx((eigenvalues(i))**(0.5_dp*exponent), threshold, kind=dp)
304 ELSE
305 IF (PRESENT(min_eigval)) THEN
306 eigenvalues_exponent(i) = cmplx(min_eigval, 0.0_dp, kind=dp)
307 ELSE
308 eigenvalues_exponent(i) = z_zero
309 END IF
310 END IF
311 END DO
312
313 CALL cp_cfm_column_scale(cfm_work, eigenvalues_exponent)
314
315 CALL parallel_gemm("N", "C", nrow_global, nrow_global, nrow_global, z_one, &
316 cfm_work, cfm_work, z_zero, matrix)
317
318 DEALLOCATE (eigenvalues, eigenvalues_exponent)
319
320 CALL cp_cfm_release(cfm_work)
321
322 CALL timestop(handle)
323
324 END SUBROUTINE cp_cfm_power
325
326! **************************************************************************************************
327!> \brief ...
328!> \param cfm_mat_Q ...
329!> \param mat_P_omega_kp ...
330!> \param fm_mat_L_re ...
331!> \param fm_mat_L_im ...
332!> \param fm_mat_RI_global_work ...
333!> \param dimen_RI ...
334!> \param ikp ...
335!> \param nkp ...
336!> \param ikp_local ...
337!> \param para_env ...
338!> \param make_chi_pos_definite ...
339! **************************************************************************************************
340 SUBROUTINE compute_q_kp_rpa(cfm_mat_Q, mat_P_omega_kp, fm_mat_L_re, fm_mat_L_im, &
341 fm_mat_RI_global_work, dimen_RI, ikp, nkp, ikp_local, para_env, &
342 make_chi_pos_definite)
343
344 TYPE(cp_cfm_type) :: cfm_mat_q
345 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega_kp
346 TYPE(cp_fm_type) :: fm_mat_l_re, fm_mat_l_im, &
347 fm_mat_ri_global_work
348 INTEGER, INTENT(IN) :: dimen_ri, ikp, nkp, ikp_local
349 TYPE(mp_para_env_type), POINTER :: para_env
350 LOGICAL, INTENT(IN) :: make_chi_pos_definite
351
352 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_Q_kp_RPA'
353
354 INTEGER :: handle
355 TYPE(cp_cfm_type) :: cfm_mat_l, cfm_mat_work
356 TYPE(cp_fm_type) :: fm_mat_work
357
358 CALL timeset(routinen, handle)
359
360 CALL cp_cfm_create(cfm_mat_work, fm_mat_l_re%matrix_struct)
361 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
362
363 CALL cp_cfm_create(cfm_mat_l, fm_mat_l_re%matrix_struct)
364 CALL cp_cfm_set_all(cfm_mat_l, z_zero)
365
366 CALL cp_fm_create(fm_mat_work, fm_mat_l_re%matrix_struct)
367 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
368
369 ! 1. Convert the dbcsr matrix mat_P_omega_kp (that is chi(k,iw)) to a full matrix and
370 ! distribute it to subgroups
371 CALL mat_p_to_subgroup(mat_p_omega_kp, fm_mat_ri_global_work, &
372 fm_mat_work, cfm_mat_q, ikp, nkp, ikp_local, para_env)
373
374 ! 2. Remove all negative eigenvalues from chi(k,iw)
375 IF (make_chi_pos_definite) THEN
376 CALL cp_cfm_power(cfm_mat_q, threshold=0.0_dp, exponent=1.0_dp)
377 END IF
378
379 ! 3. Copy fm_mat_L_re and fm_mat_L_re to cfm_mat_L
380 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_l, z_one, fm_mat_l_re)
381 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_l, gaussi, fm_mat_l_im)
382
383 ! 4. work = P(iw,k)*L(k)
384 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_q, cfm_mat_l, &
385 z_zero, cfm_mat_work)
386
387 ! 5. Q(iw,k) = L^H(k)*work
388 CALL parallel_gemm('C', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_l, cfm_mat_work, &
389 z_zero, cfm_mat_q)
390
391 CALL cp_cfm_release(cfm_mat_work)
392 CALL cp_cfm_release(cfm_mat_l)
393 CALL cp_fm_release(fm_mat_work)
394
395 CALL timestop(handle)
396
397 END SUBROUTINE compute_q_kp_rpa
398
399! **************************************************************************************************
400!> \brief ...
401!> \param mat_P_omega_kp ...
402!> \param fm_mat_RI_global_work ...
403!> \param fm_mat_work ...
404!> \param cfm_mat_Q ...
405!> \param ikp ...
406!> \param nkp ...
407!> \param ikp_local ...
408!> \param para_env ...
409! **************************************************************************************************
410 SUBROUTINE mat_p_to_subgroup(mat_P_omega_kp, fm_mat_RI_global_work, &
411 fm_mat_work, cfm_mat_Q, ikp, nkp, ikp_local, para_env)
412
413 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega_kp
414 TYPE(cp_fm_type), INTENT(INOUT) :: fm_mat_ri_global_work, fm_mat_work
415 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
416 INTEGER, INTENT(IN) :: ikp, nkp, ikp_local
417 TYPE(mp_para_env_type), POINTER :: para_env
418
419 CHARACTER(LEN=*), PARAMETER :: routinen = 'mat_P_to_subgroup'
420
421 INTEGER :: handle, jkp
422 TYPE(cp_fm_type) :: fm_dummy
423 TYPE(dbcsr_type), POINTER :: mat_p_omega_im, mat_p_omega_re
424
425 CALL timeset(routinen, handle)
426
427 IF (ikp_local == -1) THEN
428
429 mat_p_omega_re => mat_p_omega_kp(1, ikp)%matrix
430 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
431 CALL copy_dbcsr_to_fm(mat_p_omega_re, fm_mat_work)
432 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_q, z_one, fm_mat_work)
433
434 mat_p_omega_im => mat_p_omega_kp(2, ikp)%matrix
435 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
436 CALL copy_dbcsr_to_fm(mat_p_omega_im, fm_mat_work)
437 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_q, gaussi, fm_mat_work)
438
439 ELSE
440
441 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
442
443 DO jkp = 1, nkp
444
445 mat_p_omega_re => mat_p_omega_kp(1, jkp)%matrix
446
447 CALL cp_fm_set_all(fm_mat_ri_global_work, 0.0_dp)
448 CALL copy_dbcsr_to_fm(mat_p_omega_re, fm_mat_ri_global_work)
449
450 CALL para_env%sync()
451
452 IF (ikp_local == jkp) THEN
453 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_mat_work, para_env)
454 ELSE
455 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_dummy, para_env)
456 END IF
457
458 CALL para_env%sync()
459
460 END DO
461
462 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_q, z_one, fm_mat_work)
463
464 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
465
466 DO jkp = 1, nkp
467
468 mat_p_omega_im => mat_p_omega_kp(2, jkp)%matrix
469
470 CALL cp_fm_set_all(fm_mat_ri_global_work, 0.0_dp)
471 CALL copy_dbcsr_to_fm(mat_p_omega_im, fm_mat_ri_global_work)
472
473 CALL para_env%sync()
474
475 IF (ikp_local == jkp) THEN
476 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_mat_work, para_env)
477 ELSE
478 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_dummy, para_env)
479 END IF
480
481 CALL para_env%sync()
482
483 END DO
484
485 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_q, gaussi, fm_mat_work)
486
487 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
488
489 END IF
490
491 CALL para_env%sync()
492
493 CALL timestop(handle)
494
495 END SUBROUTINE mat_p_to_subgroup
496
497! **************************************************************************************************
498!> \brief ...
499!> \param cfm_mat_Q ...
500!> \param para_env ...
501!> \param trace_Qomega ...
502!> \param dimen_RI ...
503! **************************************************************************************************
504 SUBROUTINE cholesky_decomp_q(cfm_mat_Q, para_env, trace_Qomega, dimen_RI)
505
506 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
507 TYPE(mp_para_env_type), INTENT(IN) :: para_env
508 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: trace_qomega
509 INTEGER, INTENT(IN) :: dimen_ri
510
511 CHARACTER(LEN=*), PARAMETER :: routinen = 'cholesky_decomp_Q'
512
513 INTEGER :: handle, i_global, iib, info_chol, &
514 j_global, jjb, ncol_local, nrow_local
515 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
516 TYPE(cp_cfm_type) :: cfm_mat_q_tmp, cfm_mat_work
517
518 CALL timeset(routinen, handle)
519
520 CALL cp_cfm_create(cfm_mat_work, cfm_mat_q%matrix_struct)
521 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
522
523 CALL cp_cfm_create(cfm_mat_q_tmp, cfm_mat_q%matrix_struct)
524 CALL cp_cfm_set_all(cfm_mat_q_tmp, z_zero)
525
526 ! get info of fm_mat_Q
527 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
528 nrow_local=nrow_local, &
529 ncol_local=ncol_local, &
530 row_indices=row_indices, &
531 col_indices=col_indices)
532
533 ! calculate the trace of Q and add 1 on the diagonal
534 trace_qomega = 0.0_dp
535!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
536!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,trace_Qomega,cfm_mat_Q,dimen_RI)
537 DO jjb = 1, ncol_local
538 j_global = col_indices(jjb)
539 DO iib = 1, nrow_local
540 i_global = row_indices(iib)
541 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
542 trace_qomega(i_global) = real(cfm_mat_q%local_data(iib, jjb))
543 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) + z_one
544 END IF
545 END DO
546 END DO
547 CALL para_env%sum(trace_qomega)
548
549 CALL cp_cfm_to_cfm(cfm_mat_q, cfm_mat_q_tmp)
550
551 CALL cp_cfm_cholesky_decompose(matrix=cfm_mat_q, n=dimen_ri, info_out=info_chol)
552
553 cpassert(info_chol == 0)
554
555 CALL cp_cfm_release(cfm_mat_work)
556 CALL cp_cfm_release(cfm_mat_q_tmp)
557
558 CALL timestop(handle)
559
560 END SUBROUTINE cholesky_decomp_q
561
562! **************************************************************************************************
563!> \brief ...
564!> \param Erpa ...
565!> \param cfm_mat_Q ...
566!> \param para_env ...
567!> \param trace_Qomega ...
568!> \param dimen_RI ...
569!> \param freq_weight ...
570!> \param kp_weight ...
571! **************************************************************************************************
572 SUBROUTINE frequency_and_kpoint_integration(Erpa, cfm_mat_Q, para_env, trace_Qomega, &
573 dimen_RI, freq_weight, kp_weight)
574
575 REAL(kind=dp), INTENT(INOUT) :: erpa
576 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
577 TYPE(mp_para_env_type), INTENT(IN) :: para_env
578 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: trace_qomega
579 INTEGER, INTENT(IN) :: dimen_ri
580 REAL(kind=dp), INTENT(IN) :: freq_weight, kp_weight
581
582 CHARACTER(LEN=*), PARAMETER :: routinen = 'frequency_and_kpoint_integration'
583
584 INTEGER :: handle, i_global, iib, j_global, jjb, &
585 ncol_local, nrow_local
586 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
587 REAL(kind=dp) :: fcomega
588 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: q_log
589
590 CALL timeset(routinen, handle)
591
592 ! get info of cholesky_decomposed(fm_mat_Q)
593 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
594 nrow_local=nrow_local, &
595 ncol_local=ncol_local, &
596 row_indices=row_indices, &
597 col_indices=col_indices)
598
599 ALLOCATE (q_log(dimen_ri))
600 q_log = 0.0_dp
601!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
602!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,Q_log,cfm_mat_Q,dimen_RI)
603 DO jjb = 1, ncol_local
604 j_global = col_indices(jjb)
605 DO iib = 1, nrow_local
606 i_global = row_indices(iib)
607 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
608 q_log(i_global) = 2.0_dp*log(real(cfm_mat_q%local_data(iib, jjb)))
609 END IF
610 END DO
611 END DO
612 CALL para_env%sum(q_log)
613
614 fcomega = 0.0_dp
615 DO iib = 1, dimen_ri
616 IF (modulo(iib, para_env%num_pe) /= para_env%mepos) cycle
617 ! FComega=FComega+(LOG(Q_log(iiB))-trace_Qomega(iiB))/2.0_dp
618 fcomega = fcomega + (q_log(iib) - trace_qomega(iib))/2.0_dp
619 END DO
620
621 erpa = erpa + fcomega*freq_weight*kp_weight
622
623 DEALLOCATE (q_log)
624
625 CALL timestop(handle)
626
627 END SUBROUTINE frequency_and_kpoint_integration
628
629! **************************************************************************************************
630!> \brief ...
631!> \param tj_dummy ...
632!> \param tau_tj_dummy ...
633!> \param weights_cos_tf_w_to_t_dummy ...
634! **************************************************************************************************
635 SUBROUTINE get_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
636
637 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
638 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
639 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
640 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
641
642 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_dummys'
643
644 INTEGER :: handle
645
646 CALL timeset(routinen, handle)
647
648 ALLOCATE (weights_cos_tf_w_to_t_dummy(1, 1))
649 ALLOCATE (tj_dummy(1))
650 ALLOCATE (tau_tj_dummy(1))
651
652 tj_dummy(1) = 0.0_dp
653 tau_tj_dummy(1) = 0.0_dp
654 weights_cos_tf_w_to_t_dummy(1, 1) = 1.0_dp
655
656 CALL timestop(handle)
657
658 END SUBROUTINE
659
660! **************************************************************************************************
661!> \brief ...
662!> \param tj_dummy ...
663!> \param tau_tj_dummy ...
664!> \param weights_cos_tf_w_to_t_dummy ...
665! **************************************************************************************************
666 SUBROUTINE release_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
667
668 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
669 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
670 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
671 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
672
673 CHARACTER(LEN=*), PARAMETER :: routinen = 'release_dummys'
674
675 INTEGER :: handle
676
677 CALL timeset(routinen, handle)
678
679 DEALLOCATE (weights_cos_tf_w_to_t_dummy)
680 DEALLOCATE (tj_dummy)
681 DEALLOCATE (tau_tj_dummy)
682
683 CALL timestop(handle)
684
685 END SUBROUTINE
686
687! **************************************************************************************************
688!> \brief ...
689!> \param mat_P_omega ...
690!> \param qs_env ...
691!> \param kpoints ...
692!> \param jquad ...
693!> \param unit_nr ...
694! **************************************************************************************************
695 SUBROUTINE get_mat_cell_t_from_mat_gamma(mat_P_omega, qs_env, kpoints, jquad, unit_nr)
696 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN) :: mat_p_omega
697 TYPE(qs_environment_type), POINTER :: qs_env
698 TYPE(kpoint_type), POINTER :: kpoints
699 INTEGER, INTENT(IN) :: jquad, unit_nr
700
701 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_mat_cell_T_from_mat_gamma'
702
703 INTEGER :: col, handle, i_cell, i_dim, j_cell, &
704 num_cells_p, num_integ_points, row
705 INTEGER, DIMENSION(3) :: cell_grid_p, periodic
706 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell_p
707 LOGICAL :: i_cell_is_the_minimum_image_cell
708 REAL(kind=dp) :: abs_rab_cell_i, abs_rab_cell_j
709 REAL(kind=dp), DIMENSION(3) :: cell_vector, cell_vector_j, rab_cell_i, &
710 rab_cell_j
711 REAL(kind=dp), DIMENSION(3, 3) :: hmat
712 REAL(kind=dp), DIMENSION(:, :), POINTER :: data_block
713 TYPE(cell_type), POINTER :: cell
714 TYPE(dbcsr_iterator_type) :: iter
715 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
716
717 CALL timeset(routinen, handle)
718
719 NULLIFY (cell, particle_set)
720 CALL get_qs_env(qs_env, cell=cell, &
721 particle_set=particle_set)
722 CALL get_cell(cell=cell, h=hmat, periodic=periodic)
723
724 DO i_dim = 1, 3
725 ! we have at most 3 neigboring cells per dimension and at least one because
726 ! the density response at Gamma is only divided to neighboring
727 IF (periodic(i_dim) == 1) THEN
728 cell_grid_p(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
729 ELSE
730 cell_grid_p(i_dim) = 1
731 END IF
732 END DO
733
734 ! overwrite the cell indices in kpoints
735 CALL init_cell_index_rpa(cell_grid_p, kpoints%cell_to_index, kpoints%index_to_cell, cell)
736
737 index_to_cell_p => kpoints%index_to_cell
738
739 num_cells_p = SIZE(index_to_cell_p, 2)
740
741 num_integ_points = SIZE(mat_p_omega, 1)
742
743 ! first, copy the Gamma-only result from mat_P_omega(1) into all other matrices and
744 ! remove the blocks later which do not belong to the cell index
745 DO i_cell = 2, num_cells_p
746 CALL dbcsr_copy(mat_p_omega(i_cell)%matrix, &
747 mat_p_omega(1)%matrix)
748 END DO
749
750 IF (jquad == 1 .AND. unit_nr > 0) THEN
751 WRITE (unit_nr, '(T3,A,T66,ES15.2)') 'GW_INFO| RI regularization parameter: ', &
752 qs_env%mp2_env%ri_rpa_im_time%regularization_RI
753 WRITE (unit_nr, '(T3,A,T66,ES15.2)') 'GW_INFO| eps_eigval_S: ', &
754 qs_env%mp2_env%ri_rpa_im_time%eps_eigval_S
755 IF (qs_env%mp2_env%ri_rpa_im_time%make_chi_pos_definite) THEN
756 WRITE (unit_nr, '(T3,A,T81)') &
757 'GW_INFO| Make chi(iw,k) positive definite? TRUE'
758 ELSE
759 WRITE (unit_nr, '(T3,A,T81)') &
760 'GW_INFO| Make chi(iw,k) positive definite? FALSE'
761 END IF
762
763 END IF
764
765 DO i_cell = 1, num_cells_p
766
767 CALL dbcsr_iterator_start(iter, mat_p_omega(i_cell)%matrix)
768 DO WHILE (dbcsr_iterator_blocks_left(iter))
769 CALL dbcsr_iterator_next_block(iter, row, col, data_block)
770
771 cell_vector(1:3) = matmul(hmat, real(index_to_cell_p(1:3, i_cell), dp))
772 rab_cell_i(1:3) = pbc(particle_set(row)%r(1:3), cell) - &
773 (pbc(particle_set(col)%r(1:3), cell) + cell_vector(1:3))
774 abs_rab_cell_i = sqrt(rab_cell_i(1)**2 + rab_cell_i(2)**2 + rab_cell_i(3)**2)
775
776 ! minimum image convention
777 i_cell_is_the_minimum_image_cell = .true.
778 DO j_cell = 1, num_cells_p
779 cell_vector_j(1:3) = matmul(hmat, real(index_to_cell_p(1:3, j_cell), dp))
780 rab_cell_j(1:3) = pbc(particle_set(row)%r(1:3), cell) - &
781 (pbc(particle_set(col)%r(1:3), cell) + cell_vector_j(1:3))
782 abs_rab_cell_j = sqrt(rab_cell_j(1)**2 + rab_cell_j(2)**2 + rab_cell_j(3)**2)
783
784 IF (abs_rab_cell_i > abs_rab_cell_j + 1.0e-6_dp) THEN
785 i_cell_is_the_minimum_image_cell = .false.
786 END IF
787 END DO
788
789 IF (.NOT. i_cell_is_the_minimum_image_cell) THEN
790 data_block(:, :) = data_block(:, :)*0.0_dp
791 END IF
792
793 END DO
794 CALL dbcsr_iterator_stop(iter)
795
796 END DO
797
798 CALL timestop(handle)
799
800 END SUBROUTINE get_mat_cell_t_from_mat_gamma
801
802! **************************************************************************************************
803!> \brief ...
804!> \param mat_P_omega ...
805!> \param mat_P_omega_kp ...
806!> \param kpoints ...
807!> \param eps_filter_im_time ...
808!> \param jquad ...
809! **************************************************************************************************
810 SUBROUTINE transform_p_from_real_space_to_kpoints(mat_P_omega, mat_P_omega_kp, &
811 kpoints, eps_filter_im_time, jquad)
812
813 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
814 TYPE(kpoint_type), POINTER :: kpoints
815 REAL(kind=dp), INTENT(IN) :: eps_filter_im_time
816 INTEGER, INTENT(IN) :: jquad
817
818 CHARACTER(LEN=*), PARAMETER :: routinen = 'transform_P_from_real_space_to_kpoints'
819
820 INTEGER :: handle, icell, nkp, num_integ_points
821
822 CALL timeset(routinen, handle)
823
824 num_integ_points = SIZE(mat_p_omega, 1)
825 nkp = SIZE(mat_p_omega, 2)
826
827 CALL real_space_to_kpoint_transform_rpa(mat_p_omega_kp(1, :), mat_p_omega_kp(2, :), mat_p_omega(jquad, :), &
828 kpoints, eps_filter_im_time)
829
830 DO icell = 1, SIZE(mat_p_omega, 2)
831 CALL dbcsr_set(mat_p_omega(jquad, icell)%matrix, 0.0_dp)
832 CALL dbcsr_filter(mat_p_omega(jquad, icell)%matrix, 1.0_dp)
833 END DO
834
835 CALL timestop(handle)
836
837 END SUBROUTINE transform_p_from_real_space_to_kpoints
838
839! **************************************************************************************************
840!> \brief ...
841!> \param real_mat_kp ...
842!> \param imag_mat_kp ...
843!> \param mat_real_space ...
844!> \param kpoints ...
845!> \param eps_filter_im_time ...
846!> \param real_mat_real_space ...
847! **************************************************************************************************
848 SUBROUTINE real_space_to_kpoint_transform_rpa(real_mat_kp, imag_mat_kp, mat_real_space, &
849 kpoints, eps_filter_im_time, real_mat_real_space)
850
851 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT) :: real_mat_kp, imag_mat_kp, mat_real_space
852 TYPE(kpoint_type), POINTER :: kpoints
853 REAL(kind=dp), INTENT(IN) :: eps_filter_im_time
854 LOGICAL, INTENT(IN), OPTIONAL :: real_mat_real_space
855
856 CHARACTER(LEN=*), PARAMETER :: routinen = 'real_space_to_kpoint_transform_rpa'
857
858 INTEGER :: handle, i_cell, ik, nkp, num_cells
859 INTEGER, DIMENSION(3) :: cell
860 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell
861 LOGICAL :: my_real_mat_real_space
862 REAL(kind=dp) :: arg, coskl, sinkl
863 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
864 TYPE(dbcsr_type) :: mat_work
865
866 CALL timeset(routinen, handle)
867
868 my_real_mat_real_space = .true.
869 IF (PRESENT(real_mat_real_space)) my_real_mat_real_space = real_mat_real_space
870
871 CALL dbcsr_create(matrix=mat_work, &
872 template=real_mat_kp(1)%matrix, &
873 matrix_type=dbcsr_type_no_symmetry)
874 CALL dbcsr_reserve_all_blocks(mat_work)
875 CALL dbcsr_set(mat_work, 0.0_dp)
876
877 ! this kpoint environme t should be the kpoints for D(it) and X(it) created in init_cell_index_rpa
878 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp)
879
880 NULLIFY (index_to_cell)
881 index_to_cell => kpoints%index_to_cell
882
883 num_cells = SIZE(index_to_cell, 2)
884
885 cpassert(SIZE(mat_real_space) >= num_cells/2 + 1)
886
887 DO ik = 1, nkp
888
889 CALL dbcsr_reserve_all_blocks(real_mat_kp(ik)%matrix)
890 CALL dbcsr_reserve_all_blocks(imag_mat_kp(ik)%matrix)
891
892 CALL dbcsr_set(real_mat_kp(ik)%matrix, 0.0_dp)
893 CALL dbcsr_set(imag_mat_kp(ik)%matrix, 0.0_dp)
894
895 DO i_cell = 1, num_cells/2 + 1
896
897 cell(:) = index_to_cell(:, i_cell)
898
899 arg = real(cell(1), dp)*xkp(1, ik) + real(cell(2), dp)*xkp(2, ik) + real(cell(3), dp)*xkp(3, ik)
900 coskl = cos(twopi*arg)
901 sinkl = sin(twopi*arg)
902
903 IF (my_real_mat_real_space) THEN
904 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
905 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, sinkl)
906 ELSE
907 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, -sinkl)
908 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
909 END IF
910
911 IF (.NOT. (cell(1) == 0 .AND. cell(2) == 0 .AND. cell(3) == 0)) THEN
912
913 CALL dbcsr_transposed(mat_work, mat_real_space(i_cell)%matrix)
914
915 IF (my_real_mat_real_space) THEN
916 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, coskl)
917 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
918 ELSE
919 ! for an imaginary real-space matrix, we need to consider the imaginary unit
920 ! and we need to take into account that the transposed gives an extra "-" sign
921 ! because the transposed is actually Hermitian conjugate
922 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
923 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -coskl)
924 END IF
925
926 CALL dbcsr_set(mat_work, 0.0_dp)
927
928 END IF
929
930 END DO
931
932 CALL dbcsr_filter(real_mat_kp(ik)%matrix, eps_filter_im_time)
933 CALL dbcsr_filter(imag_mat_kp(ik)%matrix, eps_filter_im_time)
934
935 END DO
936
937 CALL dbcsr_release(mat_work)
938
939 CALL timestop(handle)
940
942
943! **************************************************************************************************
944!> \brief ...
945!> \param mat_a ...
946!> \param mat_b ...
947!> \param alpha ...
948!> \param beta ...
949! **************************************************************************************************
950 SUBROUTINE dbcsr_add_local(mat_a, mat_b, alpha, beta)
951 TYPE(dbcsr_type), INTENT(INOUT) :: mat_a, mat_b
952 REAL(kind=dp), INTENT(IN) :: alpha, beta
953
954 INTEGER :: col, row
955 LOGICAL :: found
956 REAL(kind=dp), DIMENSION(:, :), POINTER :: block_to_compute, data_block
957 TYPE(dbcsr_iterator_type) :: iter
958
959 CALL dbcsr_iterator_start(iter, mat_b)
960 DO WHILE (dbcsr_iterator_blocks_left(iter))
961 CALL dbcsr_iterator_next_block(iter, row, col, data_block)
962
963 NULLIFY (block_to_compute)
964 CALL dbcsr_get_block_p(matrix=mat_a, &
965 row=row, col=col, block=block_to_compute, found=found)
966
967 cpassert(found)
968
969 block_to_compute(:, :) = alpha*block_to_compute(:, :) + beta*data_block(:, :)
970
971 END DO
972 CALL dbcsr_iterator_stop(iter)
973
974 END SUBROUTINE dbcsr_add_local
975
976! **************************************************************************************************
977!> \brief ...
978!> \param fm_mat_W_tau ...
979!> \param cfm_mat_Q ...
980!> \param fm_mat_L_re ...
981!> \param fm_mat_L_im ...
982!> \param dimen_RI ...
983!> \param num_integ_points ...
984!> \param jquad ...
985!> \param ikp ...
986!> \param tj ...
987!> \param tau_tj ...
988!> \param weights_cos_tf_w_to_t ...
989!> \param ikp_local ...
990!> \param para_env ...
991!> \param kpoints ...
992!> \param qs_env ...
993!> \param wkp_W ...
994! **************************************************************************************************
995 SUBROUTINE compute_wc_real_space_tau_gw(fm_mat_W_tau, cfm_mat_Q, fm_mat_L_re, fm_mat_L_im, &
996 dimen_RI, num_integ_points, jquad, &
997 ikp, tj, tau_tj, weights_cos_tf_w_to_t, ikp_local, &
998 para_env, kpoints, qs_env, wkp_W)
999
1000 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_w_tau
1001 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
1002 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_l_re, fm_mat_l_im
1003 INTEGER, INTENT(IN) :: dimen_ri, num_integ_points, jquad, ikp
1004 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: tj
1005 REAL(kind=dp), DIMENSION(num_integ_points), &
1006 INTENT(IN) :: tau_tj
1007 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: weights_cos_tf_w_to_t
1008 INTEGER, INTENT(IN) :: ikp_local
1009 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1010 TYPE(kpoint_type), INTENT(IN), POINTER :: kpoints
1011 TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
1012 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: wkp_w
1013
1014 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_Wc_real_space_tau_GW'
1015
1016 INTEGER :: handle, handle2, i_global, iatom, iatom_old, iib, iquad, irow, j_global, jatom, &
1017 jatom_old, jcol, jjb, jkp, ncol_local, nkp, nrow_local, num_cells
1018 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_from_ri_index
1019 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
1020 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell
1021 REAL(kind=dp) :: contribution, omega, tau, weight, &
1022 weight_im, weight_re
1023 REAL(kind=dp), DIMENSION(3, 3) :: hmat
1024 REAL(kind=dp), DIMENSION(:), POINTER :: wkp
1025 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1026 TYPE(cell_type), POINTER :: cell
1027 TYPE(cp_cfm_type) :: cfm_mat_l, cfm_mat_work, cfm_mat_work_2
1028 TYPE(cp_fm_type) :: fm_dummy, fm_mat_work_global, &
1029 fm_mat_work_local
1030 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1031
1032 CALL timeset(routinen, handle)
1033
1034 CALL timeset(routinen//"_1", handle2)
1035
1036 CALL cp_cfm_create(cfm_mat_work, cfm_mat_q%matrix_struct)
1037 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
1038
1039 CALL cp_cfm_create(cfm_mat_work_2, cfm_mat_q%matrix_struct)
1040 CALL cp_cfm_set_all(cfm_mat_work_2, z_zero)
1041
1042 CALL cp_cfm_create(cfm_mat_l, cfm_mat_q%matrix_struct)
1043 CALL cp_cfm_set_all(cfm_mat_l, z_zero)
1044
1045 ! Copy fm_mat_L_re and fm_mat_L_re to cfm_mat_L
1046 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_l, z_one, fm_mat_l_re)
1047 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_l, gaussi, fm_mat_l_im)
1048
1049 CALL cp_fm_create(fm_mat_work_global, fm_mat_w_tau(1)%matrix_struct)
1050 CALL cp_fm_set_all(fm_mat_work_global, 0.0_dp)
1051
1052 CALL cp_fm_create(fm_mat_work_local, cfm_mat_q%matrix_struct)
1053 CALL cp_fm_set_all(fm_mat_work_local, 0.0_dp)
1054
1055 CALL timestop(handle2)
1056
1057 CALL timeset(routinen//"_2", handle2)
1058
1059 ! calculate [1+Q(iw')]^-1
1060 CALL cp_cfm_cholesky_invert(cfm_mat_q)
1061
1062 ! symmetrize the result
1063 CALL cp_cfm_uplo_to_full(cfm_mat_q)
1064
1065 ! subtract exchange part by subtracing identity matrix from epsilon
1066 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
1067 nrow_local=nrow_local, &
1068 ncol_local=ncol_local, &
1069 row_indices=row_indices, &
1070 col_indices=col_indices)
1071
1072 DO jjb = 1, ncol_local
1073 j_global = col_indices(jjb)
1074 DO iib = 1, nrow_local
1075 i_global = row_indices(iib)
1076 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
1077 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) - z_one
1078 END IF
1079 END DO
1080 END DO
1081
1082 CALL timestop(handle2)
1083
1084 CALL timeset(routinen//"_3", handle2)
1085
1086 ! work = epsilon(iw,k)*V^1/2(k)
1087 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_q, cfm_mat_l, &
1088 z_zero, cfm_mat_work)
1089
1090 ! W(iw,k) = V^1/2(k)*work
1091 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_l, cfm_mat_work, &
1092 z_zero, cfm_mat_work_2)
1093
1094 CALL timestop(handle2)
1095
1096 CALL timeset(routinen//"_4", handle2)
1097
1098 CALL get_kpoint_info(kpoints, xkp=xkp, wkp=wkp, nkp=nkp)
1099 index_to_cell => kpoints%index_to_cell
1100 num_cells = SIZE(index_to_cell, 2)
1101
1102 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
1103
1104 ALLOCATE (atom_from_ri_index(dimen_ri))
1105
1106 CALL get_atom_index_from_basis_function_index(qs_env, atom_from_ri_index, dimen_ri, "RI_AUX")
1107
1108 NULLIFY (cell, particle_set)
1109 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1110 CALL get_cell(cell=cell, h=hmat)
1111 iatom_old = 0
1112 jatom_old = 0
1113
1114 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
1115 nrow_local=nrow_local, &
1116 ncol_local=ncol_local, &
1117 row_indices=row_indices, &
1118 col_indices=col_indices)
1119
1120 DO irow = 1, nrow_local
1121 DO jcol = 1, ncol_local
1122
1123 iatom = atom_from_ri_index(row_indices(irow))
1124 jatom = atom_from_ri_index(col_indices(jcol))
1125
1126 IF (iatom .NE. iatom_old .OR. jatom .NE. jatom_old) THEN
1127
1128 ! symmetrize=.FALSE. necessary since we already have a symmetrized index_to_cell
1129 CALL compute_weight_re_im(weight_re, weight_im, &
1130 num_cells, iatom, jatom, xkp(1:3, ikp), wkp_w(ikp), &
1131 cell, index_to_cell, hmat, particle_set)
1132
1133 iatom_old = iatom
1134 jatom_old = jatom
1135
1136 END IF
1137
1138 contribution = weight_re*real(cfm_mat_work_2%local_data(irow, jcol)) + &
1139 weight_im*aimag(cfm_mat_work_2%local_data(irow, jcol))
1140
1141 fm_mat_work_local%local_data(irow, jcol) = fm_mat_work_local%local_data(irow, jcol) + contribution
1142
1143 END DO
1144 END DO
1145
1146 CALL timestop(handle2)
1147
1148 CALL timeset(routinen//"_5", handle2)
1149
1150 IF (ikp_local == -1) THEN
1151
1152 CALL cp_fm_copy_general(fm_mat_work_local, fm_mat_work_global, para_env)
1153
1154 DO iquad = 1, num_integ_points
1155
1156 omega = tj(jquad)
1157 tau = tau_tj(iquad)
1158 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1159
1160 IF (jquad == 1 .AND. ikp == 1) THEN
1161 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1162 END IF
1163
1164 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_w_tau(iquad), beta=weight, matrix_b=fm_mat_work_global)
1165
1166 END DO
1167
1168 ELSE
1169
1170 DO jkp = 1, nkp
1171
1172 CALL para_env%sync()
1173
1174 IF (ikp_local == jkp) THEN
1175 CALL cp_fm_copy_general(fm_mat_work_local, fm_mat_work_global, para_env)
1176 ELSE
1177 CALL cp_fm_copy_general(fm_dummy, fm_mat_work_global, para_env)
1178 END IF
1179
1180 CALL para_env%sync()
1181
1182 DO iquad = 1, num_integ_points
1183
1184 omega = tj(jquad)
1185 tau = tau_tj(iquad)
1186 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1187
1188 IF (jquad == 1 .AND. jkp == 1) THEN
1189 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1190 END IF
1191
1192 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_w_tau(iquad), beta=weight, &
1193 matrix_b=fm_mat_work_global)
1194
1195 END DO
1196
1197 END DO
1198
1199 END IF
1200
1201 CALL cp_cfm_release(cfm_mat_work)
1202 CALL cp_cfm_release(cfm_mat_work_2)
1203 CALL cp_cfm_release(cfm_mat_l)
1204 CALL cp_fm_release(fm_mat_work_global)
1205 CALL cp_fm_release(fm_mat_work_local)
1206
1207 DEALLOCATE (atom_from_ri_index)
1208
1209 CALL timestop(handle2)
1210
1211 CALL timestop(handle)
1212
1213 END SUBROUTINE compute_wc_real_space_tau_gw
1214
1215! **************************************************************************************************
1216!> \brief ...
1217!> \param fm_mat_W ...
1218!> \param fm_matrix_Minv ...
1219!> \param para_env ...
1220!> \param dimen_RI ...
1221!> \param num_integ_points ...
1222! **************************************************************************************************
1223 SUBROUTINE wc_to_minv_wc_minv(fm_mat_W, fm_matrix_Minv, para_env, dimen_RI, num_integ_points)
1224 TYPE(cp_fm_type), DIMENSION(:) :: fm_mat_w
1225 TYPE(cp_fm_type), DIMENSION(:, :) :: fm_matrix_minv
1226 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1227 INTEGER :: dimen_ri, num_integ_points
1228
1229 CHARACTER(LEN=*), PARAMETER :: routinen = 'Wc_to_Minv_Wc_Minv'
1230
1231 INTEGER :: handle, jquad
1232 TYPE(cp_fm_type) :: fm_work_minv, fm_work_minv_w
1233
1234 CALL timeset(routinen, handle)
1235
1236 CALL cp_fm_create(fm_work_minv, fm_mat_w(1)%matrix_struct)
1237 CALL cp_fm_copy_general(fm_matrix_minv(1, 1), fm_work_minv, para_env)
1238
1239 CALL cp_fm_create(fm_work_minv_w, fm_mat_w(1)%matrix_struct)
1240
1241 DO jquad = 1, num_integ_points
1242
1243 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv, fm_mat_w(jquad), &
1244 0.0_dp, fm_work_minv_w)
1245 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv_w, fm_work_minv, &
1246 0.0_dp, fm_mat_w(jquad))
1247
1248 END DO
1249
1250 CALL cp_fm_release(fm_work_minv)
1251
1252 CALL cp_fm_release(fm_work_minv_w)
1253
1254 CALL timestop(handle)
1255
1256 END SUBROUTINE wc_to_minv_wc_minv
1257
1258! **************************************************************************************************
1259!> \brief ...
1260!> \param qs_env ...
1261!> \param wkp_W ...
1262!> \param wkp_V ...
1263!> \param kpoints ...
1264!> \param h_inv ...
1265!> \param periodic ...
1266! **************************************************************************************************
1267 SUBROUTINE compute_wkp_w(qs_env, wkp_W, wkp_V, kpoints, h_inv, periodic)
1268
1269 TYPE(qs_environment_type), POINTER :: qs_env
1270 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
1271 INTENT(OUT) :: wkp_w, wkp_v
1272 TYPE(kpoint_type), POINTER :: kpoints
1273 REAL(kind=dp), DIMENSION(3, 3) :: h_inv
1274 INTEGER, DIMENSION(3) :: periodic
1275
1276 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_wkp_W'
1277
1278 INTEGER :: handle, i_x, ikp, info, j_y, k_z, &
1279 kpoint_weights_w_method, n_x, n_y, &
1280 n_z, nkp, nsuperfine, num_lin_eqs
1281 REAL(kind=dp) :: exp_kpoints, integral, k_sq, weight
1282 REAL(kind=dp), DIMENSION(3) :: k_vec, x_vec
1283 REAL(kind=dp), DIMENSION(:), POINTER :: right_side, wkp, wkp_tmp
1284 REAL(kind=dp), DIMENSION(:, :), POINTER :: matrix_lin_eqs, xkp
1285
1286 CALL timeset(routinen, handle)
1287
1288 kpoint_weights_w_method = qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method
1289
1290 CALL get_kpoint_info(kpoints, xkp=xkp, wkp=wkp, nkp=nkp)
1291
1292 ! we determine the kpoint weights of the Monkhors Pack mesh new
1293 ! such that the functions 1/k^2, 1/k and const are integrated exactly
1294 ! in the Brillouin zone
1295 ! this is done by minimizing sum_i |w_i|^2 where w_i are the weights of
1296 ! the i-th kpoint under the following constraints:
1297 ! 1) 1/k^2, 1/k and const are integrated exactly
1298 ! 2) the kpoint weights of kpoints with identical absolute value are
1299 ! the same, of e.g. (1/8,3/8,3/8) same weight as for (3/8,1/8,3/8)
1300 ! for 1d and 2d materials: we use ordinary Monkhorst-Pack weights, checked
1301 ! by SUM(periodic) == 3
1302 ALLOCATE (wkp_v(nkp), wkp_w(nkp))
1303
1304 ! for exchange part of self-energy, we use truncated Coulomb operator that should be fine
1305 ! with uniform weights (without k-point extrapolation)
1306 IF (ALLOCATED(qs_env%mp2_env%ri_rpa_im_time%wkp_V)) THEN
1307 wkp_v(:) = qs_env%mp2_env%ri_rpa_im_time%wkp_V(:)
1308 ELSE
1309 wkp_v(:) = wkp(:)
1310 END IF
1311
1312 IF (kpoint_weights_w_method == kp_weights_w_uniform) THEN
1313
1314 ! in the k-point weights wkp, there might be k-point extrapolation included
1315 wkp_w(:) = wkp(:)
1316
1317 ELSE IF (kpoint_weights_w_method == kp_weights_w_tailored .OR. &
1318 kpoint_weights_w_method == kp_weights_w_auto) THEN
1319
1320 IF (kpoint_weights_w_method == kp_weights_w_tailored) &
1321 exp_kpoints = qs_env%mp2_env%ri_rpa_im_time%exp_tailored_weights
1322
1323 IF (kpoint_weights_w_method == kp_weights_w_auto) THEN
1324 IF (sum(periodic) == 2) exp_kpoints = -1.0_dp
1325 END IF
1326
1327 ! first, compute the integral of f(k)=1/k^2 and 1/k on super fine grid
1328 nsuperfine = 500
1329 integral = 0.0_dp
1330
1331 IF (periodic(1) == 1) THEN
1332 n_x = nsuperfine
1333 ELSE
1334 n_x = 1
1335 END IF
1336 IF (periodic(2) == 1) THEN
1337 n_y = nsuperfine
1338 ELSE
1339 n_y = 1
1340 END IF
1341 IF (periodic(3) == 1) THEN
1342 n_z = nsuperfine
1343 ELSE
1344 n_z = 1
1345 END IF
1346
1347 ! actually, there is the factor *det_3x3(h_inv) missing to account for the
1348 ! integration volume but for wkp det_3x3(h_inv) is needed
1349 weight = 1.0_dp/(real(n_x, dp)*real(n_y, dp)*real(n_z, dp))
1350 DO i_x = 1, n_x
1351 DO j_y = 1, n_y
1352 DO k_z = 1, n_z
1353
1354 IF (periodic(1) == 1) THEN
1355 x_vec(1) = (real(i_x - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1356 ELSE
1357 x_vec(1) = 0.0_dp
1358 END IF
1359 IF (periodic(2) == 1) THEN
1360 x_vec(2) = (real(j_y - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1361 ELSE
1362 x_vec(2) = 0.0_dp
1363 END IF
1364 IF (periodic(3) == 1) THEN
1365 x_vec(3) = (real(k_z - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1366 ELSE
1367 x_vec(3) = 0.0_dp
1368 END IF
1369
1370 k_vec = matmul(h_inv(1:3, 1:3), x_vec)
1371 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1372 integral = integral + weight*k_sq**(exp_kpoints*0.5_dp)
1373
1374 END DO
1375 END DO
1376 END DO
1377
1378 num_lin_eqs = nkp + 2
1379
1380 ALLOCATE (matrix_lin_eqs(num_lin_eqs, num_lin_eqs))
1381 matrix_lin_eqs(:, :) = 0.0_dp
1382
1383 DO ikp = 1, nkp
1384
1385 k_vec = matmul(h_inv(1:3, 1:3), xkp(1:3, ikp))
1386 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1387
1388 matrix_lin_eqs(ikp, ikp) = 2.0_dp
1389 matrix_lin_eqs(ikp, nkp + 1) = 1.0_dp
1390 matrix_lin_eqs(nkp + 1, ikp) = 1.0_dp
1391
1392 matrix_lin_eqs(ikp, nkp + 2) = k_sq**(exp_kpoints*0.5_dp)
1393 matrix_lin_eqs(nkp + 2, ikp) = k_sq**(exp_kpoints*0.5_dp)
1394
1395 END DO
1396
1397 CALL invmat(matrix_lin_eqs, info)
1398 ! check whether inversion was successful
1399 cpassert(info == 0)
1400
1401 ALLOCATE (right_side(num_lin_eqs))
1402 right_side = 0.0_dp
1403 right_side(nkp + 1) = 1.0_dp
1404 ! divide integral by two because CP2K k-mesh already considers symmetry k <-> -k
1405 right_side(nkp + 2) = integral
1406
1407 ALLOCATE (wkp_tmp(num_lin_eqs))
1408
1409 wkp_tmp(1:num_lin_eqs) = matmul(matrix_lin_eqs, right_side)
1410
1411 wkp_w(1:nkp) = wkp_tmp(1:nkp)
1412
1413 DEALLOCATE (matrix_lin_eqs, right_side, wkp_tmp)
1414
1415 END IF
1416
1417 CALL timestop(handle)
1418
1419 END SUBROUTINE
1420
1421! **************************************************************************************************
1422!> \brief ...
1423!> \param qs_env ...
1424!> \param Eigenval_kp ...
1425! **************************************************************************************************
1426 SUBROUTINE get_bandstruc_and_k_dependent_mos(qs_env, Eigenval_kp)
1427 TYPE(qs_environment_type), POINTER :: qs_env
1428 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: eigenval_kp
1429
1430 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_bandstruc_and_k_dependent_MOs'
1431
1432 INTEGER :: handle, ikp, ispin, nmo, nspins
1433 INTEGER, DIMENSION(3) :: nkp_grid_g
1434 REAL(kind=dp), DIMENSION(:), POINTER :: ev
1435 REAL(kind=dp), DIMENSION(:, :), POINTER :: kpgeneral
1436 TYPE(kpoint_type), POINTER :: kpoints_sigma
1437 TYPE(mp_para_env_type), POINTER :: para_env
1438
1439 CALL timeset(routinen, handle)
1440
1441 NULLIFY (qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1442 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1443 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1444 para_env)
1445
1446 nkp_grid_g(1:3) = (/1, 1, 1/)
1447
1448 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
1449
1450 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1451 "MONKHORST-PACK", para_env%num_pe, &
1452 mp_grid=nkp_grid_g(1:3))
1453
1454 IF (qs_env%mp2_env%ri_g0w0%do_kpoints_Sigma) THEN
1455
1456 ! set up k-points for GW band structure calculation, will be completed later
1457 CALL get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1458
1459 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1460 "GENERAL", para_env%num_pe, &
1461 kpgeneral=kpgeneral)
1462
1463 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1464 "GENERAL", para_env%num_pe, &
1465 kpgeneral=kpgeneral, with_xc_terms=.false.)
1466
1467 kpoints_sigma => qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma
1468 nmo = SIZE(eigenval_kp, 1)
1469 nspins = SIZE(eigenval_kp, 3)
1470
1471 ALLOCATE (qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(nmo))
1472 qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(:) = eigenval_kp(:, 1, 1)
1473
1474 DEALLOCATE (eigenval_kp)
1475
1476 ALLOCATE (eigenval_kp(nmo, kpoints_sigma%nkp, nspins))
1477
1478 DO ikp = 1, kpoints_sigma%nkp
1479
1480 DO ispin = 1, nspins
1481
1482 ev => kpoints_sigma%kp_env(ikp)%kpoint_env%mos(1, ispin)%eigenvalues
1483
1484 eigenval_kp(:, ikp, ispin) = ev(:)
1485
1486 END DO
1487
1488 END DO
1489
1490 DEALLOCATE (kpgeneral)
1491
1492 END IF
1493
1494 CALL release_hfx_stuff(qs_env)
1495
1496 CALL timestop(handle)
1497
1499
1500! **************************************************************************************************
1501!> \brief releases part of the given qs_env in order to save memory
1502!> \param qs_env the object to release
1503! **************************************************************************************************
1504 SUBROUTINE release_hfx_stuff(qs_env)
1505 TYPE(qs_environment_type), POINTER :: qs_env
1506
1507 IF (ASSOCIATED(qs_env%x_data) .AND. .NOT. qs_env%mp2_env%ri_g0w0%do_ri_Sigma_x) THEN
1508 CALL hfx_release(qs_env%x_data)
1509 END IF
1510
1511 END SUBROUTINE release_hfx_stuff
1512
1513! **************************************************************************************************
1514!> \brief ...
1515!> \param qs_env ...
1516!> \param kpoints ...
1517!> \param scheme ...
1518!> \param group_size_ext ...
1519!> \param mp_grid ...
1520!> \param kpgeneral ...
1521!> \param with_xc_terms ...
1522! **************************************************************************************************
1523 SUBROUTINE create_kp_and_calc_kp_orbitals(qs_env, kpoints, scheme, &
1524 group_size_ext, mp_grid, kpgeneral, with_xc_terms)
1525
1526 TYPE(qs_environment_type), POINTER :: qs_env
1527 TYPE(kpoint_type), POINTER :: kpoints
1528 CHARACTER(LEN=*), INTENT(IN) :: scheme
1529 INTEGER :: group_size_ext
1530 INTEGER, DIMENSION(3), INTENT(IN), OPTIONAL :: mp_grid
1531 REAL(kind=dp), DIMENSION(:, :), INTENT(IN), &
1532 OPTIONAL :: kpgeneral
1533 LOGICAL, OPTIONAL :: with_xc_terms
1534
1535 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_kp_and_calc_kp_orbitals'
1536
1537 INTEGER :: handle, i_dim, i_re_im, ikp, ispin, nkp, &
1538 nspins
1539 INTEGER, DIMENSION(3) :: cell_grid, periodic
1540 LOGICAL :: my_with_xc_terms
1541 REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues
1542 TYPE(cell_type), POINTER :: cell
1543 TYPE(cp_blacs_env_type), POINTER :: blacs_env
1544 TYPE(cp_cfm_type) :: cksmat, cmos, csmat, cwork
1545 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1546 TYPE(cp_fm_type) :: fm_work
1547 TYPE(cp_fm_type), POINTER :: imos, rmos
1548 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, matrix_s_desymm
1549 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_ks_kp, mat_s_kp
1550 TYPE(dft_control_type), POINTER :: dft_control
1551 TYPE(kpoint_env_type), POINTER :: kp
1552 TYPE(mp_para_env_type), POINTER :: para_env
1553 TYPE(qs_scf_env_type), POINTER :: scf_env
1554 TYPE(scf_control_type), POINTER :: scf_control
1555
1556 CALL timeset(routinen, handle)
1557
1558 my_with_xc_terms = .true.
1559 IF (PRESENT(with_xc_terms)) my_with_xc_terms = with_xc_terms
1560
1561 CALL get_qs_env(qs_env, &
1562 para_env=para_env, &
1563 blacs_env=blacs_env, &
1564 matrix_s=matrix_s, &
1565 scf_env=scf_env, &
1566 scf_control=scf_control, &
1567 cell=cell)
1568
1569 ! get kpoints
1570 CALL calculate_kpoints_for_bs(kpoints, scheme, kpgeneral=kpgeneral, mp_grid=mp_grid, &
1571 group_size_ext=group_size_ext)
1572
1573 CALL kpoint_env_initialize(kpoints, para_env, blacs_env)
1574
1575 ! calculate all MOs that are accessible in the given
1576 ! Gaussian AO basis, therefore nadd=1E10
1577 CALL kpoint_initialize_mos(kpoints, qs_env%mos, 2000000000)
1578 CALL kpoint_initialize_mo_set(kpoints)
1579
1580 CALL get_cell(cell=cell, periodic=periodic)
1581
1582 DO i_dim = 1, 3
1583 ! we have at most 3 neigboring cells per dimension and at least one because
1584 ! the density response at Gamma is only divided to neighboring
1585 IF (periodic(i_dim) == 1) THEN
1586 cell_grid(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
1587 ELSE
1588 cell_grid(i_dim) = 1
1589 END IF
1590 END DO
1591 CALL init_cell_index_rpa(cell_grid, kpoints%cell_to_index, kpoints%index_to_cell, cell)
1592
1593 ! get S(k)
1594 CALL get_qs_env(qs_env, matrix_s=matrix_s, scf_env=scf_env, scf_control=scf_control, dft_control=dft_control)
1595
1596 NULLIFY (matrix_s_desymm)
1597 CALL dbcsr_allocate_matrix_set(matrix_s_desymm, 1)
1598 ALLOCATE (matrix_s_desymm(1)%matrix)
1599 CALL dbcsr_create(matrix=matrix_s_desymm(1)%matrix, template=matrix_s(1)%matrix, &
1600 matrix_type=dbcsr_type_no_symmetry)
1601 CALL dbcsr_desymmetrize(matrix_s(1)%matrix, matrix_s_desymm(1)%matrix)
1602
1603 CALL mat_kp_from_mat_gamma(qs_env, mat_s_kp, matrix_s_desymm(1)%matrix, kpoints, 1)
1604
1605 CALL get_kpoint_info(kpoints, nkp=nkp)
1606
1607 matrix_struct => kpoints%kp_env(1)%kpoint_env%wmat(1, 1)%matrix_struct
1608
1609 CALL cp_cfm_create(cksmat, matrix_struct)
1610 CALL cp_cfm_create(csmat, matrix_struct)
1611 CALL cp_cfm_create(cmos, matrix_struct)
1612 CALL cp_cfm_create(cwork, matrix_struct)
1613 CALL cp_fm_create(fm_work, matrix_struct)
1614
1615 nspins = dft_control%nspins
1616
1617 DO ispin = 1, nspins
1618
1619 ! get H(k)
1620 IF (my_with_xc_terms) THEN
1621 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_ks(ispin)%matrix, kpoints, ispin)
1622 ELSE
1623 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_sigma_x_minus_vxc(ispin)%matrix, &
1624 kpoints, ispin)
1625 END IF
1626
1627 DO ikp = 1, nkp
1628
1629 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 1)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1630 CALL cp_cfm_scale_and_add_fm(z_zero, cksmat, z_one, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1631
1632 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 2)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1633 CALL cp_cfm_scale_and_add_fm(z_one, cksmat, gaussi, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1634
1635 CALL copy_dbcsr_to_fm(mat_s_kp(ikp, 1)%matrix, fm_work)
1636 CALL cp_cfm_scale_and_add_fm(z_zero, csmat, z_one, fm_work)
1637
1638 CALL copy_dbcsr_to_fm(mat_s_kp(ikp, 2)%matrix, fm_work)
1639 CALL cp_cfm_scale_and_add_fm(z_one, csmat, gaussi, fm_work)
1640
1641 kp => kpoints%kp_env(ikp)%kpoint_env
1642
1643 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos, eigenvalues=eigenvalues)
1644 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1645
1646 IF (scf_env%cholesky_method == cholesky_off .OR. &
1647 qs_env%mp2_env%ri_rpa_im_time%make_overlap_mat_ao_pos_definite) THEN
1648 CALL cp_cfm_geeig_canon(cksmat, csmat, cmos, eigenvalues, cwork, scf_control%eps_eigval)
1649 ELSE
1650 CALL cp_cfm_geeig(cksmat, csmat, cmos, eigenvalues, cwork)
1651 END IF
1652
1653 CALL cp_cfm_to_fm(cmos, rmos, imos)
1654
1655 kp%mos(2, ispin)%eigenvalues = eigenvalues
1656
1657 END DO
1658
1659 END DO
1660
1661 DO ikp = 1, nkp
1662 DO i_re_im = 1, 2
1663 CALL dbcsr_deallocate_matrix(mat_ks_kp(ikp, i_re_im)%matrix)
1664 END DO
1665 END DO
1666 DEALLOCATE (mat_ks_kp)
1667
1668 DO ikp = 1, nkp
1669 DO i_re_im = 1, 2
1670 CALL dbcsr_deallocate_matrix(mat_s_kp(ikp, i_re_im)%matrix)
1671 END DO
1672 END DO
1673 DEALLOCATE (mat_s_kp)
1674
1675 CALL dbcsr_deallocate_matrix(matrix_s_desymm(1)%matrix)
1676 DEALLOCATE (matrix_s_desymm)
1677
1678 CALL cp_cfm_release(cksmat)
1679 CALL cp_cfm_release(csmat)
1680 CALL cp_cfm_release(cwork)
1681 CALL cp_cfm_release(cmos)
1682 CALL cp_fm_release(fm_work)
1683
1684 CALL timestop(handle)
1685
1686 END SUBROUTINE create_kp_and_calc_kp_orbitals
1687
1688! **************************************************************************************************
1689!> \brief ...
1690!> \param qs_env ...
1691!> \param mat_kp ...
1692!> \param mat_gamma ...
1693!> \param kpoints ...
1694!> \param ispin ...
1695!> \param real_mat_real_space ...
1696! **************************************************************************************************
1697 SUBROUTINE mat_kp_from_mat_gamma(qs_env, mat_kp, mat_gamma, kpoints, ispin, real_mat_real_space)
1698
1699 TYPE(qs_environment_type), POINTER :: qs_env
1700 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_kp
1701 TYPE(dbcsr_type) :: mat_gamma
1702 TYPE(kpoint_type), POINTER :: kpoints
1703 INTEGER :: ispin
1704 LOGICAL, INTENT(IN), OPTIONAL :: real_mat_real_space
1705
1706 CHARACTER(LEN=*), PARAMETER :: routinen = 'mat_kp_from_mat_gamma'
1707
1708 INTEGER :: handle, i_cell, i_re_im, ikp, nkp, &
1709 num_cells
1710 INTEGER, DIMENSION(3) :: periodic
1711 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1712 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1713 TYPE(cell_type), POINTER :: cell
1714 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mat_real_space
1715
1716 CALL timeset(routinen, handle)
1717
1718 CALL get_qs_env(qs_env, cell=cell)
1719 CALL get_cell(cell=cell, periodic=periodic)
1720 num_cells = 3**(periodic(1) + periodic(2) + periodic(3))
1721
1722 NULLIFY (mat_real_space)
1723 CALL dbcsr_allocate_matrix_set(mat_real_space, num_cells)
1724 DO i_cell = 1, num_cells
1725 ALLOCATE (mat_real_space(i_cell)%matrix)
1726 CALL dbcsr_create(matrix=mat_real_space(i_cell)%matrix, &
1727 template=mat_gamma)
1728 CALL dbcsr_reserve_all_blocks(mat_real_space(i_cell)%matrix)
1729 CALL dbcsr_set(mat_real_space(i_cell)%matrix, 0.0_dp)
1730 END DO
1731
1732 CALL dbcsr_copy(mat_real_space(1)%matrix, mat_gamma)
1733
1734 CALL get_mat_cell_t_from_mat_gamma(mat_real_space, qs_env, kpoints, 2, 0)
1735
1736 NULLIFY (xkp, cell_to_index)
1737 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp, cell_to_index=cell_to_index)
1738
1739 IF (ispin == 1) THEN
1740 NULLIFY (mat_kp)
1741 CALL dbcsr_allocate_matrix_set(mat_kp, nkp, 2)
1742 DO ikp = 1, nkp
1743 DO i_re_im = 1, 2
1744 ALLOCATE (mat_kp(ikp, i_re_im)%matrix)
1745 CALL dbcsr_create(matrix=mat_kp(ikp, i_re_im)%matrix, template=mat_gamma)
1746 CALL dbcsr_reserve_all_blocks(mat_kp(ikp, i_re_im)%matrix)
1747 CALL dbcsr_set(mat_kp(ikp, i_re_im)%matrix, 0.0_dp)
1748 END DO
1749 END DO
1750 END IF
1751
1752 IF (PRESENT(real_mat_real_space)) THEN
1753 CALL real_space_to_kpoint_transform_rpa(mat_kp(:, 1), mat_kp(:, 2), mat_real_space, kpoints, 0.0_dp, &
1754 real_mat_real_space)
1755 ELSE
1756 CALL real_space_to_kpoint_transform_rpa(mat_kp(:, 1), mat_kp(:, 2), mat_real_space, kpoints, 0.0_dp)
1757 END IF
1758
1759 DO i_cell = 1, num_cells
1760 CALL dbcsr_deallocate_matrix(mat_real_space(i_cell)%matrix)
1761 END DO
1762 DEALLOCATE (mat_real_space)
1763
1764 CALL timestop(handle)
1765
1766 END SUBROUTINE mat_kp_from_mat_gamma
1767
1768! **************************************************************************************************
1769!> \brief ...
1770!> \param qs_env ...
1771!> \param kpgeneral ...
1772! **************************************************************************************************
1773 SUBROUTINE get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1774 TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
1775 REAL(kind=dp), DIMENSION(:, :), POINTER :: kpgeneral
1776
1777 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_kpgeneral_for_Sigma_kpoints'
1778
1779 INTEGER :: handle, i_kp_in_kp_line, i_special_kp, &
1780 i_x, ikk, j_y, k_z, n_kp_in_kp_line, &
1781 n_special_kp
1782 INTEGER, DIMENSION(:), POINTER :: nkp_grid
1783
1784 CALL timeset(routinen, handle)
1785
1786 n_special_kp = qs_env%mp2_env%ri_g0w0%n_special_kp
1787 n_kp_in_kp_line = qs_env%mp2_env%ri_g0w0%n_kp_in_kp_line
1788 IF (n_special_kp > 0) THEN
1789 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = n_kp_in_kp_line*(n_special_kp - 1) + 1
1790 ELSE
1791 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = 0
1792 END IF
1793
1794 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack = qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(1)* &
1795 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(2)* &
1796 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(3)
1797
1798 qs_env%mp2_env%ri_g0w0%nkp_self_energy = qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp + &
1799 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack
1800
1801 ALLOCATE (kpgeneral(3, qs_env%mp2_env%ri_g0w0%nkp_self_energy))
1802
1803 IF (n_special_kp > 0) THEN
1804
1805 kpgeneral(1:3, 1) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, 1)
1806
1807 ikk = 1
1808
1809 DO i_special_kp = 2, n_special_kp
1810 DO i_kp_in_kp_line = 1, n_kp_in_kp_line
1811
1812 ikk = ikk + 1
1813 kpgeneral(1:3, ikk) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1) + &
1814 REAL(i_kp_in_kp_line, kind=dp)/real(n_kp_in_kp_line, kind=dp)* &
1815 (qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp) - &
1816 qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1))
1817
1818 END DO
1819 END DO
1820
1821 ELSE
1822
1823 ikk = 0
1824
1825 END IF
1826
1827 nkp_grid => qs_env%mp2_env%ri_g0w0%kp_grid_Sigma
1828
1829 DO i_x = 1, nkp_grid(1)
1830 DO j_y = 1, nkp_grid(2)
1831 DO k_z = 1, nkp_grid(3)
1832 ikk = ikk + 1
1833 kpgeneral(1, ikk) = real(2*i_x - nkp_grid(1) - 1, kind=dp)/(2._dp*real(nkp_grid(1), kind=dp))
1834 kpgeneral(2, ikk) = real(2*j_y - nkp_grid(2) - 1, kind=dp)/(2._dp*real(nkp_grid(2), kind=dp))
1835 kpgeneral(3, ikk) = real(2*k_z - nkp_grid(3) - 1, kind=dp)/(2._dp*real(nkp_grid(3), kind=dp))
1836 END DO
1837 END DO
1838 END DO
1839
1840 CALL timestop(handle)
1841
1842 END SUBROUTINE get_kpgeneral_for_sigma_kpoints
1843
1844END MODULE rpa_gw_kpoints_util
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
methods related to the blacs parallel environment
Basic linear algebra operations for complex full matrices.
subroutine, public cp_cfm_scale_and_add_fm(alpha, matrix_a, beta, matrix_b)
Scale and add two BLACS matrices (a = alpha*a + beta*b). where b is a real matrix (adapted from cp_cf...
subroutine, public cp_cfm_uplo_to_full(matrix, workspace, uplo)
...
subroutine, public cp_cfm_column_scale(matrix_a, scaling)
Scales columns of the full matrix by corresponding factors.
various cholesky decomposition related routines
subroutine, public cp_cfm_cholesky_decompose(matrix, n, info_out)
Used to replace a symmetric positive definite matrix M with its Cholesky decomposition U: M = U^T * U...
subroutine, public cp_cfm_cholesky_invert(matrix, n, info_out)
Used to replace Cholesky decomposition by the inverse.
used for collecting diagonalization schemes available for cp_cfm_type
Definition cp_cfm_diag.F:14
subroutine, public cp_cfm_geeig(amatrix, bmatrix, eigenvectors, eigenvalues, work)
General Eigenvalue Problem AX = BXE Single option version: Cholesky decomposition of B.
subroutine, public cp_cfm_heevd(matrix, eigenvectors, eigenvalues)
Perform a diagonalisation of a complex matrix.
Definition cp_cfm_diag.F:60
subroutine, public cp_cfm_geeig_canon(amatrix, bmatrix, eigenvectors, eigenvalues, work, epseig)
General Eigenvalue Problem AX = BXE Use canonical orthogonalization.
Represents a complex full matrix distributed on many processors.
subroutine, public cp_cfm_create(matrix, matrix_struct, name)
Creates a new full matrix with the given structure.
subroutine, public cp_cfm_release(matrix)
Releases a full matrix.
subroutine, public cp_cfm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, matrix_struct, para_env)
Returns information about a full matrix.
subroutine, public cp_cfm_set_all(matrix, alpha, beta)
Set all elements of the full matrix to alpha. Besides, set all diagonal matrix elements to beta (if g...
subroutine, public cp_cfm_to_fm(msource, mtargetr, mtargeti)
Copy real and imaginary parts of a complex full matrix into separate real-value full matrices.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_transposed(transposed, normal, shallow_data_copy, transpose_distribution, use_distribution)
...
subroutine, public dbcsr_deallocate_matrix(matrix)
...
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_get_block_p(matrix, row, col, block, found, row_size, col_size)
...
subroutine, public dbcsr_filter(matrix, eps)
...
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_reserve_all_blocks(matrix)
Reserves all blocks.
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
subroutine, public copy_fm_to_dbcsr(fm, matrix, keep_sparsity)
Copy a BLACS matrix to a dbcsr matrix.
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....
represent the structure of a full matrix
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_copy_general(source, destination, para_env)
General copy of a fm matrix to another fm matrix. Uses non-blocking MPI rather than ScaLAPACK.
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
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
Types and set/get functions for HFX.
Definition hfx_types.F:15
subroutine, public hfx_release(x_data)
This routine deallocates all data structures
Definition hfx_types.F:1910
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public kp_weights_w_auto
integer, parameter, public kp_weights_w_uniform
integer, parameter, public cholesky_off
integer, parameter, public kp_weights_w_tailored
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize_mo_set(kpoint)
...
subroutine, public kpoint_initialize_mos(kpoint, mos, added_mos, for_aux_fit)
Initialize a set of MOs and density matrix for each kpoint (kpoint group)
subroutine, public kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
Initialize the kpoint environment.
Types and basic routines needed for a kpoint calculation.
subroutine, public get_kpoint_info(kpoint, kp_scheme, nkp_grid, kp_shift, symmetry, verbose, full_grid, use_real_wfn, eps_geo, parallel_group_size, kp_range, nkp, xkp, wkp, para_env, blacs_env_all, para_env_kp, para_env_inter_kp, blacs_env, kp_env, kp_aux_env, mpools, iogrp, nkp_groups, kp_dist, cell_to_index, index_to_cell, sab_nl, sab_nl_nosym)
Retrieve information from a kpoint environment.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:148
Definition of mathematical constants and functions.
complex(kind=dp), parameter, public z_one
complex(kind=dp), parameter, public gaussi
real(kind=dp), parameter, public twopi
complex(kind=dp), parameter, public z_zero
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
subroutine, public invmat(a, info)
returns inverse of matrix using the lapack routines DGETRF and DGETRI
Definition mathlib.F:543
Interface to the message passing library MPI.
basic linear algebra operations for full matrixes
Define the data structure for the particle information.
Calculation of band structures.
subroutine, public calculate_kpoints_for_bs(kpoint, scheme, group_size_ext, mp_grid, kpgeneral)
...
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, 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, 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.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
module that contains the definitions of the scf types
Utility routines for GW with imaginary time.
subroutine, public compute_weight_re_im(weight_re, weight_im, num_cells, iatom, jatom, xkp, wkp_w, cell, index_to_cell, hmat, particle_set)
...
subroutine, public get_atom_index_from_basis_function_index(qs_env, atom_from_basis_index, basis_size, basis_type, first_bf_from_atom)
...
Routines treating GW and RPA calculations with kpoints.
subroutine, public get_mat_cell_t_from_mat_gamma(mat_p_omega, qs_env, kpoints, jquad, unit_nr)
...
subroutine, public invert_eps_compute_w_and_erpa_kp(dimen_ri, num_integ_points, jquad, nkp, count_ev_sc_gw, para_env, erpa, tau_tj, tj, wj, weights_cos_tf_w_to_t, wkp_w, do_gw_im_time, do_ri_sigma_x, do_kpoints_from_gamma, cfm_mat_q, ikp_local, mat_p_omega, mat_p_omega_kp, qs_env, eps_filter_im_time, unit_nr, kpoints, fm_mat_minv_l_kpoints, fm_matrix_l_kpoints, fm_mat_w, fm_mat_ri_global_work, mat_minvvminv, fm_matrix_minv, fm_matrix_minv_vtrunc_minv)
...
subroutine, public real_space_to_kpoint_transform_rpa(real_mat_kp, imag_mat_kp, mat_real_space, kpoints, eps_filter_im_time, real_mat_real_space)
...
subroutine, public compute_wkp_w(qs_env, wkp_w, wkp_v, kpoints, h_inv, periodic)
...
subroutine, public get_bandstruc_and_k_dependent_mos(qs_env, eigenval_kp)
...
subroutine, public cp_cfm_power(matrix, threshold, exponent, min_eigval)
...
subroutine, public mat_kp_from_mat_gamma(qs_env, mat_kp, mat_gamma, kpoints, ispin, real_mat_real_space)
...
Routines for low-scaling RPA/GW with imaginary time.
Definition rpa_im_time.F:13
subroutine, public init_cell_index_rpa(cell_grid, cell_to_index, index_to_cell, cell)
...
parameters that control an scf iteration
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Represent a complex full matrix.
keeps the information about the structure of a full matrix
represent a full matrix
Keeps information about a specific k-point.
Contains information about kpoints.
stores all the informations relevant to an mpi environment