80#include "./base/base_uses.f90"
86 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'rpa_gw_kpoints_util'
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)
136 INTEGER,
INTENT(IN) :: dimen_ri, num_integ_points, jquad, nkp, &
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
147 INTEGER,
INTENT(IN) :: ikp_local
148 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
150 REAL(kind=
dp),
INTENT(IN) :: eps_filter_im_time
151 INTEGER,
INTENT(IN) :: unit_nr
153 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_mat_minv_l_kpoints, &
155 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: fm_mat_w
158 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_matrix_minv, &
159 fm_matrix_minv_vtrunc_minv
161 CHARACTER(LEN=*),
PARAMETER :: routinen =
'invert_eps_compute_W_and_Erpa_kp'
163 INTEGER :: handle, ikp
164 LOGICAL :: do_this_ikp
165 REAL(kind=
dp) :: t1, t2
166 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: trace_qomega
168 CALL timeset(routinen, handle)
172 IF (do_kpoints_from_gamma)
THEN
176 CALL transform_p_from_real_space_to_kpoints(mat_p_omega, mat_p_omega_kp, &
177 kpoints, eps_filter_im_time, jquad)
179 ALLOCATE (trace_qomega(dimen_ri))
181 IF (unit_nr > 0)
WRITE (unit_nr,
'(/T3,A,1X,I3)') &
182 'GW_INFO| Computing chi and W frequency point', jquad
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
191 CALL compute_q_kp_rpa(cfm_mat_q, &
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)
200 CALL cholesky_decomp_q(cfm_mat_q, para_env, trace_qomega, dimen_ri)
203 CALL frequency_and_kpoint_integration(erpa, cfm_mat_q, para_env, trace_qomega, &
204 dimen_ri, wj(jquad), kpoints%wkp(ikp))
206 IF (do_gw_im_time)
THEN
209 IF (do_ri_sigma_x .AND. jquad == 1 .AND. count_ev_sc_gw == 1 .AND. do_kpoints_from_gamma)
THEN
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.)
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)
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)
233 DEALLOCATE (trace_qomega)
237 IF (unit_nr > 0)
WRITE (unit_nr,
'(T6,A,T56,F25.1)')
'Execution time (s):', t2 - t1
239 CALL timestop(handle)
248 SUBROUTINE deallocate_kp_matrices(fm_matrix_L_kpoints, fm_mat_Minv_L_kpoints)
250 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_matrix_l_kpoints, &
251 fm_mat_minv_l_kpoints
253 CHARACTER(LEN=*),
PARAMETER :: routinen =
'deallocate_kp_matrices'
257 CALL timeset(routinen, handle)
262 CALL timestop(handle)
264 END SUBROUTINE deallocate_kp_matrices
275 REAL(kind=
dp) :: threshold, exponent
276 REAL(kind=
dp),
OPTIONAL :: min_eigval
278 CHARACTER(LEN=*),
PARAMETER :: routinen =
'cp_cfm_power'
279 COMPLEX(KIND=dp),
PARAMETER :: czero = cmplx(0.0_dp, 0.0_dp, kind=
dp)
281 COMPLEX(KIND=dp),
ALLOCATABLE,
DIMENSION(:) :: eigenvalues_exponent
282 INTEGER :: handle, i, ncol_global, nrow_global
283 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eigenvalues
286 CALL timeset(routinen, handle)
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
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)
306 IF (
PRESENT(min_eigval))
THEN
307 eigenvalues_exponent(i) = cmplx(min_eigval, 0.0_dp, kind=
dp)
309 eigenvalues_exponent(i) = czero
317 cfm_work, cfm_work,
z_zero, matrix)
319 DEALLOCATE (eigenvalues, eigenvalues_exponent)
323 CALL timestop(handle)
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)
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
351 LOGICAL,
INTENT(IN) :: make_chi_pos_definite
353 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_Q_kp_RPA'
359 CALL timeset(routinen, handle)
367 CALL cp_fm_create(fm_mat_work, fm_mat_l_re%matrix_struct)
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)
376 IF (make_chi_pos_definite)
THEN
377 CALL cp_cfm_power(cfm_mat_q, threshold=0.0_dp, exponent=1.0_dp)
385 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_q, cfm_mat_l, &
389 CALL parallel_gemm(
'C',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_l, cfm_mat_work, &
396 CALL timestop(handle)
398 END SUBROUTINE compute_q_kp_rpa
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)
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
417 INTEGER,
INTENT(IN) :: ikp, nkp, ikp_local
420 CHARACTER(LEN=*),
PARAMETER :: routinen =
'mat_P_to_subgroup'
422 INTEGER :: handle, jkp
424 TYPE(
dbcsr_type),
POINTER :: mat_p_omega_im, mat_p_omega_re
426 CALL timeset(routinen, handle)
428 IF (ikp_local == -1)
THEN
430 mat_p_omega_re => mat_p_omega_kp(1, ikp)%matrix
435 mat_p_omega_im => mat_p_omega_kp(2, ikp)%matrix
446 mat_p_omega_re => mat_p_omega_kp(1, jkp)%matrix
453 IF (ikp_local == jkp)
THEN
469 mat_p_omega_im => mat_p_omega_kp(2, jkp)%matrix
476 IF (ikp_local == jkp)
THEN
494 CALL timestop(handle)
496 END SUBROUTINE mat_p_to_subgroup
505 SUBROUTINE cholesky_decomp_q(cfm_mat_Q, para_env, trace_Qomega, dimen_RI)
509 REAL(kind=
dp),
DIMENSION(:),
INTENT(OUT) :: trace_qomega
510 INTEGER,
INTENT(IN) :: dimen_ri
512 CHARACTER(LEN=*),
PARAMETER :: routinen =
'cholesky_decomp_Q'
514 INTEGER :: handle, i_global, iib, info_chol, &
515 j_global, jjb, ncol_local, nrow_local
516 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
519 CALL timeset(routinen, handle)
529 nrow_local=nrow_local, &
530 ncol_local=ncol_local, &
531 row_indices=row_indices, &
532 col_indices=col_indices)
535 trace_qomega = 0.0_dp
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
548 CALL para_env%sum(trace_qomega)
554 cpassert(info_chol == 0)
559 CALL timestop(handle)
561 END SUBROUTINE cholesky_decomp_q
573 SUBROUTINE frequency_and_kpoint_integration(Erpa, cfm_mat_Q, para_env, trace_Qomega, &
574 dimen_RI, freq_weight, kp_weight)
576 REAL(kind=
dp),
INTENT(INOUT) :: erpa
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
583 CHARACTER(LEN=*),
PARAMETER :: routinen =
'frequency_and_kpoint_integration'
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
591 CALL timeset(routinen, handle)
595 nrow_local=nrow_local, &
596 ncol_local=ncol_local, &
597 row_indices=row_indices, &
598 col_indices=col_indices)
600 ALLOCATE (q_log(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)))
613 CALL para_env%sum(q_log)
617 IF (
modulo(iib, para_env%num_pe) /= para_env%mepos) cycle
619 fcomega = fcomega + (q_log(iib) - trace_qomega(iib))/2.0_dp
622 erpa = erpa + fcomega*freq_weight*kp_weight
626 CALL timestop(handle)
628 END SUBROUTINE frequency_and_kpoint_integration
636 SUBROUTINE get_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
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
643 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_dummys'
647 CALL timeset(routinen, handle)
649 ALLOCATE (weights_cos_tf_w_to_t_dummy(1, 1))
650 ALLOCATE (tj_dummy(1))
651 ALLOCATE (tau_tj_dummy(1))
654 tau_tj_dummy(1) = 0.0_dp
655 weights_cos_tf_w_to_t_dummy(1, 1) = 1.0_dp
657 CALL timestop(handle)
667 SUBROUTINE release_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
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
674 CHARACTER(LEN=*),
PARAMETER :: routinen =
'release_dummys'
678 CALL timeset(routinen, handle)
680 DEALLOCATE (weights_cos_tf_w_to_t_dummy)
681 DEALLOCATE (tj_dummy)
682 DEALLOCATE (tau_tj_dummy)
684 CALL timestop(handle)
697 TYPE(
dbcsr_p_type),
DIMENSION(:),
INTENT(IN) :: mat_p_omega
700 INTEGER,
INTENT(IN) :: jquad, unit_nr
702 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_mat_cell_T_from_mat_gamma'
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, &
712 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
713 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
718 CALL timeset(routinen, handle)
720 NULLIFY (cell, particle_set)
722 particle_set=particle_set)
723 CALL get_cell(cell=cell, h=hmat, periodic=periodic)
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)
731 cell_grid_p(i_dim) = 1
738 index_to_cell_p => kpoints%index_to_cell
740 num_cells_p =
SIZE(index_to_cell_p, 2)
742 num_integ_points =
SIZE(mat_p_omega, 1)
746 DO i_cell = 2, num_cells_p
748 mat_p_omega(1)%matrix)
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'
760 WRITE (unit_nr,
'(T3,A,T81)') &
761 'GW_INFO| Make chi(iw,k) positive definite? FALSE'
766 DO i_cell = 1, num_cells_p
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)
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)
785 IF (abs_rab_cell_i > abs_rab_cell_j + 1.0e-6_dp)
THEN
786 i_cell_is_the_minimum_image_cell = .false.
790 IF (.NOT. i_cell_is_the_minimum_image_cell)
THEN
791 data_block(:, :) = data_block(:, :)*0.0_dp
799 CALL timestop(handle)
811 SUBROUTINE transform_p_from_real_space_to_kpoints(mat_P_omega, mat_P_omega_kp, &
812 kpoints, eps_filter_im_time, jquad)
814 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
816 REAL(kind=
dp),
INTENT(IN) :: eps_filter_im_time
817 INTEGER,
INTENT(IN) :: jquad
819 CHARACTER(LEN=*),
PARAMETER :: routinen =
'transform_P_from_real_space_to_kpoints'
821 INTEGER :: handle, icell, nkp, num_integ_points
823 CALL timeset(routinen, handle)
825 num_integ_points =
SIZE(mat_p_omega, 1)
826 nkp =
SIZE(mat_p_omega, 2)
829 kpoints, eps_filter_im_time)
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)
836 CALL timestop(handle)
838 END SUBROUTINE transform_p_from_real_space_to_kpoints
850 kpoints, eps_filter_im_time, real_mat_real_space)
852 TYPE(
dbcsr_p_type),
DIMENSION(:),
INTENT(INOUT) :: real_mat_kp, imag_mat_kp, mat_real_space
854 REAL(kind=
dp),
INTENT(IN) :: eps_filter_im_time
855 LOGICAL,
INTENT(IN),
OPTIONAL :: real_mat_real_space
857 CHARACTER(LEN=*),
PARAMETER :: routinen =
'real_space_to_kpoint_transform_rpa'
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
867 CALL timeset(routinen, handle)
869 my_real_mat_real_space = .true.
870 IF (
PRESENT(real_mat_real_space)) my_real_mat_real_space = real_mat_real_space
873 template=real_mat_kp(1)%matrix, &
874 matrix_type=dbcsr_type_no_symmetry)
881 NULLIFY (index_to_cell)
882 index_to_cell => kpoints%index_to_cell
884 num_cells =
SIZE(index_to_cell, 2)
886 cpassert(
SIZE(mat_real_space) >= num_cells/2 + 1)
893 CALL dbcsr_set(real_mat_kp(ik)%matrix, 0.0_dp)
894 CALL dbcsr_set(imag_mat_kp(ik)%matrix, 0.0_dp)
896 DO i_cell = 1, num_cells/2 + 1
898 cell(:) = index_to_cell(:, i_cell)
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)
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)
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)
912 IF (.NOT. (cell(1) == 0 .AND. cell(2) == 0 .AND. cell(3) == 0))
THEN
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)
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)
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)
940 CALL timestop(handle)
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
957 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: block_to_compute, data_block
964 NULLIFY (block_to_compute)
966 row=row, col=col, block=block_to_compute, found=found)
970 block_to_compute(:, :) = alpha*block_to_compute(:, :) + beta*data_block(:, :)
975 END SUBROUTINE dbcsr_add_local
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)
1001 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: fm_mat_w_tau
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
1013 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: wkp_w
1015 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_Wc_real_space_tau_GW'
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
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, &
1033 CALL timeset(routinen, handle)
1035 CALL timeset(routinen//
"_1", handle2)
1050 CALL cp_fm_create(fm_mat_work_global, fm_mat_w_tau(1)%matrix_struct)
1053 CALL cp_fm_create(fm_mat_work_local, cfm_mat_q%matrix_struct)
1056 CALL timestop(handle2)
1058 CALL timeset(routinen//
"_2", handle2)
1068 nrow_local=nrow_local, &
1069 ncol_local=ncol_local, &
1070 row_indices=row_indices, &
1071 col_indices=col_indices)
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
1083 CALL timestop(handle2)
1085 CALL timeset(routinen//
"_3", handle2)
1088 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_q, cfm_mat_l, &
1092 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_l, cfm_mat_work, &
1095 CALL timestop(handle2)
1097 CALL timeset(routinen//
"_4", handle2)
1100 index_to_cell => kpoints%index_to_cell
1101 num_cells =
SIZE(index_to_cell, 2)
1105 ALLOCATE (atom_from_ri_index(dimen_ri))
1109 NULLIFY (cell, particle_set)
1110 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1116 nrow_local=nrow_local, &
1117 ncol_local=ncol_local, &
1118 row_indices=row_indices, &
1119 col_indices=col_indices)
1121 DO irow = 1, nrow_local
1122 DO jcol = 1, ncol_local
1124 iatom = atom_from_ri_index(row_indices(irow))
1125 jatom = atom_from_ri_index(col_indices(jcol))
1127 IF (iatom .NE. iatom_old .OR. jatom .NE. jatom_old)
THEN
1131 num_cells, iatom, jatom, xkp(1:3, ikp), wkp_w(ikp), &
1132 cell, index_to_cell, hmat, particle_set)
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))
1142 fm_mat_work_local%local_data(irow, jcol) = fm_mat_work_local%local_data(irow, jcol) + contribution
1147 CALL timestop(handle2)
1149 CALL timeset(routinen//
"_5", handle2)
1151 IF (ikp_local == -1)
THEN
1155 DO iquad = 1, num_integ_points
1159 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1161 IF (jquad == 1 .AND. ikp == 1)
THEN
1162 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
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)
1173 CALL para_env%sync()
1175 IF (ikp_local == jkp)
THEN
1181 CALL para_env%sync()
1183 DO iquad = 1, num_integ_points
1187 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1189 IF (jquad == 1 .AND. jkp == 1)
THEN
1190 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1194 matrix_b=fm_mat_work_global)
1208 DEALLOCATE (atom_from_ri_index)
1210 CALL timestop(handle2)
1212 CALL timestop(handle)
1214 END SUBROUTINE compute_wc_real_space_tau_gw
1224 SUBROUTINE wc_to_minv_wc_minv(fm_mat_W, fm_matrix_Minv, para_env, dimen_RI, num_integ_points)
1226 TYPE(
cp_fm_type),
DIMENSION(:, :) :: fm_matrix_minv
1228 INTEGER :: dimen_ri, num_integ_points
1230 CHARACTER(LEN=*),
PARAMETER :: routinen =
'Wc_to_Minv_Wc_Minv'
1232 INTEGER :: handle, jquad
1233 TYPE(
cp_fm_type) :: fm_work_minv, fm_work_minv_w
1235 CALL timeset(routinen, handle)
1237 CALL cp_fm_create(fm_work_minv, fm_mat_w(1)%matrix_struct)
1240 CALL cp_fm_create(fm_work_minv_w, fm_mat_w(1)%matrix_struct)
1242 DO jquad = 1, num_integ_points
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))
1255 CALL timestop(handle)
1257 END SUBROUTINE wc_to_minv_wc_minv
1271 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:), &
1272 INTENT(OUT) :: wkp_w, wkp_v
1274 REAL(kind=
dp),
DIMENSION(3, 3) :: h_inv
1275 INTEGER,
DIMENSION(3) :: periodic
1277 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_wkp_W'
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
1287 CALL timeset(routinen, handle)
1289 kpoint_weights_w_method = qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method
1303 ALLOCATE (wkp_v(nkp), wkp_w(nkp))
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(:)
1322 exp_kpoints = qs_env%mp2_env%ri_rpa_im_time%exp_tailored_weights
1325 IF (sum(periodic) == 2) exp_kpoints = -1.0_dp
1332 IF (periodic(1) == 1)
THEN
1337 IF (periodic(2) == 1)
THEN
1342 IF (periodic(3) == 1)
THEN
1350 weight = 1.0_dp/(real(n_x,
dp)*real(n_y,
dp)*real(n_z,
dp))
1355 IF (periodic(1) == 1)
THEN
1356 x_vec(1) = (real(i_x - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
1360 IF (periodic(2) == 1)
THEN
1361 x_vec(2) = (real(j_y - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
1365 IF (periodic(3) == 1)
THEN
1366 x_vec(3) = (real(k_z - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
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)
1379 num_lin_eqs = nkp + 2
1381 ALLOCATE (matrix_lin_eqs(num_lin_eqs, num_lin_eqs))
1382 matrix_lin_eqs(:, :) = 0.0_dp
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
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
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)
1398 CALL invmat(matrix_lin_eqs, info)
1402 ALLOCATE (right_side(num_lin_eqs))
1404 right_side(nkp + 1) = 1.0_dp
1406 right_side(nkp + 2) = integral
1408 ALLOCATE (wkp_tmp(num_lin_eqs))
1410 wkp_tmp(1:num_lin_eqs) = matmul(matrix_lin_eqs, right_side)
1412 wkp_w(1:nkp) = wkp_tmp(1:nkp)
1414 DEALLOCATE (matrix_lin_eqs, right_side, wkp_tmp)
1418 CALL timestop(handle)
1429 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: eigenval_kp
1431 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_bandstruc_and_k_dependent_MOs'
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
1440 CALL timeset(routinen, handle)
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, &
1447 nkp_grid_g(1:3) = (/1, 1, 1/)
1449 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
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))
1455 IF (qs_env%mp2_env%ri_g0w0%do_kpoints_Sigma)
THEN
1458 CALL get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
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)
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.)
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)
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)
1475 DEALLOCATE (eigenval_kp)
1477 ALLOCATE (eigenval_kp(nmo, kpoints_sigma%nkp, nspins))
1479 DO ikp = 1, kpoints_sigma%nkp
1481 DO ispin = 1, nspins
1483 ev => kpoints_sigma%kp_env(ikp)%kpoint_env%mos(1, ispin)%eigenvalues
1485 eigenval_kp(:, ikp, ispin) = ev(:)
1491 DEALLOCATE (kpgeneral)
1495 CALL release_hfx_stuff(qs_env)
1497 CALL timestop(handle)
1505 SUBROUTINE release_hfx_stuff(qs_env)
1508 IF (
ASSOCIATED(qs_env%x_data) .AND. .NOT. qs_env%mp2_env%ri_g0w0%do_ri_Sigma_x)
THEN
1512 END SUBROUTINE release_hfx_stuff
1524 SUBROUTINE create_kp_and_calc_kp_orbitals(qs_env, kpoints, scheme, &
1525 group_size_ext, mp_grid, kpgeneral, with_xc_terms)
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
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)
1540 INTEGER :: handle, i_dim, i_re_im, ikp, ispin, nkp, &
1542 INTEGER,
DIMENSION(3) :: cell_grid, periodic
1543 LOGICAL :: my_with_xc_terms
1544 REAL(kind=
dp),
DIMENSION(:),
POINTER :: eigenvalues
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
1559 CALL timeset(routinen, handle)
1561 my_with_xc_terms = .true.
1562 IF (
PRESENT(with_xc_terms)) my_with_xc_terms = with_xc_terms
1565 para_env=para_env, &
1566 blacs_env=blacs_env, &
1567 matrix_s=matrix_s, &
1569 scf_control=scf_control, &
1574 group_size_ext=group_size_ext)
1583 CALL get_cell(cell=cell, periodic=periodic)
1588 IF (periodic(i_dim) == 1)
THEN
1589 cell_grid(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
1591 cell_grid(i_dim) = 1
1597 CALL get_qs_env(qs_env, matrix_s=matrix_s, scf_env=scf_env, scf_control=scf_control, dft_control=dft_control)
1599 NULLIFY (matrix_s_desymm)
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)
1610 matrix_struct => kpoints%kp_env(1)%kpoint_env%wmat(1, 1)%matrix_struct
1618 nspins = dft_control%nspins
1620 DO ispin = 1, nspins
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)
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, &
1632 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 1)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1635 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 2)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1644 kp => kpoints%kp_env(ikp)%kpoint_env
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)
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)
1653 CALL cp_cfm_geeig(cksmat, csmat, cmos, eigenvalues, cwork)
1658 kp%mos(2, ispin)%eigenvalues = eigenvalues
1669 DEALLOCATE (mat_ks_kp)
1676 DEALLOCATE (mat_s_kp)
1679 DEALLOCATE (matrix_s_desymm)
1687 CALL timestop(handle)
1689 END SUBROUTINE create_kp_and_calc_kp_orbitals
1707 LOGICAL,
INTENT(IN),
OPTIONAL :: real_mat_real_space
1709 CHARACTER(LEN=*),
PARAMETER :: routinen =
'mat_kp_from_mat_gamma'
1711 INTEGER :: handle, i_cell, i_re_im, ikp, nkp, &
1713 INTEGER,
DIMENSION(3) :: periodic
1714 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
1715 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: xkp
1717 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: mat_real_space
1719 CALL timeset(routinen, handle)
1722 CALL get_cell(cell=cell, periodic=periodic)
1723 num_cells = 3**(periodic(1) + periodic(2) + periodic(3))
1725 NULLIFY (mat_real_space)
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, &
1732 CALL dbcsr_set(mat_real_space(i_cell)%matrix, 0.0_dp)
1735 CALL dbcsr_copy(mat_real_space(1)%matrix, mat_gamma)
1739 NULLIFY (xkp, cell_to_index)
1740 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp, cell_to_index=cell_to_index)
1742 IF (ispin == 1)
THEN
1747 ALLOCATE (mat_kp(ikp, i_re_im)%matrix)
1748 CALL dbcsr_create(matrix=mat_kp(ikp, i_re_im)%matrix, template=mat_gamma)
1750 CALL dbcsr_set(mat_kp(ikp, i_re_im)%matrix, 0.0_dp)
1755 IF (
PRESENT(real_mat_real_space))
THEN
1757 real_mat_real_space)
1762 DO i_cell = 1, num_cells
1765 DEALLOCATE (mat_real_space)
1767 CALL timestop(handle)
1776 SUBROUTINE get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1778 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: kpgeneral
1780 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_kpgeneral_for_Sigma_kpoints'
1782 INTEGER :: handle, i_kp_in_kp_line, i_special_kp, &
1783 i_x, ikk, j_y, k_z, n_kp_in_kp_line, &
1785 INTEGER,
DIMENSION(:),
POINTER :: nkp_grid
1787 CALL timeset(routinen, handle)
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
1794 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = 0
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)
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
1804 ALLOCATE (kpgeneral(3, qs_env%mp2_env%ri_g0w0%nkp_self_energy))
1806 IF (n_special_kp > 0)
THEN
1808 kpgeneral(1:3, 1) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, 1)
1812 DO i_special_kp = 2, n_special_kp
1813 DO i_kp_in_kp_line = 1, n_kp_in_kp_line
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))
1830 nkp_grid => qs_env%mp2_env%ri_g0w0%kp_grid_Sigma
1832 DO i_x = 1, nkp_grid(1)
1833 DO j_y = 1, nkp_grid(2)
1834 DO k_z = 1, nkp_grid(3)
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))
1843 CALL timestop(handle)
1845 END SUBROUTINE get_kpgeneral_for_sigma_kpoints
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.
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
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
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.
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
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.
subroutine, public hfx_release(x_data)
This routine deallocates all data structures
Defines the basic variable types.
integer, parameter, public dp
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.
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
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.
subroutine, public invmat(a, info)
returns inverse of matrix using the lapack routines DGETRF and DGETRI
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.
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.
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.
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
Keeps information about a specific k-point.
Contains information about kpoints.
stores all the informations relevant to an mpi environment