(git:ed6f26b)
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 COMPLEX(KIND=dp), PARAMETER :: czero = cmplx(0.0_dp, 0.0_dp, kind=dp)
280
281 COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues_exponent
282 INTEGER :: handle, i, ncol_global, nrow_global
283 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues
284 TYPE(cp_cfm_type) :: cfm_work
285
286 CALL timeset(routinen, handle)
287
288 CALL cp_cfm_create(cfm_work, matrix%matrix_struct)
289 CALL cp_cfm_set_all(cfm_work, z_zero)
290
291 ! Test that matrix is square
292 CALL cp_cfm_get_info(matrix, nrow_global=nrow_global, ncol_global=ncol_global)
293 cpassert(nrow_global == ncol_global)
294 ALLOCATE (eigenvalues(nrow_global))
295 eigenvalues(:) = 0.0_dp
296 ALLOCATE (eigenvalues_exponent(nrow_global))
297 eigenvalues_exponent(:) = czero
298
299 ! Diagonalize matrix: get eigenvectors and eigenvalues
300 CALL cp_cfm_heevd(matrix, cfm_work, eigenvalues)
301
302 DO i = 1, nrow_global
303 IF (eigenvalues(i) > threshold) THEN
304 eigenvalues_exponent(i) = cmplx((eigenvalues(i))**(0.5_dp*exponent), threshold, kind=dp)
305 ELSE
306 IF (PRESENT(min_eigval)) THEN
307 eigenvalues_exponent(i) = cmplx(min_eigval, 0.0_dp, kind=dp)
308 ELSE
309 eigenvalues_exponent(i) = czero
310 END IF
311 END IF
312 END DO
313
314 CALL cp_cfm_column_scale(cfm_work, eigenvalues_exponent)
315
316 CALL parallel_gemm("N", "C", nrow_global, nrow_global, nrow_global, z_one, &
317 cfm_work, cfm_work, z_zero, matrix)
318
319 DEALLOCATE (eigenvalues, eigenvalues_exponent)
320
321 CALL cp_cfm_release(cfm_work)
322
323 CALL timestop(handle)
324
325 END SUBROUTINE cp_cfm_power
326
327! **************************************************************************************************
328!> \brief ...
329!> \param cfm_mat_Q ...
330!> \param mat_P_omega_kp ...
331!> \param fm_mat_L_re ...
332!> \param fm_mat_L_im ...
333!> \param fm_mat_RI_global_work ...
334!> \param dimen_RI ...
335!> \param ikp ...
336!> \param nkp ...
337!> \param ikp_local ...
338!> \param para_env ...
339!> \param make_chi_pos_definite ...
340! **************************************************************************************************
341 SUBROUTINE compute_q_kp_rpa(cfm_mat_Q, mat_P_omega_kp, fm_mat_L_re, fm_mat_L_im, &
342 fm_mat_RI_global_work, dimen_RI, ikp, nkp, ikp_local, para_env, &
343 make_chi_pos_definite)
344
345 TYPE(cp_cfm_type) :: cfm_mat_q
346 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega_kp
347 TYPE(cp_fm_type) :: fm_mat_l_re, fm_mat_l_im, &
348 fm_mat_ri_global_work
349 INTEGER, INTENT(IN) :: dimen_ri, ikp, nkp, ikp_local
350 TYPE(mp_para_env_type), POINTER :: para_env
351 LOGICAL, INTENT(IN) :: make_chi_pos_definite
352
353 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_Q_kp_RPA'
354
355 INTEGER :: handle
356 TYPE(cp_cfm_type) :: cfm_mat_l, cfm_mat_work
357 TYPE(cp_fm_type) :: fm_mat_work
358
359 CALL timeset(routinen, handle)
360
361 CALL cp_cfm_create(cfm_mat_work, fm_mat_l_re%matrix_struct)
362 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
363
364 CALL cp_cfm_create(cfm_mat_l, fm_mat_l_re%matrix_struct)
365 CALL cp_cfm_set_all(cfm_mat_l, z_zero)
366
367 CALL cp_fm_create(fm_mat_work, fm_mat_l_re%matrix_struct)
368 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
369
370 ! 1. Convert the dbcsr matrix mat_P_omega_kp (that is chi(k,iw)) to a full matrix and
371 ! distribute it to subgroups
372 CALL mat_p_to_subgroup(mat_p_omega_kp, fm_mat_ri_global_work, &
373 fm_mat_work, cfm_mat_q, ikp, nkp, ikp_local, para_env)
374
375 ! 2. Remove all negative eigenvalues from chi(k,iw)
376 IF (make_chi_pos_definite) THEN
377 CALL cp_cfm_power(cfm_mat_q, threshold=0.0_dp, exponent=1.0_dp)
378 END IF
379
380 ! 3. Copy fm_mat_L_re and fm_mat_L_re to cfm_mat_L
381 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_l, z_one, fm_mat_l_re)
382 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_l, gaussi, fm_mat_l_im)
383
384 ! 4. work = P(iw,k)*L(k)
385 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_q, cfm_mat_l, &
386 z_zero, cfm_mat_work)
387
388 ! 5. Q(iw,k) = L^H(k)*work
389 CALL parallel_gemm('C', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_l, cfm_mat_work, &
390 z_zero, cfm_mat_q)
391
392 CALL cp_cfm_release(cfm_mat_work)
393 CALL cp_cfm_release(cfm_mat_l)
394 CALL cp_fm_release(fm_mat_work)
395
396 CALL timestop(handle)
397
398 END SUBROUTINE compute_q_kp_rpa
399
400! **************************************************************************************************
401!> \brief ...
402!> \param mat_P_omega_kp ...
403!> \param fm_mat_RI_global_work ...
404!> \param fm_mat_work ...
405!> \param cfm_mat_Q ...
406!> \param ikp ...
407!> \param nkp ...
408!> \param ikp_local ...
409!> \param para_env ...
410! **************************************************************************************************
411 SUBROUTINE mat_p_to_subgroup(mat_P_omega_kp, fm_mat_RI_global_work, &
412 fm_mat_work, cfm_mat_Q, ikp, nkp, ikp_local, para_env)
413
414 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega_kp
415 TYPE(cp_fm_type), INTENT(INOUT) :: fm_mat_ri_global_work, fm_mat_work
416 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
417 INTEGER, INTENT(IN) :: ikp, nkp, ikp_local
418 TYPE(mp_para_env_type), POINTER :: para_env
419
420 CHARACTER(LEN=*), PARAMETER :: routinen = 'mat_P_to_subgroup'
421
422 INTEGER :: handle, jkp
423 TYPE(cp_fm_type) :: fm_dummy
424 TYPE(dbcsr_type), POINTER :: mat_p_omega_im, mat_p_omega_re
425
426 CALL timeset(routinen, handle)
427
428 IF (ikp_local == -1) THEN
429
430 mat_p_omega_re => mat_p_omega_kp(1, ikp)%matrix
431 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
432 CALL copy_dbcsr_to_fm(mat_p_omega_re, fm_mat_work)
433 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_q, z_one, fm_mat_work)
434
435 mat_p_omega_im => mat_p_omega_kp(2, ikp)%matrix
436 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
437 CALL copy_dbcsr_to_fm(mat_p_omega_im, fm_mat_work)
438 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_q, gaussi, fm_mat_work)
439
440 ELSE
441
442 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
443
444 DO jkp = 1, nkp
445
446 mat_p_omega_re => mat_p_omega_kp(1, jkp)%matrix
447
448 CALL cp_fm_set_all(fm_mat_ri_global_work, 0.0_dp)
449 CALL copy_dbcsr_to_fm(mat_p_omega_re, fm_mat_ri_global_work)
450
451 CALL para_env%sync()
452
453 IF (ikp_local == jkp) THEN
454 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_mat_work, para_env)
455 ELSE
456 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_dummy, para_env)
457 END IF
458
459 CALL para_env%sync()
460
461 END DO
462
463 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_q, z_one, fm_mat_work)
464
465 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
466
467 DO jkp = 1, nkp
468
469 mat_p_omega_im => mat_p_omega_kp(2, jkp)%matrix
470
471 CALL cp_fm_set_all(fm_mat_ri_global_work, 0.0_dp)
472 CALL copy_dbcsr_to_fm(mat_p_omega_im, fm_mat_ri_global_work)
473
474 CALL para_env%sync()
475
476 IF (ikp_local == jkp) THEN
477 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_mat_work, para_env)
478 ELSE
479 CALL cp_fm_copy_general(fm_mat_ri_global_work, fm_dummy, para_env)
480 END IF
481
482 CALL para_env%sync()
483
484 END DO
485
486 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_q, gaussi, fm_mat_work)
487
488 CALL cp_fm_set_all(fm_mat_work, 0.0_dp)
489
490 END IF
491
492 CALL para_env%sync()
493
494 CALL timestop(handle)
495
496 END SUBROUTINE mat_p_to_subgroup
497
498! **************************************************************************************************
499!> \brief ...
500!> \param cfm_mat_Q ...
501!> \param para_env ...
502!> \param trace_Qomega ...
503!> \param dimen_RI ...
504! **************************************************************************************************
505 SUBROUTINE cholesky_decomp_q(cfm_mat_Q, para_env, trace_Qomega, dimen_RI)
506
507 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
508 TYPE(mp_para_env_type), INTENT(IN) :: para_env
509 REAL(kind=dp), DIMENSION(:), INTENT(OUT) :: trace_qomega
510 INTEGER, INTENT(IN) :: dimen_ri
511
512 CHARACTER(LEN=*), PARAMETER :: routinen = 'cholesky_decomp_Q'
513
514 INTEGER :: handle, i_global, iib, info_chol, &
515 j_global, jjb, ncol_local, nrow_local
516 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
517 TYPE(cp_cfm_type) :: cfm_mat_q_tmp, cfm_mat_work
518
519 CALL timeset(routinen, handle)
520
521 CALL cp_cfm_create(cfm_mat_work, cfm_mat_q%matrix_struct)
522 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
523
524 CALL cp_cfm_create(cfm_mat_q_tmp, cfm_mat_q%matrix_struct)
525 CALL cp_cfm_set_all(cfm_mat_q_tmp, z_zero)
526
527 ! get info of fm_mat_Q
528 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
529 nrow_local=nrow_local, &
530 ncol_local=ncol_local, &
531 row_indices=row_indices, &
532 col_indices=col_indices)
533
534 ! calculate the trace of Q and add 1 on the diagonal
535 trace_qomega = 0.0_dp
536!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
537!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,trace_Qomega,cfm_mat_Q,dimen_RI)
538 DO jjb = 1, ncol_local
539 j_global = col_indices(jjb)
540 DO iib = 1, nrow_local
541 i_global = row_indices(iib)
542 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
543 trace_qomega(i_global) = real(cfm_mat_q%local_data(iib, jjb))
544 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) + z_one
545 END IF
546 END DO
547 END DO
548 CALL para_env%sum(trace_qomega)
549
550 CALL cp_cfm_to_cfm(cfm_mat_q, cfm_mat_q_tmp)
551
552 CALL cp_cfm_cholesky_decompose(matrix=cfm_mat_q, n=dimen_ri, info_out=info_chol)
553
554 cpassert(info_chol == 0)
555
556 CALL cp_cfm_release(cfm_mat_work)
557 CALL cp_cfm_release(cfm_mat_q_tmp)
558
559 CALL timestop(handle)
560
561 END SUBROUTINE cholesky_decomp_q
562
563! **************************************************************************************************
564!> \brief ...
565!> \param Erpa ...
566!> \param cfm_mat_Q ...
567!> \param para_env ...
568!> \param trace_Qomega ...
569!> \param dimen_RI ...
570!> \param freq_weight ...
571!> \param kp_weight ...
572! **************************************************************************************************
573 SUBROUTINE frequency_and_kpoint_integration(Erpa, cfm_mat_Q, para_env, trace_Qomega, &
574 dimen_RI, freq_weight, kp_weight)
575
576 REAL(kind=dp), INTENT(INOUT) :: erpa
577 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
578 TYPE(mp_para_env_type), INTENT(IN) :: para_env
579 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: trace_qomega
580 INTEGER, INTENT(IN) :: dimen_ri
581 REAL(kind=dp), INTENT(IN) :: freq_weight, kp_weight
582
583 CHARACTER(LEN=*), PARAMETER :: routinen = 'frequency_and_kpoint_integration'
584
585 INTEGER :: handle, i_global, iib, j_global, jjb, &
586 ncol_local, nrow_local
587 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
588 REAL(kind=dp) :: fcomega
589 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: q_log
590
591 CALL timeset(routinen, handle)
592
593 ! get info of cholesky_decomposed(fm_mat_Q)
594 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
595 nrow_local=nrow_local, &
596 ncol_local=ncol_local, &
597 row_indices=row_indices, &
598 col_indices=col_indices)
599
600 ALLOCATE (q_log(dimen_ri))
601 q_log = 0.0_dp
602!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
603!$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,Q_log,cfm_mat_Q,dimen_RI)
604 DO jjb = 1, ncol_local
605 j_global = col_indices(jjb)
606 DO iib = 1, nrow_local
607 i_global = row_indices(iib)
608 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
609 q_log(i_global) = 2.0_dp*log(real(cfm_mat_q%local_data(iib, jjb)))
610 END IF
611 END DO
612 END DO
613 CALL para_env%sum(q_log)
614
615 fcomega = 0.0_dp
616 DO iib = 1, dimen_ri
617 IF (modulo(iib, para_env%num_pe) /= para_env%mepos) cycle
618 ! FComega=FComega+(LOG(Q_log(iiB))-trace_Qomega(iiB))/2.0_dp
619 fcomega = fcomega + (q_log(iib) - trace_qomega(iib))/2.0_dp
620 END DO
621
622 erpa = erpa + fcomega*freq_weight*kp_weight
623
624 DEALLOCATE (q_log)
625
626 CALL timestop(handle)
627
628 END SUBROUTINE frequency_and_kpoint_integration
629
630! **************************************************************************************************
631!> \brief ...
632!> \param tj_dummy ...
633!> \param tau_tj_dummy ...
634!> \param weights_cos_tf_w_to_t_dummy ...
635! **************************************************************************************************
636 SUBROUTINE get_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
637
638 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
639 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
640 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
641 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
642
643 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_dummys'
644
645 INTEGER :: handle
646
647 CALL timeset(routinen, handle)
648
649 ALLOCATE (weights_cos_tf_w_to_t_dummy(1, 1))
650 ALLOCATE (tj_dummy(1))
651 ALLOCATE (tau_tj_dummy(1))
652
653 tj_dummy(1) = 0.0_dp
654 tau_tj_dummy(1) = 0.0_dp
655 weights_cos_tf_w_to_t_dummy(1, 1) = 1.0_dp
656
657 CALL timestop(handle)
658
659 END SUBROUTINE
660
661! **************************************************************************************************
662!> \brief ...
663!> \param tj_dummy ...
664!> \param tau_tj_dummy ...
665!> \param weights_cos_tf_w_to_t_dummy ...
666! **************************************************************************************************
667 SUBROUTINE release_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
668
669 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
670 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
671 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
672 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
673
674 CHARACTER(LEN=*), PARAMETER :: routinen = 'release_dummys'
675
676 INTEGER :: handle
677
678 CALL timeset(routinen, handle)
679
680 DEALLOCATE (weights_cos_tf_w_to_t_dummy)
681 DEALLOCATE (tj_dummy)
682 DEALLOCATE (tau_tj_dummy)
683
684 CALL timestop(handle)
685
686 END SUBROUTINE
687
688! **************************************************************************************************
689!> \brief ...
690!> \param mat_P_omega ...
691!> \param qs_env ...
692!> \param kpoints ...
693!> \param jquad ...
694!> \param unit_nr ...
695! **************************************************************************************************
696 SUBROUTINE get_mat_cell_t_from_mat_gamma(mat_P_omega, qs_env, kpoints, jquad, unit_nr)
697 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN) :: mat_p_omega
698 TYPE(qs_environment_type), POINTER :: qs_env
699 TYPE(kpoint_type), POINTER :: kpoints
700 INTEGER, INTENT(IN) :: jquad, unit_nr
701
702 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_mat_cell_T_from_mat_gamma'
703
704 INTEGER :: col, handle, i_cell, i_dim, j_cell, &
705 num_cells_p, num_integ_points, row
706 INTEGER, DIMENSION(3) :: cell_grid_p, periodic
707 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell_p
708 LOGICAL :: i_cell_is_the_minimum_image_cell
709 REAL(kind=dp) :: abs_rab_cell_i, abs_rab_cell_j
710 REAL(kind=dp), DIMENSION(3) :: cell_vector, cell_vector_j, rab_cell_i, &
711 rab_cell_j
712 REAL(kind=dp), DIMENSION(3, 3) :: hmat
713 REAL(kind=dp), DIMENSION(:, :), POINTER :: data_block
714 TYPE(cell_type), POINTER :: cell
715 TYPE(dbcsr_iterator_type) :: iter
716 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
717
718 CALL timeset(routinen, handle)
719
720 NULLIFY (cell, particle_set)
721 CALL get_qs_env(qs_env, cell=cell, &
722 particle_set=particle_set)
723 CALL get_cell(cell=cell, h=hmat, periodic=periodic)
724
725 DO i_dim = 1, 3
726 ! we have at most 3 neigboring cells per dimension and at least one because
727 ! the density response at Gamma is only divided to neighboring
728 IF (periodic(i_dim) == 1) THEN
729 cell_grid_p(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
730 ELSE
731 cell_grid_p(i_dim) = 1
732 END IF
733 END DO
734
735 ! overwrite the cell indices in kpoints
736 CALL init_cell_index_rpa(cell_grid_p, kpoints%cell_to_index, kpoints%index_to_cell, cell)
737
738 index_to_cell_p => kpoints%index_to_cell
739
740 num_cells_p = SIZE(index_to_cell_p, 2)
741
742 num_integ_points = SIZE(mat_p_omega, 1)
743
744 ! first, copy the Gamma-only result from mat_P_omega(1) into all other matrices and
745 ! remove the blocks later which do not belong to the cell index
746 DO i_cell = 2, num_cells_p
747 CALL dbcsr_copy(mat_p_omega(i_cell)%matrix, &
748 mat_p_omega(1)%matrix)
749 END DO
750
751 IF (jquad == 1 .AND. unit_nr > 0) THEN
752 WRITE (unit_nr, '(T3,A,T66,ES15.2)') 'GW_INFO| RI regularization parameter: ', &
753 qs_env%mp2_env%ri_rpa_im_time%regularization_RI
754 WRITE (unit_nr, '(T3,A,T66,ES15.2)') 'GW_INFO| eps_eigval_S: ', &
755 qs_env%mp2_env%ri_rpa_im_time%eps_eigval_S
756 IF (qs_env%mp2_env%ri_rpa_im_time%make_chi_pos_definite) THEN
757 WRITE (unit_nr, '(T3,A,T81)') &
758 'GW_INFO| Make chi(iw,k) positive definite? TRUE'
759 ELSE
760 WRITE (unit_nr, '(T3,A,T81)') &
761 'GW_INFO| Make chi(iw,k) positive definite? FALSE'
762 END IF
763
764 END IF
765
766 DO i_cell = 1, num_cells_p
767
768 CALL dbcsr_iterator_start(iter, mat_p_omega(i_cell)%matrix)
769 DO WHILE (dbcsr_iterator_blocks_left(iter))
770 CALL dbcsr_iterator_next_block(iter, row, col, data_block)
771
772 cell_vector(1:3) = matmul(hmat, real(index_to_cell_p(1:3, i_cell), dp))
773 rab_cell_i(1:3) = pbc(particle_set(row)%r(1:3), cell) - &
774 (pbc(particle_set(col)%r(1:3), cell) + cell_vector(1:3))
775 abs_rab_cell_i = sqrt(rab_cell_i(1)**2 + rab_cell_i(2)**2 + rab_cell_i(3)**2)
776
777 ! minimum image convention
778 i_cell_is_the_minimum_image_cell = .true.
779 DO j_cell = 1, num_cells_p
780 cell_vector_j(1:3) = matmul(hmat, real(index_to_cell_p(1:3, j_cell), dp))
781 rab_cell_j(1:3) = pbc(particle_set(row)%r(1:3), cell) - &
782 (pbc(particle_set(col)%r(1:3), cell) + cell_vector_j(1:3))
783 abs_rab_cell_j = sqrt(rab_cell_j(1)**2 + rab_cell_j(2)**2 + rab_cell_j(3)**2)
784
785 IF (abs_rab_cell_i > abs_rab_cell_j + 1.0e-6_dp) THEN
786 i_cell_is_the_minimum_image_cell = .false.
787 END IF
788 END DO
789
790 IF (.NOT. i_cell_is_the_minimum_image_cell) THEN
791 data_block(:, :) = data_block(:, :)*0.0_dp
792 END IF
793
794 END DO
795 CALL dbcsr_iterator_stop(iter)
796
797 END DO
798
799 CALL timestop(handle)
800
801 END SUBROUTINE get_mat_cell_t_from_mat_gamma
802
803! **************************************************************************************************
804!> \brief ...
805!> \param mat_P_omega ...
806!> \param mat_P_omega_kp ...
807!> \param kpoints ...
808!> \param eps_filter_im_time ...
809!> \param jquad ...
810! **************************************************************************************************
811 SUBROUTINE transform_p_from_real_space_to_kpoints(mat_P_omega, mat_P_omega_kp, &
812 kpoints, eps_filter_im_time, jquad)
813
814 TYPE(dbcsr_p_type), DIMENSION(:, :), INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
815 TYPE(kpoint_type), POINTER :: kpoints
816 REAL(kind=dp), INTENT(IN) :: eps_filter_im_time
817 INTEGER, INTENT(IN) :: jquad
818
819 CHARACTER(LEN=*), PARAMETER :: routinen = 'transform_P_from_real_space_to_kpoints'
820
821 INTEGER :: handle, icell, nkp, num_integ_points
822
823 CALL timeset(routinen, handle)
824
825 num_integ_points = SIZE(mat_p_omega, 1)
826 nkp = SIZE(mat_p_omega, 2)
827
828 CALL real_space_to_kpoint_transform_rpa(mat_p_omega_kp(1, :), mat_p_omega_kp(2, :), mat_p_omega(jquad, :), &
829 kpoints, eps_filter_im_time)
830
831 DO icell = 1, SIZE(mat_p_omega, 2)
832 CALL dbcsr_set(mat_p_omega(jquad, icell)%matrix, 0.0_dp)
833 CALL dbcsr_filter(mat_p_omega(jquad, icell)%matrix, 1.0_dp)
834 END DO
835
836 CALL timestop(handle)
837
838 END SUBROUTINE transform_p_from_real_space_to_kpoints
839
840! **************************************************************************************************
841!> \brief ...
842!> \param real_mat_kp ...
843!> \param imag_mat_kp ...
844!> \param mat_real_space ...
845!> \param kpoints ...
846!> \param eps_filter_im_time ...
847!> \param real_mat_real_space ...
848! **************************************************************************************************
849 SUBROUTINE real_space_to_kpoint_transform_rpa(real_mat_kp, imag_mat_kp, mat_real_space, &
850 kpoints, eps_filter_im_time, real_mat_real_space)
851
852 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(INOUT) :: real_mat_kp, imag_mat_kp, mat_real_space
853 TYPE(kpoint_type), POINTER :: kpoints
854 REAL(kind=dp), INTENT(IN) :: eps_filter_im_time
855 LOGICAL, INTENT(IN), OPTIONAL :: real_mat_real_space
856
857 CHARACTER(LEN=*), PARAMETER :: routinen = 'real_space_to_kpoint_transform_rpa'
858
859 INTEGER :: handle, i_cell, ik, nkp, num_cells
860 INTEGER, DIMENSION(3) :: cell
861 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell
862 LOGICAL :: my_real_mat_real_space
863 REAL(kind=dp) :: arg, coskl, sinkl
864 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
865 TYPE(dbcsr_type) :: mat_work
866
867 CALL timeset(routinen, handle)
868
869 my_real_mat_real_space = .true.
870 IF (PRESENT(real_mat_real_space)) my_real_mat_real_space = real_mat_real_space
871
872 CALL dbcsr_create(matrix=mat_work, &
873 template=real_mat_kp(1)%matrix, &
874 matrix_type=dbcsr_type_no_symmetry)
875 CALL dbcsr_reserve_all_blocks(mat_work)
876 CALL dbcsr_set(mat_work, 0.0_dp)
877
878 ! this kpoint environme t should be the kpoints for D(it) and X(it) created in init_cell_index_rpa
879 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp)
880
881 NULLIFY (index_to_cell)
882 index_to_cell => kpoints%index_to_cell
883
884 num_cells = SIZE(index_to_cell, 2)
885
886 cpassert(SIZE(mat_real_space) >= num_cells/2 + 1)
887
888 DO ik = 1, nkp
889
890 CALL dbcsr_reserve_all_blocks(real_mat_kp(ik)%matrix)
891 CALL dbcsr_reserve_all_blocks(imag_mat_kp(ik)%matrix)
892
893 CALL dbcsr_set(real_mat_kp(ik)%matrix, 0.0_dp)
894 CALL dbcsr_set(imag_mat_kp(ik)%matrix, 0.0_dp)
895
896 DO i_cell = 1, num_cells/2 + 1
897
898 cell(:) = index_to_cell(:, i_cell)
899
900 arg = real(cell(1), dp)*xkp(1, ik) + real(cell(2), dp)*xkp(2, ik) + real(cell(3), dp)*xkp(3, ik)
901 coskl = cos(twopi*arg)
902 sinkl = sin(twopi*arg)
903
904 IF (my_real_mat_real_space) THEN
905 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
906 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, sinkl)
907 ELSE
908 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, -sinkl)
909 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
910 END IF
911
912 IF (.NOT. (cell(1) == 0 .AND. cell(2) == 0 .AND. cell(3) == 0)) THEN
913
914 CALL dbcsr_transposed(mat_work, mat_real_space(i_cell)%matrix)
915
916 IF (my_real_mat_real_space) THEN
917 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, coskl)
918 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
919 ELSE
920 ! for an imaginary real-space matrix, we need to consider the imaginary unit
921 ! and we need to take into account that the transposed gives an extra "-" sign
922 ! because the transposed is actually Hermitian conjugate
923 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
924 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -coskl)
925 END IF
926
927 CALL dbcsr_set(mat_work, 0.0_dp)
928
929 END IF
930
931 END DO
932
933 CALL dbcsr_filter(real_mat_kp(ik)%matrix, eps_filter_im_time)
934 CALL dbcsr_filter(imag_mat_kp(ik)%matrix, eps_filter_im_time)
935
936 END DO
937
938 CALL dbcsr_release(mat_work)
939
940 CALL timestop(handle)
941
943
944! **************************************************************************************************
945!> \brief ...
946!> \param mat_a ...
947!> \param mat_b ...
948!> \param alpha ...
949!> \param beta ...
950! **************************************************************************************************
951 SUBROUTINE dbcsr_add_local(mat_a, mat_b, alpha, beta)
952 TYPE(dbcsr_type), INTENT(INOUT) :: mat_a, mat_b
953 REAL(kind=dp), INTENT(IN) :: alpha, beta
954
955 INTEGER :: col, row
956 LOGICAL :: found
957 REAL(kind=dp), DIMENSION(:, :), POINTER :: block_to_compute, data_block
958 TYPE(dbcsr_iterator_type) :: iter
959
960 CALL dbcsr_iterator_start(iter, mat_b)
961 DO WHILE (dbcsr_iterator_blocks_left(iter))
962 CALL dbcsr_iterator_next_block(iter, row, col, data_block)
963
964 NULLIFY (block_to_compute)
965 CALL dbcsr_get_block_p(matrix=mat_a, &
966 row=row, col=col, block=block_to_compute, found=found)
967
968 cpassert(found)
969
970 block_to_compute(:, :) = alpha*block_to_compute(:, :) + beta*data_block(:, :)
971
972 END DO
973 CALL dbcsr_iterator_stop(iter)
974
975 END SUBROUTINE dbcsr_add_local
976
977! **************************************************************************************************
978!> \brief ...
979!> \param fm_mat_W_tau ...
980!> \param cfm_mat_Q ...
981!> \param fm_mat_L_re ...
982!> \param fm_mat_L_im ...
983!> \param dimen_RI ...
984!> \param num_integ_points ...
985!> \param jquad ...
986!> \param ikp ...
987!> \param tj ...
988!> \param tau_tj ...
989!> \param weights_cos_tf_w_to_t ...
990!> \param ikp_local ...
991!> \param para_env ...
992!> \param kpoints ...
993!> \param qs_env ...
994!> \param wkp_W ...
995! **************************************************************************************************
996 SUBROUTINE compute_wc_real_space_tau_gw(fm_mat_W_tau, cfm_mat_Q, fm_mat_L_re, fm_mat_L_im, &
997 dimen_RI, num_integ_points, jquad, &
998 ikp, tj, tau_tj, weights_cos_tf_w_to_t, ikp_local, &
999 para_env, kpoints, qs_env, wkp_W)
1000
1001 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_w_tau
1002 TYPE(cp_cfm_type), INTENT(IN) :: cfm_mat_q
1003 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_l_re, fm_mat_l_im
1004 INTEGER, INTENT(IN) :: dimen_ri, num_integ_points, jquad, ikp
1005 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: tj
1006 REAL(kind=dp), DIMENSION(num_integ_points), &
1007 INTENT(IN) :: tau_tj
1008 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: weights_cos_tf_w_to_t
1009 INTEGER, INTENT(IN) :: ikp_local
1010 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1011 TYPE(kpoint_type), INTENT(IN), POINTER :: kpoints
1012 TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
1013 REAL(kind=dp), DIMENSION(:), INTENT(IN) :: wkp_w
1014
1015 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_Wc_real_space_tau_GW'
1016
1017 INTEGER :: handle, handle2, i_global, iatom, iatom_old, iib, iquad, irow, j_global, jatom, &
1018 jatom_old, jcol, jjb, jkp, ncol_local, nkp, nrow_local, num_cells
1019 INTEGER, ALLOCATABLE, DIMENSION(:) :: atom_from_ri_index
1020 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
1021 INTEGER, DIMENSION(:, :), POINTER :: index_to_cell
1022 REAL(kind=dp) :: contribution, omega, tau, weight, &
1023 weight_im, weight_re
1024 REAL(kind=dp), DIMENSION(3, 3) :: hmat
1025 REAL(kind=dp), DIMENSION(:), POINTER :: wkp
1026 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1027 TYPE(cell_type), POINTER :: cell
1028 TYPE(cp_cfm_type) :: cfm_mat_l, cfm_mat_work, cfm_mat_work_2
1029 TYPE(cp_fm_type) :: fm_dummy, fm_mat_work_global, &
1030 fm_mat_work_local
1031 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1032
1033 CALL timeset(routinen, handle)
1034
1035 CALL timeset(routinen//"_1", handle2)
1036
1037 CALL cp_cfm_create(cfm_mat_work, cfm_mat_q%matrix_struct)
1038 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
1039
1040 CALL cp_cfm_create(cfm_mat_work_2, cfm_mat_q%matrix_struct)
1041 CALL cp_cfm_set_all(cfm_mat_work_2, z_zero)
1042
1043 CALL cp_cfm_create(cfm_mat_l, cfm_mat_q%matrix_struct)
1044 CALL cp_cfm_set_all(cfm_mat_l, z_zero)
1045
1046 ! Copy fm_mat_L_re and fm_mat_L_re to cfm_mat_L
1047 CALL cp_cfm_scale_and_add_fm(z_zero, cfm_mat_l, z_one, fm_mat_l_re)
1048 CALL cp_cfm_scale_and_add_fm(z_one, cfm_mat_l, gaussi, fm_mat_l_im)
1049
1050 CALL cp_fm_create(fm_mat_work_global, fm_mat_w_tau(1)%matrix_struct)
1051 CALL cp_fm_set_all(fm_mat_work_global, 0.0_dp)
1052
1053 CALL cp_fm_create(fm_mat_work_local, cfm_mat_q%matrix_struct)
1054 CALL cp_fm_set_all(fm_mat_work_local, 0.0_dp)
1055
1056 CALL timestop(handle2)
1057
1058 CALL timeset(routinen//"_2", handle2)
1059
1060 ! calculate [1+Q(iw')]^-1
1061 CALL cp_cfm_cholesky_invert(cfm_mat_q)
1062
1063 ! symmetrize the result
1064 CALL cp_cfm_uplo_to_full(cfm_mat_q)
1065
1066 ! subtract exchange part by subtracing identity matrix from epsilon
1067 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
1068 nrow_local=nrow_local, &
1069 ncol_local=ncol_local, &
1070 row_indices=row_indices, &
1071 col_indices=col_indices)
1072
1073 DO jjb = 1, ncol_local
1074 j_global = col_indices(jjb)
1075 DO iib = 1, nrow_local
1076 i_global = row_indices(iib)
1077 IF (j_global == i_global .AND. i_global <= dimen_ri) THEN
1078 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) - z_one
1079 END IF
1080 END DO
1081 END DO
1082
1083 CALL timestop(handle2)
1084
1085 CALL timeset(routinen//"_3", handle2)
1086
1087 ! work = epsilon(iw,k)*V^1/2(k)
1088 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_q, cfm_mat_l, &
1089 z_zero, cfm_mat_work)
1090
1091 ! W(iw,k) = V^1/2(k)*work
1092 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, z_one, cfm_mat_l, cfm_mat_work, &
1093 z_zero, cfm_mat_work_2)
1094
1095 CALL timestop(handle2)
1096
1097 CALL timeset(routinen//"_4", handle2)
1098
1099 CALL get_kpoint_info(kpoints, xkp=xkp, wkp=wkp, nkp=nkp)
1100 index_to_cell => kpoints%index_to_cell
1101 num_cells = SIZE(index_to_cell, 2)
1102
1103 CALL cp_cfm_set_all(cfm_mat_work, z_zero)
1104
1105 ALLOCATE (atom_from_ri_index(dimen_ri))
1106
1107 CALL get_atom_index_from_basis_function_index(qs_env, atom_from_ri_index, dimen_ri, "RI_AUX")
1108
1109 NULLIFY (cell, particle_set)
1110 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1111 CALL get_cell(cell=cell, h=hmat)
1112 iatom_old = 0
1113 jatom_old = 0
1114
1115 CALL cp_cfm_get_info(matrix=cfm_mat_q, &
1116 nrow_local=nrow_local, &
1117 ncol_local=ncol_local, &
1118 row_indices=row_indices, &
1119 col_indices=col_indices)
1120
1121 DO irow = 1, nrow_local
1122 DO jcol = 1, ncol_local
1123
1124 iatom = atom_from_ri_index(row_indices(irow))
1125 jatom = atom_from_ri_index(col_indices(jcol))
1126
1127 IF (iatom .NE. iatom_old .OR. jatom .NE. jatom_old) THEN
1128
1129 ! symmetrize=.FALSE. necessary since we already have a symmetrized index_to_cell
1130 CALL compute_weight_re_im(weight_re, weight_im, &
1131 num_cells, iatom, jatom, xkp(1:3, ikp), wkp_w(ikp), &
1132 cell, index_to_cell, hmat, particle_set)
1133
1134 iatom_old = iatom
1135 jatom_old = jatom
1136
1137 END IF
1138
1139 contribution = weight_re*real(cfm_mat_work_2%local_data(irow, jcol)) + &
1140 weight_im*aimag(cfm_mat_work_2%local_data(irow, jcol))
1141
1142 fm_mat_work_local%local_data(irow, jcol) = fm_mat_work_local%local_data(irow, jcol) + contribution
1143
1144 END DO
1145 END DO
1146
1147 CALL timestop(handle2)
1148
1149 CALL timeset(routinen//"_5", handle2)
1150
1151 IF (ikp_local == -1) THEN
1152
1153 CALL cp_fm_copy_general(fm_mat_work_local, fm_mat_work_global, para_env)
1154
1155 DO iquad = 1, num_integ_points
1156
1157 omega = tj(jquad)
1158 tau = tau_tj(iquad)
1159 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1160
1161 IF (jquad == 1 .AND. ikp == 1) THEN
1162 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1163 END IF
1164
1165 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)
1166
1167 END DO
1168
1169 ELSE
1170
1171 DO jkp = 1, nkp
1172
1173 CALL para_env%sync()
1174
1175 IF (ikp_local == jkp) THEN
1176 CALL cp_fm_copy_general(fm_mat_work_local, fm_mat_work_global, para_env)
1177 ELSE
1178 CALL cp_fm_copy_general(fm_dummy, fm_mat_work_global, para_env)
1179 END IF
1180
1181 CALL para_env%sync()
1182
1183 DO iquad = 1, num_integ_points
1184
1185 omega = tj(jquad)
1186 tau = tau_tj(iquad)
1187 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1188
1189 IF (jquad == 1 .AND. jkp == 1) THEN
1190 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1191 END IF
1192
1193 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_w_tau(iquad), beta=weight, &
1194 matrix_b=fm_mat_work_global)
1195
1196 END DO
1197
1198 END DO
1199
1200 END IF
1201
1202 CALL cp_cfm_release(cfm_mat_work)
1203 CALL cp_cfm_release(cfm_mat_work_2)
1204 CALL cp_cfm_release(cfm_mat_l)
1205 CALL cp_fm_release(fm_mat_work_global)
1206 CALL cp_fm_release(fm_mat_work_local)
1207
1208 DEALLOCATE (atom_from_ri_index)
1209
1210 CALL timestop(handle2)
1211
1212 CALL timestop(handle)
1213
1214 END SUBROUTINE compute_wc_real_space_tau_gw
1215
1216! **************************************************************************************************
1217!> \brief ...
1218!> \param fm_mat_W ...
1219!> \param fm_matrix_Minv ...
1220!> \param para_env ...
1221!> \param dimen_RI ...
1222!> \param num_integ_points ...
1223! **************************************************************************************************
1224 SUBROUTINE wc_to_minv_wc_minv(fm_mat_W, fm_matrix_Minv, para_env, dimen_RI, num_integ_points)
1225 TYPE(cp_fm_type), DIMENSION(:) :: fm_mat_w
1226 TYPE(cp_fm_type), DIMENSION(:, :) :: fm_matrix_minv
1227 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
1228 INTEGER :: dimen_ri, num_integ_points
1229
1230 CHARACTER(LEN=*), PARAMETER :: routinen = 'Wc_to_Minv_Wc_Minv'
1231
1232 INTEGER :: handle, jquad
1233 TYPE(cp_fm_type) :: fm_work_minv, fm_work_minv_w
1234
1235 CALL timeset(routinen, handle)
1236
1237 CALL cp_fm_create(fm_work_minv, fm_mat_w(1)%matrix_struct)
1238 CALL cp_fm_copy_general(fm_matrix_minv(1, 1), fm_work_minv, para_env)
1239
1240 CALL cp_fm_create(fm_work_minv_w, fm_mat_w(1)%matrix_struct)
1241
1242 DO jquad = 1, num_integ_points
1243
1244 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv, fm_mat_w(jquad), &
1245 0.0_dp, fm_work_minv_w)
1246 CALL parallel_gemm('N', 'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv_w, fm_work_minv, &
1247 0.0_dp, fm_mat_w(jquad))
1248
1249 END DO
1250
1251 CALL cp_fm_release(fm_work_minv)
1252
1253 CALL cp_fm_release(fm_work_minv_w)
1254
1255 CALL timestop(handle)
1256
1257 END SUBROUTINE wc_to_minv_wc_minv
1258
1259! **************************************************************************************************
1260!> \brief ...
1261!> \param qs_env ...
1262!> \param wkp_W ...
1263!> \param wkp_V ...
1264!> \param kpoints ...
1265!> \param h_inv ...
1266!> \param periodic ...
1267! **************************************************************************************************
1268 SUBROUTINE compute_wkp_w(qs_env, wkp_W, wkp_V, kpoints, h_inv, periodic)
1269
1270 TYPE(qs_environment_type), POINTER :: qs_env
1271 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), &
1272 INTENT(OUT) :: wkp_w, wkp_v
1273 TYPE(kpoint_type), POINTER :: kpoints
1274 REAL(kind=dp), DIMENSION(3, 3) :: h_inv
1275 INTEGER, DIMENSION(3) :: periodic
1276
1277 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_wkp_W'
1278
1279 INTEGER :: handle, i_x, ikp, info, j_y, k_z, &
1280 kpoint_weights_w_method, n_x, n_y, &
1281 n_z, nkp, nsuperfine, num_lin_eqs
1282 REAL(kind=dp) :: exp_kpoints, integral, k_sq, weight
1283 REAL(kind=dp), DIMENSION(3) :: k_vec, x_vec
1284 REAL(kind=dp), DIMENSION(:), POINTER :: right_side, wkp, wkp_tmp
1285 REAL(kind=dp), DIMENSION(:, :), POINTER :: matrix_lin_eqs, xkp
1286
1287 CALL timeset(routinen, handle)
1288
1289 kpoint_weights_w_method = qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method
1290
1291 CALL get_kpoint_info(kpoints, xkp=xkp, wkp=wkp, nkp=nkp)
1292
1293 ! we determine the kpoint weights of the Monkhors Pack mesh new
1294 ! such that the functions 1/k^2, 1/k and const are integrated exactly
1295 ! in the Brillouin zone
1296 ! this is done by minimizing sum_i |w_i|^2 where w_i are the weights of
1297 ! the i-th kpoint under the following constraints:
1298 ! 1) 1/k^2, 1/k and const are integrated exactly
1299 ! 2) the kpoint weights of kpoints with identical absolute value are
1300 ! the same, of e.g. (1/8,3/8,3/8) same weight as for (3/8,1/8,3/8)
1301 ! for 1d and 2d materials: we use ordinary Monkhorst-Pack weights, checked
1302 ! by SUM(periodic) == 3
1303 ALLOCATE (wkp_v(nkp), wkp_w(nkp))
1304
1305 ! for exchange part of self-energy, we use truncated Coulomb operator that should be fine
1306 ! with uniform weights (without k-point extrapolation)
1307 IF (ALLOCATED(qs_env%mp2_env%ri_rpa_im_time%wkp_V)) THEN
1308 wkp_v(:) = qs_env%mp2_env%ri_rpa_im_time%wkp_V(:)
1309 ELSE
1310 wkp_v(:) = wkp(:)
1311 END IF
1312
1313 IF (kpoint_weights_w_method == kp_weights_w_uniform) THEN
1314
1315 ! in the k-point weights wkp, there might be k-point extrapolation included
1316 wkp_w(:) = wkp(:)
1317
1318 ELSE IF (kpoint_weights_w_method == kp_weights_w_tailored .OR. &
1319 kpoint_weights_w_method == kp_weights_w_auto) THEN
1320
1321 IF (kpoint_weights_w_method == kp_weights_w_tailored) &
1322 exp_kpoints = qs_env%mp2_env%ri_rpa_im_time%exp_tailored_weights
1323
1324 IF (kpoint_weights_w_method == kp_weights_w_auto) THEN
1325 IF (sum(periodic) == 2) exp_kpoints = -1.0_dp
1326 END IF
1327
1328 ! first, compute the integral of f(k)=1/k^2 and 1/k on super fine grid
1329 nsuperfine = 500
1330 integral = 0.0_dp
1331
1332 IF (periodic(1) == 1) THEN
1333 n_x = nsuperfine
1334 ELSE
1335 n_x = 1
1336 END IF
1337 IF (periodic(2) == 1) THEN
1338 n_y = nsuperfine
1339 ELSE
1340 n_y = 1
1341 END IF
1342 IF (periodic(3) == 1) THEN
1343 n_z = nsuperfine
1344 ELSE
1345 n_z = 1
1346 END IF
1347
1348 ! actually, there is the factor *det_3x3(h_inv) missing to account for the
1349 ! integration volume but for wkp det_3x3(h_inv) is needed
1350 weight = 1.0_dp/(real(n_x, dp)*real(n_y, dp)*real(n_z, dp))
1351 DO i_x = 1, n_x
1352 DO j_y = 1, n_y
1353 DO k_z = 1, n_z
1354
1355 IF (periodic(1) == 1) THEN
1356 x_vec(1) = (real(i_x - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1357 ELSE
1358 x_vec(1) = 0.0_dp
1359 END IF
1360 IF (periodic(2) == 1) THEN
1361 x_vec(2) = (real(j_y - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1362 ELSE
1363 x_vec(2) = 0.0_dp
1364 END IF
1365 IF (periodic(3) == 1) THEN
1366 x_vec(3) = (real(k_z - nsuperfine/2, dp) - 0.5_dp)/real(nsuperfine, dp)
1367 ELSE
1368 x_vec(3) = 0.0_dp
1369 END IF
1370
1371 k_vec = matmul(h_inv(1:3, 1:3), x_vec)
1372 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1373 integral = integral + weight*k_sq**(exp_kpoints*0.5_dp)
1374
1375 END DO
1376 END DO
1377 END DO
1378
1379 num_lin_eqs = nkp + 2
1380
1381 ALLOCATE (matrix_lin_eqs(num_lin_eqs, num_lin_eqs))
1382 matrix_lin_eqs(:, :) = 0.0_dp
1383
1384 DO ikp = 1, nkp
1385
1386 k_vec = matmul(h_inv(1:3, 1:3), xkp(1:3, ikp))
1387 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1388
1389 matrix_lin_eqs(ikp, ikp) = 2.0_dp
1390 matrix_lin_eqs(ikp, nkp + 1) = 1.0_dp
1391 matrix_lin_eqs(nkp + 1, ikp) = 1.0_dp
1392
1393 matrix_lin_eqs(ikp, nkp + 2) = k_sq**(exp_kpoints*0.5_dp)
1394 matrix_lin_eqs(nkp + 2, ikp) = k_sq**(exp_kpoints*0.5_dp)
1395
1396 END DO
1397
1398 CALL invmat(matrix_lin_eqs, info)
1399 ! check whether inversion was successful
1400 cpassert(info == 0)
1401
1402 ALLOCATE (right_side(num_lin_eqs))
1403 right_side = 0.0_dp
1404 right_side(nkp + 1) = 1.0_dp
1405 ! divide integral by two because CP2K k-mesh already considers symmetry k <-> -k
1406 right_side(nkp + 2) = integral
1407
1408 ALLOCATE (wkp_tmp(num_lin_eqs))
1409
1410 wkp_tmp(1:num_lin_eqs) = matmul(matrix_lin_eqs, right_side)
1411
1412 wkp_w(1:nkp) = wkp_tmp(1:nkp)
1413
1414 DEALLOCATE (matrix_lin_eqs, right_side, wkp_tmp)
1415
1416 END IF
1417
1418 CALL timestop(handle)
1419
1420 END SUBROUTINE
1421
1422! **************************************************************************************************
1423!> \brief ...
1424!> \param qs_env ...
1425!> \param Eigenval_kp ...
1426! **************************************************************************************************
1427 SUBROUTINE get_bandstruc_and_k_dependent_mos(qs_env, Eigenval_kp)
1428 TYPE(qs_environment_type), POINTER :: qs_env
1429 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: eigenval_kp
1430
1431 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_bandstruc_and_k_dependent_MOs'
1432
1433 INTEGER :: handle, ikp, ispin, nmo, nspins
1434 INTEGER, DIMENSION(3) :: nkp_grid_g
1435 REAL(kind=dp), DIMENSION(:), POINTER :: ev
1436 REAL(kind=dp), DIMENSION(:, :), POINTER :: kpgeneral
1437 TYPE(kpoint_type), POINTER :: kpoints_sigma
1438 TYPE(mp_para_env_type), POINTER :: para_env
1439
1440 CALL timeset(routinen, handle)
1441
1442 NULLIFY (qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1443 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1444 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1445 para_env)
1446
1447 nkp_grid_g(1:3) = (/1, 1, 1/)
1448
1449 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
1450
1451 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1452 "MONKHORST-PACK", para_env%num_pe, &
1453 mp_grid=nkp_grid_g(1:3))
1454
1455 IF (qs_env%mp2_env%ri_g0w0%do_kpoints_Sigma) THEN
1456
1457 ! set up k-points for GW band structure calculation, will be completed later
1458 CALL get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1459
1460 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1461 "GENERAL", para_env%num_pe, &
1462 kpgeneral=kpgeneral)
1463
1464 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1465 "GENERAL", para_env%num_pe, &
1466 kpgeneral=kpgeneral, with_xc_terms=.false.)
1467
1468 kpoints_sigma => qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma
1469 nmo = SIZE(eigenval_kp, 1)
1470 nspins = SIZE(eigenval_kp, 3)
1471
1472 ALLOCATE (qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(nmo))
1473 qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(:) = eigenval_kp(:, 1, 1)
1474
1475 DEALLOCATE (eigenval_kp)
1476
1477 ALLOCATE (eigenval_kp(nmo, kpoints_sigma%nkp, nspins))
1478
1479 DO ikp = 1, kpoints_sigma%nkp
1480
1481 DO ispin = 1, nspins
1482
1483 ev => kpoints_sigma%kp_env(ikp)%kpoint_env%mos(1, ispin)%eigenvalues
1484
1485 eigenval_kp(:, ikp, ispin) = ev(:)
1486
1487 END DO
1488
1489 END DO
1490
1491 DEALLOCATE (kpgeneral)
1492
1493 END IF
1494
1495 CALL release_hfx_stuff(qs_env)
1496
1497 CALL timestop(handle)
1498
1500
1501! **************************************************************************************************
1502!> \brief releases part of the given qs_env in order to save memory
1503!> \param qs_env the object to release
1504! **************************************************************************************************
1505 SUBROUTINE release_hfx_stuff(qs_env)
1506 TYPE(qs_environment_type), POINTER :: qs_env
1507
1508 IF (ASSOCIATED(qs_env%x_data) .AND. .NOT. qs_env%mp2_env%ri_g0w0%do_ri_Sigma_x) THEN
1509 CALL hfx_release(qs_env%x_data)
1510 END IF
1511
1512 END SUBROUTINE release_hfx_stuff
1513
1514! **************************************************************************************************
1515!> \brief ...
1516!> \param qs_env ...
1517!> \param kpoints ...
1518!> \param scheme ...
1519!> \param group_size_ext ...
1520!> \param mp_grid ...
1521!> \param kpgeneral ...
1522!> \param with_xc_terms ...
1523! **************************************************************************************************
1524 SUBROUTINE create_kp_and_calc_kp_orbitals(qs_env, kpoints, scheme, &
1525 group_size_ext, mp_grid, kpgeneral, with_xc_terms)
1526
1527 TYPE(qs_environment_type), POINTER :: qs_env
1528 TYPE(kpoint_type), POINTER :: kpoints
1529 CHARACTER(LEN=*), INTENT(IN) :: scheme
1530 INTEGER :: group_size_ext
1531 INTEGER, DIMENSION(3), INTENT(IN), OPTIONAL :: mp_grid
1532 REAL(kind=dp), DIMENSION(:, :), INTENT(IN), &
1533 OPTIONAL :: kpgeneral
1534 LOGICAL, OPTIONAL :: with_xc_terms
1535
1536 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_kp_and_calc_kp_orbitals'
1537 COMPLEX(KIND=dp), PARAMETER :: cone = cmplx(1.0_dp, 0.0_dp, kind=dp), &
1538 czero = cmplx(0.0_dp, 0.0_dp, kind=dp), ione = cmplx(0.0_dp, 1.0_dp, kind=dp)
1539
1540 INTEGER :: handle, i_dim, i_re_im, ikp, ispin, nkp, &
1541 nspins
1542 INTEGER, DIMENSION(3) :: cell_grid, periodic
1543 LOGICAL :: my_with_xc_terms
1544 REAL(kind=dp), DIMENSION(:), POINTER :: eigenvalues
1545 TYPE(cell_type), POINTER :: cell
1546 TYPE(cp_blacs_env_type), POINTER :: blacs_env
1547 TYPE(cp_cfm_type) :: cksmat, cmos, csmat, cwork
1548 TYPE(cp_fm_struct_type), POINTER :: matrix_struct
1549 TYPE(cp_fm_type) :: fm_work
1550 TYPE(cp_fm_type), POINTER :: imos, rmos
1551 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, matrix_s_desymm
1552 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_ks_kp, mat_s_kp
1553 TYPE(dft_control_type), POINTER :: dft_control
1554 TYPE(kpoint_env_type), POINTER :: kp
1555 TYPE(mp_para_env_type), POINTER :: para_env
1556 TYPE(qs_scf_env_type), POINTER :: scf_env
1557 TYPE(scf_control_type), POINTER :: scf_control
1558
1559 CALL timeset(routinen, handle)
1560
1561 my_with_xc_terms = .true.
1562 IF (PRESENT(with_xc_terms)) my_with_xc_terms = with_xc_terms
1563
1564 CALL get_qs_env(qs_env, &
1565 para_env=para_env, &
1566 blacs_env=blacs_env, &
1567 matrix_s=matrix_s, &
1568 scf_env=scf_env, &
1569 scf_control=scf_control, &
1570 cell=cell)
1571
1572 ! get kpoints
1573 CALL calculate_kpoints_for_bs(kpoints, scheme, kpgeneral=kpgeneral, mp_grid=mp_grid, &
1574 group_size_ext=group_size_ext)
1575
1576 CALL kpoint_env_initialize(kpoints, para_env, blacs_env)
1577
1578 ! calculate all MOs that are accessible in the given
1579 ! Gaussian AO basis, therefore nadd=1E10
1580 CALL kpoint_initialize_mos(kpoints, qs_env%mos, 2000000000)
1581 CALL kpoint_initialize_mo_set(kpoints)
1582
1583 CALL get_cell(cell=cell, periodic=periodic)
1584
1585 DO i_dim = 1, 3
1586 ! we have at most 3 neigboring cells per dimension and at least one because
1587 ! the density response at Gamma is only divided to neighboring
1588 IF (periodic(i_dim) == 1) THEN
1589 cell_grid(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
1590 ELSE
1591 cell_grid(i_dim) = 1
1592 END IF
1593 END DO
1594 CALL init_cell_index_rpa(cell_grid, kpoints%cell_to_index, kpoints%index_to_cell, cell)
1595
1596 ! get S(k)
1597 CALL get_qs_env(qs_env, matrix_s=matrix_s, scf_env=scf_env, scf_control=scf_control, dft_control=dft_control)
1598
1599 NULLIFY (matrix_s_desymm)
1600 CALL dbcsr_allocate_matrix_set(matrix_s_desymm, 1)
1601 ALLOCATE (matrix_s_desymm(1)%matrix)
1602 CALL dbcsr_create(matrix=matrix_s_desymm(1)%matrix, template=matrix_s(1)%matrix, &
1603 matrix_type=dbcsr_type_no_symmetry)
1604 CALL dbcsr_desymmetrize(matrix_s(1)%matrix, matrix_s_desymm(1)%matrix)
1605
1606 CALL mat_kp_from_mat_gamma(qs_env, mat_s_kp, matrix_s_desymm(1)%matrix, kpoints, 1)
1607
1608 CALL get_kpoint_info(kpoints, nkp=nkp)
1609
1610 matrix_struct => kpoints%kp_env(1)%kpoint_env%wmat(1, 1)%matrix_struct
1611
1612 CALL cp_cfm_create(cksmat, matrix_struct)
1613 CALL cp_cfm_create(csmat, matrix_struct)
1614 CALL cp_cfm_create(cmos, matrix_struct)
1615 CALL cp_cfm_create(cwork, matrix_struct)
1616 CALL cp_fm_create(fm_work, matrix_struct)
1617
1618 nspins = dft_control%nspins
1619
1620 DO ispin = 1, nspins
1621
1622 ! get H(k)
1623 IF (my_with_xc_terms) THEN
1624 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_ks(ispin)%matrix, kpoints, ispin)
1625 ELSE
1626 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_sigma_x_minus_vxc(ispin)%matrix, &
1627 kpoints, ispin)
1628 END IF
1629
1630 DO ikp = 1, nkp
1631
1632 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 1)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1633 CALL cp_cfm_scale_and_add_fm(czero, cksmat, cone, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1634
1635 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 2)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1636 CALL cp_cfm_scale_and_add_fm(cone, cksmat, ione, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1637
1638 CALL copy_dbcsr_to_fm(mat_s_kp(ikp, 1)%matrix, fm_work)
1639 CALL cp_cfm_scale_and_add_fm(czero, csmat, cone, fm_work)
1640
1641 CALL copy_dbcsr_to_fm(mat_s_kp(ikp, 2)%matrix, fm_work)
1642 CALL cp_cfm_scale_and_add_fm(cone, csmat, ione, fm_work)
1643
1644 kp => kpoints%kp_env(ikp)%kpoint_env
1645
1646 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos, eigenvalues=eigenvalues)
1647 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1648
1649 IF (scf_env%cholesky_method == cholesky_off .OR. &
1650 qs_env%mp2_env%ri_rpa_im_time%make_overlap_mat_ao_pos_definite) THEN
1651 CALL cp_cfm_geeig_canon(cksmat, csmat, cmos, eigenvalues, cwork, scf_control%eps_eigval)
1652 ELSE
1653 CALL cp_cfm_geeig(cksmat, csmat, cmos, eigenvalues, cwork)
1654 END IF
1655
1656 CALL cp_cfm_to_fm(cmos, rmos, imos)
1657
1658 kp%mos(2, ispin)%eigenvalues = eigenvalues
1659
1660 END DO
1661
1662 END DO
1663
1664 DO ikp = 1, nkp
1665 DO i_re_im = 1, 2
1666 CALL dbcsr_deallocate_matrix(mat_ks_kp(ikp, i_re_im)%matrix)
1667 END DO
1668 END DO
1669 DEALLOCATE (mat_ks_kp)
1670
1671 DO ikp = 1, nkp
1672 DO i_re_im = 1, 2
1673 CALL dbcsr_deallocate_matrix(mat_s_kp(ikp, i_re_im)%matrix)
1674 END DO
1675 END DO
1676 DEALLOCATE (mat_s_kp)
1677
1678 CALL dbcsr_deallocate_matrix(matrix_s_desymm(1)%matrix)
1679 DEALLOCATE (matrix_s_desymm)
1680
1681 CALL cp_cfm_release(cksmat)
1682 CALL cp_cfm_release(csmat)
1683 CALL cp_cfm_release(cwork)
1684 CALL cp_cfm_release(cmos)
1685 CALL cp_fm_release(fm_work)
1686
1687 CALL timestop(handle)
1688
1689 END SUBROUTINE create_kp_and_calc_kp_orbitals
1690
1691! **************************************************************************************************
1692!> \brief ...
1693!> \param qs_env ...
1694!> \param mat_kp ...
1695!> \param mat_gamma ...
1696!> \param kpoints ...
1697!> \param ispin ...
1698!> \param real_mat_real_space ...
1699! **************************************************************************************************
1700 SUBROUTINE mat_kp_from_mat_gamma(qs_env, mat_kp, mat_gamma, kpoints, ispin, real_mat_real_space)
1701
1702 TYPE(qs_environment_type), POINTER :: qs_env
1703 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_kp
1704 TYPE(dbcsr_type) :: mat_gamma
1705 TYPE(kpoint_type), POINTER :: kpoints
1706 INTEGER :: ispin
1707 LOGICAL, INTENT(IN), OPTIONAL :: real_mat_real_space
1708
1709 CHARACTER(LEN=*), PARAMETER :: routinen = 'mat_kp_from_mat_gamma'
1710
1711 INTEGER :: handle, i_cell, i_re_im, ikp, nkp, &
1712 num_cells
1713 INTEGER, DIMENSION(3) :: periodic
1714 INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
1715 REAL(kind=dp), DIMENSION(:, :), POINTER :: xkp
1716 TYPE(cell_type), POINTER :: cell
1717 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mat_real_space
1718
1719 CALL timeset(routinen, handle)
1720
1721 CALL get_qs_env(qs_env, cell=cell)
1722 CALL get_cell(cell=cell, periodic=periodic)
1723 num_cells = 3**(periodic(1) + periodic(2) + periodic(3))
1724
1725 NULLIFY (mat_real_space)
1726 CALL dbcsr_allocate_matrix_set(mat_real_space, num_cells)
1727 DO i_cell = 1, num_cells
1728 ALLOCATE (mat_real_space(i_cell)%matrix)
1729 CALL dbcsr_create(matrix=mat_real_space(i_cell)%matrix, &
1730 template=mat_gamma)
1731 CALL dbcsr_reserve_all_blocks(mat_real_space(i_cell)%matrix)
1732 CALL dbcsr_set(mat_real_space(i_cell)%matrix, 0.0_dp)
1733 END DO
1734
1735 CALL dbcsr_copy(mat_real_space(1)%matrix, mat_gamma)
1736
1737 CALL get_mat_cell_t_from_mat_gamma(mat_real_space, qs_env, kpoints, 2, 0)
1738
1739 NULLIFY (xkp, cell_to_index)
1740 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp, cell_to_index=cell_to_index)
1741
1742 IF (ispin == 1) THEN
1743 NULLIFY (mat_kp)
1744 CALL dbcsr_allocate_matrix_set(mat_kp, nkp, 2)
1745 DO ikp = 1, nkp
1746 DO i_re_im = 1, 2
1747 ALLOCATE (mat_kp(ikp, i_re_im)%matrix)
1748 CALL dbcsr_create(matrix=mat_kp(ikp, i_re_im)%matrix, template=mat_gamma)
1749 CALL dbcsr_reserve_all_blocks(mat_kp(ikp, i_re_im)%matrix)
1750 CALL dbcsr_set(mat_kp(ikp, i_re_im)%matrix, 0.0_dp)
1751 END DO
1752 END DO
1753 END IF
1754
1755 IF (PRESENT(real_mat_real_space)) THEN
1756 CALL real_space_to_kpoint_transform_rpa(mat_kp(:, 1), mat_kp(:, 2), mat_real_space, kpoints, 0.0_dp, &
1757 real_mat_real_space)
1758 ELSE
1759 CALL real_space_to_kpoint_transform_rpa(mat_kp(:, 1), mat_kp(:, 2), mat_real_space, kpoints, 0.0_dp)
1760 END IF
1761
1762 DO i_cell = 1, num_cells
1763 CALL dbcsr_deallocate_matrix(mat_real_space(i_cell)%matrix)
1764 END DO
1765 DEALLOCATE (mat_real_space)
1766
1767 CALL timestop(handle)
1768
1769 END SUBROUTINE mat_kp_from_mat_gamma
1770
1771! **************************************************************************************************
1772!> \brief ...
1773!> \param qs_env ...
1774!> \param kpgeneral ...
1775! **************************************************************************************************
1776 SUBROUTINE get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1777 TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
1778 REAL(kind=dp), DIMENSION(:, :), POINTER :: kpgeneral
1779
1780 CHARACTER(LEN=*), PARAMETER :: routinen = 'get_kpgeneral_for_Sigma_kpoints'
1781
1782 INTEGER :: handle, i_kp_in_kp_line, i_special_kp, &
1783 i_x, ikk, j_y, k_z, n_kp_in_kp_line, &
1784 n_special_kp
1785 INTEGER, DIMENSION(:), POINTER :: nkp_grid
1786
1787 CALL timeset(routinen, handle)
1788
1789 n_special_kp = qs_env%mp2_env%ri_g0w0%n_special_kp
1790 n_kp_in_kp_line = qs_env%mp2_env%ri_g0w0%n_kp_in_kp_line
1791 IF (n_special_kp > 0) THEN
1792 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = n_kp_in_kp_line*(n_special_kp - 1) + 1
1793 ELSE
1794 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = 0
1795 END IF
1796
1797 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack = qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(1)* &
1798 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(2)* &
1799 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(3)
1800
1801 qs_env%mp2_env%ri_g0w0%nkp_self_energy = qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp + &
1802 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack
1803
1804 ALLOCATE (kpgeneral(3, qs_env%mp2_env%ri_g0w0%nkp_self_energy))
1805
1806 IF (n_special_kp > 0) THEN
1807
1808 kpgeneral(1:3, 1) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, 1)
1809
1810 ikk = 1
1811
1812 DO i_special_kp = 2, n_special_kp
1813 DO i_kp_in_kp_line = 1, n_kp_in_kp_line
1814
1815 ikk = ikk + 1
1816 kpgeneral(1:3, ikk) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1) + &
1817 REAL(i_kp_in_kp_line, kind=dp)/real(n_kp_in_kp_line, kind=dp)* &
1818 (qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp) - &
1819 qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1))
1820
1821 END DO
1822 END DO
1823
1824 ELSE
1825
1826 ikk = 0
1827
1828 END IF
1829
1830 nkp_grid => qs_env%mp2_env%ri_g0w0%kp_grid_Sigma
1831
1832 DO i_x = 1, nkp_grid(1)
1833 DO j_y = 1, nkp_grid(2)
1834 DO k_z = 1, nkp_grid(3)
1835 ikk = ikk + 1
1836 kpgeneral(1, ikk) = real(2*i_x - nkp_grid(1) - 1, kind=dp)/(2._dp*real(nkp_grid(1), kind=dp))
1837 kpgeneral(2, ikk) = real(2*j_y - nkp_grid(2) - 1, kind=dp)/(2._dp*real(nkp_grid(2), kind=dp))
1838 kpgeneral(3, ikk) = real(2*k_z - nkp_grid(3) - 1, kind=dp)/(2._dp*real(nkp_grid(3), kind=dp))
1839 END DO
1840 END DO
1841 END DO
1842
1843 CALL timestop(handle)
1844
1845 END SUBROUTINE get_kpgeneral_for_sigma_kpoints
1846
1847END 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:58
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:1905
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:147
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, 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, ecoul_1c, rho0_s_rs, rho0_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)
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