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'
280 COMPLEX(KIND=dp),
ALLOCATABLE,
DIMENSION(:) :: eigenvalues_exponent
281 INTEGER :: handle, i, ncol_global, nrow_global
282 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: eigenvalues
285 CALL timeset(routinen, handle)
291 CALL cp_cfm_get_info(matrix, nrow_global=nrow_global, ncol_global=ncol_global)
292 cpassert(nrow_global == ncol_global)
293 ALLOCATE (eigenvalues(nrow_global))
294 eigenvalues(:) = 0.0_dp
295 ALLOCATE (eigenvalues_exponent(nrow_global))
296 eigenvalues_exponent(:) =
z_zero
301 DO i = 1, nrow_global
302 IF (eigenvalues(i) > threshold)
THEN
303 eigenvalues_exponent(i) = cmplx((eigenvalues(i))**(0.5_dp*exponent), threshold, kind=
dp)
305 IF (
PRESENT(min_eigval))
THEN
306 eigenvalues_exponent(i) = cmplx(min_eigval, 0.0_dp, kind=
dp)
308 eigenvalues_exponent(i) =
z_zero
316 cfm_work, cfm_work,
z_zero, matrix)
318 DEALLOCATE (eigenvalues, eigenvalues_exponent)
322 CALL timestop(handle)
340 SUBROUTINE compute_q_kp_rpa(cfm_mat_Q, mat_P_omega_kp, fm_mat_L_re, fm_mat_L_im, &
341 fm_mat_RI_global_work, dimen_RI, ikp, nkp, ikp_local, para_env, &
342 make_chi_pos_definite)
345 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
INTENT(INOUT) :: mat_p_omega_kp
346 TYPE(
cp_fm_type) :: fm_mat_l_re, fm_mat_l_im, &
347 fm_mat_ri_global_work
348 INTEGER,
INTENT(IN) :: dimen_ri, ikp, nkp, ikp_local
350 LOGICAL,
INTENT(IN) :: make_chi_pos_definite
352 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_Q_kp_RPA'
358 CALL timeset(routinen, handle)
366 CALL cp_fm_create(fm_mat_work, fm_mat_l_re%matrix_struct)
371 CALL mat_p_to_subgroup(mat_p_omega_kp, fm_mat_ri_global_work, &
372 fm_mat_work, cfm_mat_q, ikp, nkp, ikp_local, para_env)
375 IF (make_chi_pos_definite)
THEN
376 CALL cp_cfm_power(cfm_mat_q, threshold=0.0_dp, exponent=1.0_dp)
384 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_q, cfm_mat_l, &
388 CALL parallel_gemm(
'C',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_l, cfm_mat_work, &
395 CALL timestop(handle)
397 END SUBROUTINE compute_q_kp_rpa
410 SUBROUTINE mat_p_to_subgroup(mat_P_omega_kp, fm_mat_RI_global_work, &
411 fm_mat_work, cfm_mat_Q, ikp, nkp, ikp_local, para_env)
413 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
INTENT(INOUT) :: mat_p_omega_kp
414 TYPE(
cp_fm_type),
INTENT(INOUT) :: fm_mat_ri_global_work, fm_mat_work
416 INTEGER,
INTENT(IN) :: ikp, nkp, ikp_local
419 CHARACTER(LEN=*),
PARAMETER :: routinen =
'mat_P_to_subgroup'
421 INTEGER :: handle, jkp
423 TYPE(
dbcsr_type),
POINTER :: mat_p_omega_im, mat_p_omega_re
425 CALL timeset(routinen, handle)
427 IF (ikp_local == -1)
THEN
429 mat_p_omega_re => mat_p_omega_kp(1, ikp)%matrix
434 mat_p_omega_im => mat_p_omega_kp(2, ikp)%matrix
445 mat_p_omega_re => mat_p_omega_kp(1, jkp)%matrix
452 IF (ikp_local == jkp)
THEN
468 mat_p_omega_im => mat_p_omega_kp(2, jkp)%matrix
475 IF (ikp_local == jkp)
THEN
493 CALL timestop(handle)
495 END SUBROUTINE mat_p_to_subgroup
504 SUBROUTINE cholesky_decomp_q(cfm_mat_Q, para_env, trace_Qomega, dimen_RI)
508 REAL(kind=
dp),
DIMENSION(:),
INTENT(OUT) :: trace_qomega
509 INTEGER,
INTENT(IN) :: dimen_ri
511 CHARACTER(LEN=*),
PARAMETER :: routinen =
'cholesky_decomp_Q'
513 INTEGER :: handle, i_global, iib, info_chol, &
514 j_global, jjb, ncol_local, nrow_local
515 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
518 CALL timeset(routinen, handle)
528 nrow_local=nrow_local, &
529 ncol_local=ncol_local, &
530 row_indices=row_indices, &
531 col_indices=col_indices)
534 trace_qomega = 0.0_dp
537 DO jjb = 1, ncol_local
538 j_global = col_indices(jjb)
539 DO iib = 1, nrow_local
540 i_global = row_indices(iib)
541 IF (j_global == i_global .AND. i_global <= dimen_ri)
THEN
542 trace_qomega(i_global) = real(cfm_mat_q%local_data(iib, jjb))
543 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) +
z_one
547 CALL para_env%sum(trace_qomega)
553 cpassert(info_chol == 0)
558 CALL timestop(handle)
560 END SUBROUTINE cholesky_decomp_q
572 SUBROUTINE frequency_and_kpoint_integration(Erpa, cfm_mat_Q, para_env, trace_Qomega, &
573 dimen_RI, freq_weight, kp_weight)
575 REAL(kind=
dp),
INTENT(INOUT) :: erpa
578 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: trace_qomega
579 INTEGER,
INTENT(IN) :: dimen_ri
580 REAL(kind=
dp),
INTENT(IN) :: freq_weight, kp_weight
582 CHARACTER(LEN=*),
PARAMETER :: routinen =
'frequency_and_kpoint_integration'
584 INTEGER :: handle, i_global, iib, j_global, jjb, &
585 ncol_local, nrow_local
586 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
587 REAL(kind=
dp) :: fcomega
588 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:) :: q_log
590 CALL timeset(routinen, handle)
594 nrow_local=nrow_local, &
595 ncol_local=ncol_local, &
596 row_indices=row_indices, &
597 col_indices=col_indices)
599 ALLOCATE (q_log(dimen_ri))
603 DO jjb = 1, ncol_local
604 j_global = col_indices(jjb)
605 DO iib = 1, nrow_local
606 i_global = row_indices(iib)
607 IF (j_global == i_global .AND. i_global <= dimen_ri)
THEN
608 q_log(i_global) = 2.0_dp*log(real(cfm_mat_q%local_data(iib, jjb)))
612 CALL para_env%sum(q_log)
616 IF (
modulo(iib, para_env%num_pe) /= para_env%mepos) cycle
618 fcomega = fcomega + (q_log(iib) - trace_qomega(iib))/2.0_dp
621 erpa = erpa + fcomega*freq_weight*kp_weight
625 CALL timestop(handle)
627 END SUBROUTINE frequency_and_kpoint_integration
635 SUBROUTINE get_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
637 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:), &
638 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
639 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
640 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
642 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_dummys'
646 CALL timeset(routinen, handle)
648 ALLOCATE (weights_cos_tf_w_to_t_dummy(1, 1))
649 ALLOCATE (tj_dummy(1))
650 ALLOCATE (tau_tj_dummy(1))
653 tau_tj_dummy(1) = 0.0_dp
654 weights_cos_tf_w_to_t_dummy(1, 1) = 1.0_dp
656 CALL timestop(handle)
666 SUBROUTINE release_dummys(tj_dummy, tau_tj_dummy, weights_cos_tf_w_to_t_dummy)
668 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:), &
669 INTENT(INOUT) :: tj_dummy, tau_tj_dummy
670 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
671 INTENT(INOUT) :: weights_cos_tf_w_to_t_dummy
673 CHARACTER(LEN=*),
PARAMETER :: routinen =
'release_dummys'
677 CALL timeset(routinen, handle)
679 DEALLOCATE (weights_cos_tf_w_to_t_dummy)
680 DEALLOCATE (tj_dummy)
681 DEALLOCATE (tau_tj_dummy)
683 CALL timestop(handle)
696 TYPE(
dbcsr_p_type),
DIMENSION(:),
INTENT(IN) :: mat_p_omega
699 INTEGER,
INTENT(IN) :: jquad, unit_nr
701 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_mat_cell_T_from_mat_gamma'
703 INTEGER :: col, handle, i_cell, i_dim, j_cell, &
704 num_cells_p, num_integ_points, row
705 INTEGER,
DIMENSION(3) :: cell_grid_p, periodic
706 INTEGER,
DIMENSION(:, :),
POINTER :: index_to_cell_p
707 LOGICAL :: i_cell_is_the_minimum_image_cell
708 REAL(kind=
dp) :: abs_rab_cell_i, abs_rab_cell_j
709 REAL(kind=
dp),
DIMENSION(3) :: cell_vector, cell_vector_j, rab_cell_i, &
711 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
712 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: data_block
717 CALL timeset(routinen, handle)
719 NULLIFY (cell, particle_set)
721 particle_set=particle_set)
722 CALL get_cell(cell=cell, h=hmat, periodic=periodic)
727 IF (periodic(i_dim) == 1)
THEN
728 cell_grid_p(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
730 cell_grid_p(i_dim) = 1
737 index_to_cell_p => kpoints%index_to_cell
739 num_cells_p =
SIZE(index_to_cell_p, 2)
741 num_integ_points =
SIZE(mat_p_omega, 1)
745 DO i_cell = 2, num_cells_p
747 mat_p_omega(1)%matrix)
750 IF (jquad == 1 .AND. unit_nr > 0)
THEN
751 WRITE (unit_nr,
'(T3,A,T66,ES15.2)')
'GW_INFO| RI regularization parameter: ', &
752 qs_env%mp2_env%ri_rpa_im_time%regularization_RI
753 WRITE (unit_nr,
'(T3,A,T66,ES15.2)')
'GW_INFO| eps_eigval_S: ', &
754 qs_env%mp2_env%ri_rpa_im_time%eps_eigval_S
755 IF (qs_env%mp2_env%ri_rpa_im_time%make_chi_pos_definite)
THEN
756 WRITE (unit_nr,
'(T3,A,T81)') &
757 'GW_INFO| Make chi(iw,k) positive definite? TRUE'
759 WRITE (unit_nr,
'(T3,A,T81)') &
760 'GW_INFO| Make chi(iw,k) positive definite? FALSE'
765 DO i_cell = 1, num_cells_p
771 cell_vector(1:3) = matmul(hmat, real(index_to_cell_p(1:3, i_cell),
dp))
772 rab_cell_i(1:3) =
pbc(particle_set(row)%r(1:3), cell) - &
773 (
pbc(particle_set(col)%r(1:3), cell) + cell_vector(1:3))
774 abs_rab_cell_i = sqrt(rab_cell_i(1)**2 + rab_cell_i(2)**2 + rab_cell_i(3)**2)
777 i_cell_is_the_minimum_image_cell = .true.
778 DO j_cell = 1, num_cells_p
779 cell_vector_j(1:3) = matmul(hmat, real(index_to_cell_p(1:3, j_cell),
dp))
780 rab_cell_j(1:3) =
pbc(particle_set(row)%r(1:3), cell) - &
781 (
pbc(particle_set(col)%r(1:3), cell) + cell_vector_j(1:3))
782 abs_rab_cell_j = sqrt(rab_cell_j(1)**2 + rab_cell_j(2)**2 + rab_cell_j(3)**2)
784 IF (abs_rab_cell_i > abs_rab_cell_j + 1.0e-6_dp)
THEN
785 i_cell_is_the_minimum_image_cell = .false.
789 IF (.NOT. i_cell_is_the_minimum_image_cell)
THEN
790 data_block(:, :) = data_block(:, :)*0.0_dp
798 CALL timestop(handle)
810 SUBROUTINE transform_p_from_real_space_to_kpoints(mat_P_omega, mat_P_omega_kp, &
811 kpoints, eps_filter_im_time, jquad)
813 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
INTENT(INOUT) :: mat_p_omega, mat_p_omega_kp
815 REAL(kind=
dp),
INTENT(IN) :: eps_filter_im_time
816 INTEGER,
INTENT(IN) :: jquad
818 CHARACTER(LEN=*),
PARAMETER :: routinen =
'transform_P_from_real_space_to_kpoints'
820 INTEGER :: handle, icell, nkp, num_integ_points
822 CALL timeset(routinen, handle)
824 num_integ_points =
SIZE(mat_p_omega, 1)
825 nkp =
SIZE(mat_p_omega, 2)
828 kpoints, eps_filter_im_time)
830 DO icell = 1,
SIZE(mat_p_omega, 2)
831 CALL dbcsr_set(mat_p_omega(jquad, icell)%matrix, 0.0_dp)
832 CALL dbcsr_filter(mat_p_omega(jquad, icell)%matrix, 1.0_dp)
835 CALL timestop(handle)
837 END SUBROUTINE transform_p_from_real_space_to_kpoints
849 kpoints, eps_filter_im_time, real_mat_real_space)
851 TYPE(
dbcsr_p_type),
DIMENSION(:),
INTENT(INOUT) :: real_mat_kp, imag_mat_kp, mat_real_space
853 REAL(kind=
dp),
INTENT(IN) :: eps_filter_im_time
854 LOGICAL,
INTENT(IN),
OPTIONAL :: real_mat_real_space
856 CHARACTER(LEN=*),
PARAMETER :: routinen =
'real_space_to_kpoint_transform_rpa'
858 INTEGER :: handle, i_cell, ik, nkp, num_cells
859 INTEGER,
DIMENSION(3) :: cell
860 INTEGER,
DIMENSION(:, :),
POINTER :: index_to_cell
861 LOGICAL :: my_real_mat_real_space
862 REAL(kind=
dp) :: arg, coskl, sinkl
863 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: xkp
866 CALL timeset(routinen, handle)
868 my_real_mat_real_space = .true.
869 IF (
PRESENT(real_mat_real_space)) my_real_mat_real_space = real_mat_real_space
872 template=real_mat_kp(1)%matrix, &
873 matrix_type=dbcsr_type_no_symmetry)
880 NULLIFY (index_to_cell)
881 index_to_cell => kpoints%index_to_cell
883 num_cells =
SIZE(index_to_cell, 2)
885 cpassert(
SIZE(mat_real_space) >= num_cells/2 + 1)
892 CALL dbcsr_set(real_mat_kp(ik)%matrix, 0.0_dp)
893 CALL dbcsr_set(imag_mat_kp(ik)%matrix, 0.0_dp)
895 DO i_cell = 1, num_cells/2 + 1
897 cell(:) = index_to_cell(:, i_cell)
899 arg = real(cell(1),
dp)*xkp(1, ik) + real(cell(2),
dp)*xkp(2, ik) + real(cell(3),
dp)*xkp(3, ik)
900 coskl = cos(
twopi*arg)
901 sinkl = sin(
twopi*arg)
903 IF (my_real_mat_real_space)
THEN
904 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
905 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, sinkl)
907 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, -sinkl)
908 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_real_space(i_cell)%matrix, 1.0_dp, coskl)
911 IF (.NOT. (cell(1) == 0 .AND. cell(2) == 0 .AND. cell(3) == 0))
THEN
915 IF (my_real_mat_real_space)
THEN
916 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, coskl)
917 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
922 CALL dbcsr_add_local(real_mat_kp(ik)%matrix, mat_work, 1.0_dp, -sinkl)
923 CALL dbcsr_add_local(imag_mat_kp(ik)%matrix, mat_work, 1.0_dp, -coskl)
932 CALL dbcsr_filter(real_mat_kp(ik)%matrix, eps_filter_im_time)
933 CALL dbcsr_filter(imag_mat_kp(ik)%matrix, eps_filter_im_time)
939 CALL timestop(handle)
950 SUBROUTINE dbcsr_add_local(mat_a, mat_b, alpha, beta)
951 TYPE(
dbcsr_type),
INTENT(INOUT) :: mat_a, mat_b
952 REAL(kind=
dp),
INTENT(IN) :: alpha, beta
956 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: block_to_compute, data_block
963 NULLIFY (block_to_compute)
965 row=row, col=col, block=block_to_compute, found=found)
969 block_to_compute(:, :) = alpha*block_to_compute(:, :) + beta*data_block(:, :)
974 END SUBROUTINE dbcsr_add_local
995 SUBROUTINE compute_wc_real_space_tau_gw(fm_mat_W_tau, cfm_mat_Q, fm_mat_L_re, fm_mat_L_im, &
996 dimen_RI, num_integ_points, jquad, &
997 ikp, tj, tau_tj, weights_cos_tf_w_to_t, ikp_local, &
998 para_env, kpoints, qs_env, wkp_W)
1000 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: fm_mat_w_tau
1002 TYPE(
cp_fm_type),
INTENT(IN) :: fm_mat_l_re, fm_mat_l_im
1003 INTEGER,
INTENT(IN) :: dimen_ri, num_integ_points, jquad, ikp
1004 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: tj
1005 REAL(kind=
dp),
DIMENSION(num_integ_points), &
1006 INTENT(IN) :: tau_tj
1007 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN) :: weights_cos_tf_w_to_t
1008 INTEGER,
INTENT(IN) :: ikp_local
1012 REAL(kind=
dp),
DIMENSION(:),
INTENT(IN) :: wkp_w
1014 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_Wc_real_space_tau_GW'
1016 INTEGER :: handle, handle2, i_global, iatom, iatom_old, iib, iquad, irow, j_global, jatom, &
1017 jatom_old, jcol, jjb, jkp, ncol_local, nkp, nrow_local, num_cells
1018 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: atom_from_ri_index
1019 INTEGER,
DIMENSION(:),
POINTER :: col_indices, row_indices
1020 INTEGER,
DIMENSION(:, :),
POINTER :: index_to_cell
1021 REAL(kind=
dp) :: contribution, omega, tau, weight, &
1022 weight_im, weight_re
1023 REAL(kind=
dp),
DIMENSION(3, 3) :: hmat
1024 REAL(kind=
dp),
DIMENSION(:),
POINTER :: wkp
1025 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: xkp
1027 TYPE(
cp_cfm_type) :: cfm_mat_l, cfm_mat_work, cfm_mat_work_2
1028 TYPE(
cp_fm_type) :: fm_dummy, fm_mat_work_global, &
1032 CALL timeset(routinen, handle)
1034 CALL timeset(routinen//
"_1", handle2)
1049 CALL cp_fm_create(fm_mat_work_global, fm_mat_w_tau(1)%matrix_struct)
1052 CALL cp_fm_create(fm_mat_work_local, cfm_mat_q%matrix_struct)
1055 CALL timestop(handle2)
1057 CALL timeset(routinen//
"_2", handle2)
1067 nrow_local=nrow_local, &
1068 ncol_local=ncol_local, &
1069 row_indices=row_indices, &
1070 col_indices=col_indices)
1072 DO jjb = 1, ncol_local
1073 j_global = col_indices(jjb)
1074 DO iib = 1, nrow_local
1075 i_global = row_indices(iib)
1076 IF (j_global == i_global .AND. i_global <= dimen_ri)
THEN
1077 cfm_mat_q%local_data(iib, jjb) = cfm_mat_q%local_data(iib, jjb) -
z_one
1082 CALL timestop(handle2)
1084 CALL timeset(routinen//
"_3", handle2)
1087 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_q, cfm_mat_l, &
1091 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri,
z_one, cfm_mat_l, cfm_mat_work, &
1094 CALL timestop(handle2)
1096 CALL timeset(routinen//
"_4", handle2)
1099 index_to_cell => kpoints%index_to_cell
1100 num_cells =
SIZE(index_to_cell, 2)
1104 ALLOCATE (atom_from_ri_index(dimen_ri))
1108 NULLIFY (cell, particle_set)
1109 CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set)
1115 nrow_local=nrow_local, &
1116 ncol_local=ncol_local, &
1117 row_indices=row_indices, &
1118 col_indices=col_indices)
1120 DO irow = 1, nrow_local
1121 DO jcol = 1, ncol_local
1123 iatom = atom_from_ri_index(row_indices(irow))
1124 jatom = atom_from_ri_index(col_indices(jcol))
1126 IF (iatom .NE. iatom_old .OR. jatom .NE. jatom_old)
THEN
1130 num_cells, iatom, jatom, xkp(1:3, ikp), wkp_w(ikp), &
1131 cell, index_to_cell, hmat, particle_set)
1138 contribution = weight_re*real(cfm_mat_work_2%local_data(irow, jcol)) + &
1139 weight_im*aimag(cfm_mat_work_2%local_data(irow, jcol))
1141 fm_mat_work_local%local_data(irow, jcol) = fm_mat_work_local%local_data(irow, jcol) + contribution
1146 CALL timestop(handle2)
1148 CALL timeset(routinen//
"_5", handle2)
1150 IF (ikp_local == -1)
THEN
1154 DO iquad = 1, num_integ_points
1158 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1160 IF (jquad == 1 .AND. ikp == 1)
THEN
1161 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1164 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_w_tau(iquad), beta=weight, matrix_b=fm_mat_work_global)
1172 CALL para_env%sync()
1174 IF (ikp_local == jkp)
THEN
1180 CALL para_env%sync()
1182 DO iquad = 1, num_integ_points
1186 weight = weights_cos_tf_w_to_t(iquad, jquad)*cos(tau*omega)
1188 IF (jquad == 1 .AND. jkp == 1)
THEN
1189 CALL cp_fm_set_all(matrix=fm_mat_w_tau(iquad), alpha=0.0_dp)
1193 matrix_b=fm_mat_work_global)
1207 DEALLOCATE (atom_from_ri_index)
1209 CALL timestop(handle2)
1211 CALL timestop(handle)
1213 END SUBROUTINE compute_wc_real_space_tau_gw
1223 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_matrix_minv
1227 INTEGER :: dimen_ri, num_integ_points
1229 CHARACTER(LEN=*),
PARAMETER :: routinen =
'Wc_to_Minv_Wc_Minv'
1231 INTEGER :: handle, jquad
1232 TYPE(
cp_fm_type) :: fm_work_minv, fm_work_minv_w
1234 CALL timeset(routinen, handle)
1236 CALL cp_fm_create(fm_work_minv, fm_mat_w(1)%matrix_struct)
1239 CALL cp_fm_create(fm_work_minv_w, fm_mat_w(1)%matrix_struct)
1241 DO jquad = 1, num_integ_points
1243 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv, fm_mat_w(jquad), &
1244 0.0_dp, fm_work_minv_w)
1245 CALL parallel_gemm(
'N',
'N', dimen_ri, dimen_ri, dimen_ri, 1.0_dp, fm_work_minv_w, fm_work_minv, &
1246 0.0_dp, fm_mat_w(jquad))
1254 CALL timestop(handle)
1256 END SUBROUTINE wc_to_minv_wc_minv
1270 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:), &
1271 INTENT(OUT) :: wkp_w, wkp_v
1273 REAL(kind=
dp),
DIMENSION(3, 3) :: h_inv
1274 INTEGER,
DIMENSION(3) :: periodic
1276 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_wkp_W'
1278 INTEGER :: handle, i_x, ikp, info, j_y, k_z, &
1279 kpoint_weights_w_method, n_x, n_y, &
1280 n_z, nkp, nsuperfine, num_lin_eqs
1281 REAL(kind=
dp) :: exp_kpoints, integral, k_sq, weight
1282 REAL(kind=
dp),
DIMENSION(3) :: k_vec, x_vec
1283 REAL(kind=
dp),
DIMENSION(:),
POINTER :: right_side, wkp, wkp_tmp
1284 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: matrix_lin_eqs, xkp
1286 CALL timeset(routinen, handle)
1288 kpoint_weights_w_method = qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method
1302 ALLOCATE (wkp_v(nkp), wkp_w(nkp))
1306 IF (
ALLOCATED(qs_env%mp2_env%ri_rpa_im_time%wkp_V))
THEN
1307 wkp_v(:) = qs_env%mp2_env%ri_rpa_im_time%wkp_V(:)
1321 exp_kpoints = qs_env%mp2_env%ri_rpa_im_time%exp_tailored_weights
1324 IF (sum(periodic) == 2) exp_kpoints = -1.0_dp
1331 IF (periodic(1) == 1)
THEN
1336 IF (periodic(2) == 1)
THEN
1341 IF (periodic(3) == 1)
THEN
1349 weight = 1.0_dp/(real(n_x,
dp)*real(n_y,
dp)*real(n_z,
dp))
1354 IF (periodic(1) == 1)
THEN
1355 x_vec(1) = (real(i_x - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
1359 IF (periodic(2) == 1)
THEN
1360 x_vec(2) = (real(j_y - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
1364 IF (periodic(3) == 1)
THEN
1365 x_vec(3) = (real(k_z - nsuperfine/2,
dp) - 0.5_dp)/real(nsuperfine,
dp)
1370 k_vec = matmul(h_inv(1:3, 1:3), x_vec)
1371 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1372 integral = integral + weight*k_sq**(exp_kpoints*0.5_dp)
1378 num_lin_eqs = nkp + 2
1380 ALLOCATE (matrix_lin_eqs(num_lin_eqs, num_lin_eqs))
1381 matrix_lin_eqs(:, :) = 0.0_dp
1385 k_vec = matmul(h_inv(1:3, 1:3), xkp(1:3, ikp))
1386 k_sq = k_vec(1)**2 + k_vec(2)**2 + k_vec(3)**2
1388 matrix_lin_eqs(ikp, ikp) = 2.0_dp
1389 matrix_lin_eqs(ikp, nkp + 1) = 1.0_dp
1390 matrix_lin_eqs(nkp + 1, ikp) = 1.0_dp
1392 matrix_lin_eqs(ikp, nkp + 2) = k_sq**(exp_kpoints*0.5_dp)
1393 matrix_lin_eqs(nkp + 2, ikp) = k_sq**(exp_kpoints*0.5_dp)
1397 CALL invmat(matrix_lin_eqs, info)
1401 ALLOCATE (right_side(num_lin_eqs))
1403 right_side(nkp + 1) = 1.0_dp
1405 right_side(nkp + 2) = integral
1407 ALLOCATE (wkp_tmp(num_lin_eqs))
1409 wkp_tmp(1:num_lin_eqs) = matmul(matrix_lin_eqs, right_side)
1411 wkp_w(1:nkp) = wkp_tmp(1:nkp)
1413 DEALLOCATE (matrix_lin_eqs, right_side, wkp_tmp)
1417 CALL timestop(handle)
1428 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: eigenval_kp
1430 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_bandstruc_and_k_dependent_MOs'
1432 INTEGER :: handle, ikp, ispin, nmo, nspins
1433 INTEGER,
DIMENSION(3) :: nkp_grid_g
1434 REAL(kind=
dp),
DIMENSION(:),
POINTER :: ev
1435 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: kpgeneral
1439 CALL timeset(routinen, handle)
1441 NULLIFY (qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1442 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1443 qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1446 nkp_grid_g(1:3) = (/1, 1, 1/)
1448 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
1450 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_G, &
1451 "MONKHORST-PACK", para_env%num_pe, &
1452 mp_grid=nkp_grid_g(1:3))
1454 IF (qs_env%mp2_env%ri_g0w0%do_kpoints_Sigma)
THEN
1457 CALL get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1459 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma, &
1460 "GENERAL", para_env%num_pe, &
1461 kpgeneral=kpgeneral)
1463 CALL create_kp_and_calc_kp_orbitals(qs_env, qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma_no_xc, &
1464 "GENERAL", para_env%num_pe, &
1465 kpgeneral=kpgeneral, with_xc_terms=.false.)
1467 kpoints_sigma => qs_env%mp2_env%ri_rpa_im_time%kpoints_Sigma
1468 nmo =
SIZE(eigenval_kp, 1)
1469 nspins =
SIZE(eigenval_kp, 3)
1471 ALLOCATE (qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(nmo))
1472 qs_env%mp2_env%ri_rpa_im_time%Eigenval_Gamma(:) = eigenval_kp(:, 1, 1)
1474 DEALLOCATE (eigenval_kp)
1476 ALLOCATE (eigenval_kp(nmo, kpoints_sigma%nkp, nspins))
1478 DO ikp = 1, kpoints_sigma%nkp
1480 DO ispin = 1, nspins
1482 ev => kpoints_sigma%kp_env(ikp)%kpoint_env%mos(1, ispin)%eigenvalues
1484 eigenval_kp(:, ikp, ispin) = ev(:)
1490 DEALLOCATE (kpgeneral)
1494 CALL release_hfx_stuff(qs_env)
1496 CALL timestop(handle)
1504 SUBROUTINE release_hfx_stuff(qs_env)
1507 IF (
ASSOCIATED(qs_env%x_data) .AND. .NOT. qs_env%mp2_env%ri_g0w0%do_ri_Sigma_x)
THEN
1511 END SUBROUTINE release_hfx_stuff
1523 SUBROUTINE create_kp_and_calc_kp_orbitals(qs_env, kpoints, scheme, &
1524 group_size_ext, mp_grid, kpgeneral, with_xc_terms)
1528 CHARACTER(LEN=*),
INTENT(IN) :: scheme
1529 INTEGER :: group_size_ext
1530 INTEGER,
DIMENSION(3),
INTENT(IN),
OPTIONAL :: mp_grid
1531 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(IN), &
1532 OPTIONAL :: kpgeneral
1533 LOGICAL,
OPTIONAL :: with_xc_terms
1535 CHARACTER(LEN=*),
PARAMETER :: routinen =
'create_kp_and_calc_kp_orbitals'
1537 INTEGER :: handle, i_dim, i_re_im, ikp, ispin, nkp, &
1539 INTEGER,
DIMENSION(3) :: cell_grid, periodic
1540 LOGICAL :: my_with_xc_terms
1541 REAL(kind=
dp),
DIMENSION(:),
POINTER :: eigenvalues
1548 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s, matrix_s_desymm
1549 TYPE(
dbcsr_p_type),
DIMENSION(:, :),
POINTER :: mat_ks_kp, mat_s_kp
1556 CALL timeset(routinen, handle)
1558 my_with_xc_terms = .true.
1559 IF (
PRESENT(with_xc_terms)) my_with_xc_terms = with_xc_terms
1562 para_env=para_env, &
1563 blacs_env=blacs_env, &
1564 matrix_s=matrix_s, &
1566 scf_control=scf_control, &
1571 group_size_ext=group_size_ext)
1580 CALL get_cell(cell=cell, periodic=periodic)
1585 IF (periodic(i_dim) == 1)
THEN
1586 cell_grid(i_dim) = max(min((kpoints%nkp_grid(i_dim)/2)*2 - 1, 1), 3)
1588 cell_grid(i_dim) = 1
1594 CALL get_qs_env(qs_env, matrix_s=matrix_s, scf_env=scf_env, scf_control=scf_control, dft_control=dft_control)
1596 NULLIFY (matrix_s_desymm)
1598 ALLOCATE (matrix_s_desymm(1)%matrix)
1599 CALL dbcsr_create(matrix=matrix_s_desymm(1)%matrix, template=matrix_s(1)%matrix, &
1600 matrix_type=dbcsr_type_no_symmetry)
1607 matrix_struct => kpoints%kp_env(1)%kpoint_env%wmat(1, 1)%matrix_struct
1615 nspins = dft_control%nspins
1617 DO ispin = 1, nspins
1620 IF (my_with_xc_terms)
THEN
1621 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_ks(ispin)%matrix, kpoints, ispin)
1623 CALL mat_kp_from_mat_gamma(qs_env, mat_ks_kp, qs_env%mp2_env%ri_g0w0%matrix_sigma_x_minus_vxc(ispin)%matrix, &
1629 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 1)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(1, ispin))
1632 CALL copy_dbcsr_to_fm(mat_ks_kp(ikp, 2)%matrix, kpoints%kp_env(ikp)%kpoint_env%wmat(2, ispin))
1641 kp => kpoints%kp_env(ikp)%kpoint_env
1643 CALL get_mo_set(kp%mos(1, ispin), mo_coeff=rmos, eigenvalues=eigenvalues)
1644 CALL get_mo_set(kp%mos(2, ispin), mo_coeff=imos)
1647 qs_env%mp2_env%ri_rpa_im_time%make_overlap_mat_ao_pos_definite)
THEN
1648 CALL cp_cfm_geeig_canon(cksmat, csmat, cmos, eigenvalues, cwork, scf_control%eps_eigval)
1650 CALL cp_cfm_geeig(cksmat, csmat, cmos, eigenvalues, cwork)
1655 kp%mos(2, ispin)%eigenvalues = eigenvalues
1666 DEALLOCATE (mat_ks_kp)
1673 DEALLOCATE (mat_s_kp)
1676 DEALLOCATE (matrix_s_desymm)
1684 CALL timestop(handle)
1686 END SUBROUTINE create_kp_and_calc_kp_orbitals
1704 LOGICAL,
INTENT(IN),
OPTIONAL :: real_mat_real_space
1706 CHARACTER(LEN=*),
PARAMETER :: routinen =
'mat_kp_from_mat_gamma'
1708 INTEGER :: handle, i_cell, i_re_im, ikp, nkp, &
1710 INTEGER,
DIMENSION(3) :: periodic
1711 INTEGER,
DIMENSION(:, :, :),
POINTER :: cell_to_index
1712 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: xkp
1714 TYPE(
dbcsr_p_type),
DIMENSION(:),
POINTER :: mat_real_space
1716 CALL timeset(routinen, handle)
1719 CALL get_cell(cell=cell, periodic=periodic)
1720 num_cells = 3**(periodic(1) + periodic(2) + periodic(3))
1722 NULLIFY (mat_real_space)
1724 DO i_cell = 1, num_cells
1725 ALLOCATE (mat_real_space(i_cell)%matrix)
1726 CALL dbcsr_create(matrix=mat_real_space(i_cell)%matrix, &
1729 CALL dbcsr_set(mat_real_space(i_cell)%matrix, 0.0_dp)
1732 CALL dbcsr_copy(mat_real_space(1)%matrix, mat_gamma)
1736 NULLIFY (xkp, cell_to_index)
1737 CALL get_kpoint_info(kpoints, nkp=nkp, xkp=xkp, cell_to_index=cell_to_index)
1739 IF (ispin == 1)
THEN
1744 ALLOCATE (mat_kp(ikp, i_re_im)%matrix)
1745 CALL dbcsr_create(matrix=mat_kp(ikp, i_re_im)%matrix, template=mat_gamma)
1747 CALL dbcsr_set(mat_kp(ikp, i_re_im)%matrix, 0.0_dp)
1752 IF (
PRESENT(real_mat_real_space))
THEN
1754 real_mat_real_space)
1759 DO i_cell = 1, num_cells
1762 DEALLOCATE (mat_real_space)
1764 CALL timestop(handle)
1773 SUBROUTINE get_kpgeneral_for_sigma_kpoints(qs_env, kpgeneral)
1775 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: kpgeneral
1777 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_kpgeneral_for_Sigma_kpoints'
1779 INTEGER :: handle, i_kp_in_kp_line, i_special_kp, &
1780 i_x, ikk, j_y, k_z, n_kp_in_kp_line, &
1782 INTEGER,
DIMENSION(:),
POINTER :: nkp_grid
1784 CALL timeset(routinen, handle)
1786 n_special_kp = qs_env%mp2_env%ri_g0w0%n_special_kp
1787 n_kp_in_kp_line = qs_env%mp2_env%ri_g0w0%n_kp_in_kp_line
1788 IF (n_special_kp > 0)
THEN
1789 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = n_kp_in_kp_line*(n_special_kp - 1) + 1
1791 qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp = 0
1794 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack = qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(1)* &
1795 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(2)* &
1796 qs_env%mp2_env%ri_g0w0%kp_grid_Sigma(3)
1798 qs_env%mp2_env%ri_g0w0%nkp_self_energy = qs_env%mp2_env%ri_g0w0%nkp_self_energy_special_kp + &
1799 qs_env%mp2_env%ri_g0w0%nkp_self_energy_monkh_pack
1801 ALLOCATE (kpgeneral(3, qs_env%mp2_env%ri_g0w0%nkp_self_energy))
1803 IF (n_special_kp > 0)
THEN
1805 kpgeneral(1:3, 1) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, 1)
1809 DO i_special_kp = 2, n_special_kp
1810 DO i_kp_in_kp_line = 1, n_kp_in_kp_line
1813 kpgeneral(1:3, ikk) = qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1) + &
1814 REAL(i_kp_in_kp_line, kind=
dp)/real(n_kp_in_kp_line, kind=
dp)* &
1815 (qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp) - &
1816 qs_env%mp2_env%ri_g0w0%xkp_special_kp(1:3, i_special_kp - 1))
1827 nkp_grid => qs_env%mp2_env%ri_g0w0%kp_grid_Sigma
1829 DO i_x = 1, nkp_grid(1)
1830 DO j_y = 1, nkp_grid(2)
1831 DO k_z = 1, nkp_grid(3)
1833 kpgeneral(1, ikk) = real(2*i_x - nkp_grid(1) - 1, kind=
dp)/(2._dp*real(nkp_grid(1), kind=
dp))
1834 kpgeneral(2, ikk) = real(2*j_y - nkp_grid(2) - 1, kind=
dp)/(2._dp*real(nkp_grid(2), kind=
dp))
1835 kpgeneral(3, ikk) = real(2*k_z - nkp_grid(3) - 1, kind=
dp)/(2._dp*real(nkp_grid(3), kind=
dp))
1840 CALL timestop(handle)
1842 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, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Definition and initialisation of the mo data type.
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