117#include "./base/base_uses.f90"
123 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'rpa_main'
179 SUBROUTINE rpa_ri_compute_en(qs_env, Erpa, mp2_env, BIb_C, BIb_C_gw, BIb_C_bse_ij, BIb_C_bse_ab, &
180 para_env, para_env_sub, color_sub, &
181 gd_array, gd_B_virtual, gd_B_all, gd_B_occ_bse, gd_B_virt_bse, &
182 mo_coeff, fm_matrix_PQ, fm_matrix_L_kpoints, fm_matrix_Minv_L_kpoints, &
183 fm_matrix_Minv, fm_matrix_Minv_Vtrunc_Minv, kpoints, &
184 Eigenval, nmo, homo, dimen_RI, dimen_RI_red, gw_corr_lev_occ, gw_corr_lev_virt, &
186 unit_nr, do_ri_sos_laplace_mp2, my_do_gw, do_im_time, do_bse, matrix_s, &
187 mat_munu, mat_P_global, t_3c_M, t_3c_O, t_3c_O_compressed, t_3c_O_ind, &
188 starts_array_mc, ends_array_mc, &
189 starts_array_mc_block, ends_array_mc_block, calc_forces)
192 REAL(kind=
dp),
INTENT(OUT) :: erpa
193 TYPE(
mp2_type),
INTENT(INOUT) :: mp2_env
195 INTENT(INOUT) :: bib_c, bib_c_gw
196 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
197 INTENT(INOUT) :: bib_c_bse_ij, bib_c_bse_ab
199 INTEGER,
INTENT(INOUT) :: color_sub
202 INTENT(INOUT) :: gd_b_virtual
204 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(IN) :: mo_coeff
206 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_matrix_l_kpoints, &
207 fm_matrix_minv_l_kpoints, &
209 fm_matrix_minv_vtrunc_minv
211 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
212 INTENT(INOUT) :: eigenval
213 INTEGER,
INTENT(IN) :: nmo
214 INTEGER,
DIMENSION(:),
INTENT(IN) :: homo
215 INTEGER,
INTENT(IN) :: dimen_ri, dimen_ri_red
216 INTEGER,
DIMENSION(:),
INTENT(IN) :: gw_corr_lev_occ, gw_corr_lev_virt
217 INTEGER,
INTENT(IN) :: bse_lev_virt, unit_nr
218 LOGICAL,
INTENT(IN) :: do_ri_sos_laplace_mp2, my_do_gw, &
223 TYPE(dbt_type) :: t_3c_m
224 TYPE(dbt_type),
ALLOCATABLE,
DIMENSION(:, :) :: t_3c_o
226 DIMENSION(:, :, :),
INTENT(INOUT) :: t_3c_o_compressed
228 DIMENSION(:, :, :),
INTENT(INOUT) :: t_3c_o_ind
229 INTEGER,
ALLOCATABLE,
DIMENSION(:),
INTENT(IN) :: starts_array_mc, ends_array_mc, &
230 starts_array_mc_block, &
232 LOGICAL,
INTENT(IN) :: calc_forces
234 CHARACTER(LEN=*),
PARAMETER :: routinen =
'rpa_ri_compute_en'
236 INTEGER :: best_integ_group_size, best_num_integ_point, color_rpa_group, dimen_homo_square, &
237 dimen_nm_gw, dimen_virt_square, handle, handle2, handle3, ierr, iib, &
238 input_num_integ_groups, integ_group_size, ispin, jjb, min_integ_group_size, &
239 my_ab_comb_bse_end, my_ab_comb_bse_size, my_ab_comb_bse_start, my_group_l_end, &
240 my_group_l_size, my_group_l_start, my_ij_comb_bse_end, my_ij_comb_bse_size, &
241 my_ij_comb_bse_start, my_nm_gw_end, my_nm_gw_size, my_nm_gw_start, ncol_block_mat, &
242 ngroup, nrow_block_mat, nspins, num_integ_group, num_integ_points, pos_integ_group
243 INTEGER(KIND=int_8) :: mem
244 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: dimen_ia, my_ia_end, my_ia_size, &
246 LOGICAL :: do_kpoints_from_gamma, do_minimax_quad, &
247 my_open_shell, skip_integ_group_opt
248 REAL(kind=
dp) :: allowed_memory, avail_mem, e_range, emax, emin, mem_for_iak, mem_for_qk, &
249 mem_min, mem_per_group, mem_per_rank, mem_per_repl, mem_real
250 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: eigenval_kp
251 TYPE(
cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: fm_mat_q, fm_mat_q_gemm, fm_mat_s, &
253 TYPE(
cp_fm_type),
DIMENSION(1) :: fm_mat_r_gw, fm_mat_s_ab_bse, &
257 DIMENSION(:) :: bib_c_2d, bib_c_2d_gw
260 CALL timeset(routinen, handle)
271 IF (mp2_env%ri_rpa%do_rse)
THEN
287 my_open_shell = (nspins == 2)
288 ALLOCATE (virtual(nspins), dimen_ia(nspins), my_ia_end(nspins), my_ia_start(nspins), my_ia_size(nspins))
289 virtual(:) = nmo - homo(:)
290 dimen_ia(:) = virtual(:)*homo(:)
292 ALLOCATE (eigenval_kp(nmo, 1, nspins))
293 eigenval_kp(:, 1, :) = eigenval(:, :)
295 IF (do_im_time) mp2_env%ri_rpa%minimax_quad = .true.
296 do_minimax_quad = mp2_env%ri_rpa%minimax_quad
298 IF (do_ri_sos_laplace_mp2)
THEN
299 num_integ_points = mp2_env%ri_laplace%n_quadrature
300 input_num_integ_groups = mp2_env%ri_laplace%num_integ_groups
303 e_range = mp2_env%e_range
304 IF (mp2_env%e_range <= 1.0_dp .OR. mp2_env%e_gap <= 0.0_dp)
THEN
308 IF (homo(ispin) > 0)
THEN
309 emin = min(emin, 2.0_dp*(eigenval(homo(ispin) + 1, ispin) - eigenval(homo(ispin), ispin)))
310 emax = max(emax, 2.0_dp*(maxval(eigenval(:, ispin)) - minval(eigenval(:, ispin))))
315 IF (e_range < 2.0_dp) e_range = 2.0_dp
319 jjb = num_integ_points - 1
321 num_integ_points = num_integ_points - 1
327 cpassert(num_integ_points >= 1)
329 num_integ_points = mp2_env%ri_rpa%rpa_num_quad_points
330 input_num_integ_groups = mp2_env%ri_rpa%rpa_num_integ_groups
331 IF (my_do_gw .AND. do_minimax_quad)
THEN
332 IF (num_integ_points > 34)
THEN
334 CALL cp_warn(__location__, &
335 "The required number of quadrature point exceeds the maximum possible in the "// &
336 "Minimax quadrature scheme. The number of quadrature point has been reset to 30.")
337 num_integ_points = 30
340 IF (do_minimax_quad .AND. num_integ_points > 20)
THEN
342 CALL cp_warn(__location__, &
343 "The required number of quadrature point exceeds the maximum possible in the "// &
344 "Minimax quadrature scheme. The number of quadrature point has been reset to 20.")
345 num_integ_points = 20
349 allowed_memory = mp2_env%mp2_memory
351 CALL get_group_dist(gd_array, color_sub, my_group_l_start, my_group_l_end, my_group_l_size)
353 ngroup = para_env%num_pe/para_env_sub%num_pe
356 IF (do_im_time .OR. mp2_env%ri_g0w0%do_periodic .OR. do_bse)
THEN
358 integ_group_size = ngroup
359 best_num_integ_point = num_integ_points
365 mem_for_iak = real(sum(dimen_ia), kind=
dp)*dimen_ri_red*8.0_dp/(1024_dp**2)
366 mem_for_qk = real(dimen_ri_red, kind=
dp)*nspins*dimen_ri_red*8.0_dp/(1024_dp**2)
369 mem_real = (mem + 1024*1024 - 1)/(1024*1024)
370 CALL para_env%min(mem_real)
372 mem_per_rank = 0.0_dp
375 mem_per_repl = mem_for_iak
377 mem_per_repl = mem_per_repl + 2.0_dp*mem_for_qk
379 IF (calc_forces)
CALL rpa_grad_needed_mem(homo, virtual, dimen_ri_red, mem_per_rank, mem_per_repl, do_ri_sos_laplace_mp2)
382 mem_min = mem_per_repl/para_env%num_pe + mem_per_rank
384 IF (unit_nr > 0)
THEN
385 WRITE (unit_nr,
'(T3,A,T68,F9.2,A4)')
'RI_INFO| Minimum required memory per MPI process:', mem_min,
' MiB'
386 WRITE (unit_nr,
'(T3,A,T68,F9.2,A4)')
'RI_INFO| Available memory per MPI process:', mem_real,
' MiB'
390 mem_real = min(mem_real, allowed_memory)
392 mem_real = mem_real - mem_per_rank
394 mem_per_group = mem_real*para_env_sub%num_pe
397 skip_integ_group_opt = .false.
400 IF (input_num_integ_groups > 0)
THEN
401 IF (mod(num_integ_points, input_num_integ_groups) == 0)
THEN
402 IF (mod(ngroup, input_num_integ_groups) == 0)
THEN
403 best_integ_group_size = ngroup/input_num_integ_groups
404 best_num_integ_point = num_integ_points/input_num_integ_groups
405 skip_integ_group_opt = .true.
407 IF (unit_nr > 0)
WRITE (unit_nr,
'(T3,A)')
'Total number of groups not multiple of NUM_INTEG_GROUPS'
410 IF (unit_nr > 0)
WRITE (unit_nr,
'(T3,A)')
'NUM_QUAD_POINTS not multiple of NUM_INTEG_GROUPS'
414 IF (.NOT. skip_integ_group_opt)
THEN
415 best_integ_group_size = ngroup
416 best_num_integ_point = num_integ_points
418 min_integ_group_size = max(1, ngroup/num_integ_points)
420 integ_group_size = min_integ_group_size - 1
421 DO iib = min_integ_group_size + 1, ngroup
422 integ_group_size = integ_group_size + 1
425 IF (mod(ngroup, integ_group_size) /= 0) cycle
428 avail_mem = integ_group_size*mem_per_group
429 IF (avail_mem < mem_per_repl) cycle
432 num_integ_group = ngroup/integ_group_size
433 IF (mod(num_integ_points, num_integ_group) /= 0) cycle
435 best_num_integ_point = num_integ_points/num_integ_group
436 best_integ_group_size = integ_group_size
443 integ_group_size = best_integ_group_size
447 IF (unit_nr > 0 .AND. .NOT. do_im_time)
THEN
448 IF (do_ri_sos_laplace_mp2)
THEN
449 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
450 "RI_INFO| Group size for laplace numerical integration:", integ_group_size*para_env_sub%num_pe
451 WRITE (unit=unit_nr, fmt=
"(T3,A)") &
452 "INTEG_INFO| MINIMAX approximation"
453 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
454 "INTEG_INFO| Number of integration points:", num_integ_points
455 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
456 "INTEG_INFO| Number of integration points per Laplace group:", best_num_integ_point
458 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
459 "RI_INFO| Group size for frequency integration:", integ_group_size*para_env_sub%num_pe
460 IF (do_minimax_quad)
THEN
461 WRITE (unit=unit_nr, fmt=
"(T3,A)") &
462 "INTEG_INFO| MINIMAX quadrature"
464 WRITE (unit=unit_nr, fmt=
"(T3,A)") &
465 "INTEG_INFO| Clenshaw-Curtius quadrature"
467 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
468 "INTEG_INFO| Number of integration points:", num_integ_points
469 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
470 "INTEG_INFO| Number of integration points per RPA group:", best_num_integ_point
475 num_integ_group = ngroup/integ_group_size
477 pos_integ_group = mod(color_sub, integ_group_size)
478 color_rpa_group = color_sub/integ_group_size
480 CALL timeset(routinen//
"_reorder", handle2)
484 ALLOCATE (bib_c_2d(nspins))
486 IF (.NOT. do_im_time)
THEN
491 CALL calculate_bib_c_2d(bib_c_2d(ispin)%array, bib_c(ispin)%array, para_env_sub, dimen_ia(ispin), &
492 homo(ispin), virtual(ispin), gd_b_virtual(ispin), &
493 my_ia_size(ispin), my_ia_start(ispin), my_ia_end(ispin), my_group_l_size)
495 DEALLOCATE (bib_c(ispin)%array)
502 ALLOCATE (bib_c_2d_gw(nspins))
504 CALL timeset(routinen//
"_reorder_gw", handle3)
506 dimen_nm_gw = nmo*(gw_corr_lev_occ(1) + gw_corr_lev_virt(1))
510 CALL calculate_bib_c_2d(bib_c_2d_gw(ispin)%array, bib_c_gw(ispin)%array, para_env_sub, dimen_nm_gw, &
511 gw_corr_lev_occ(ispin) + gw_corr_lev_virt(ispin), nmo, gd_b_all, &
512 my_nm_gw_size, my_nm_gw_start, my_nm_gw_end, my_group_l_size)
513 DEALLOCATE (bib_c_gw(ispin)%array)
518 CALL timestop(handle3)
525 CALL timeset(routinen//
"_reorder_bse1", handle3)
527 dimen_homo_square = homo(1)**2
530 CALL calculate_bib_c_2d(bib_c_2d_bse_ij(1)%array, bib_c_bse_ij, para_env_sub, dimen_homo_square, &
531 homo(1), homo(1), gd_b_occ_bse, &
532 my_ij_comb_bse_size, my_ij_comb_bse_start, my_ij_comb_bse_end, my_group_l_size)
534 DEALLOCATE (bib_c_bse_ij)
537 CALL timestop(handle3)
539 CALL timeset(routinen//
"_reorder_bse2", handle3)
541 dimen_virt_square = bse_lev_virt**2
543 CALL calculate_bib_c_2d(bib_c_2d_bse_ab(1)%array, bib_c_bse_ab, para_env_sub, dimen_virt_square, &
544 bse_lev_virt, bse_lev_virt, gd_b_virt_bse, &
545 my_ab_comb_bse_size, my_ab_comb_bse_start, my_ab_comb_bse_end, my_group_l_size)
547 DEALLOCATE (bib_c_bse_ab)
550 CALL timestop(handle3)
554 CALL timestop(handle2)
556 IF (num_integ_group > 1)
THEN
557 ALLOCATE (para_env_rpa)
558 CALL para_env_rpa%from_split(para_env, color_rpa_group)
560 para_env_rpa => para_env
567 ALLOCATE (fm_mat_q(nspins), fm_mat_q_gemm(1), fm_mat_s(1))
569 ALLOCATE (fm_mat_q(nspins), fm_mat_q_gemm(nspins), fm_mat_s(nspins))
572 CALL create_integ_mat(bib_c_2d, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
573 dimen_ri_red, dimen_ia, color_rpa_group, &
574 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
575 my_ia_size, my_ia_start, my_ia_end, &
576 my_group_l_size, my_group_l_start, my_group_l_end, &
577 para_env_rpa, fm_mat_s, nrow_block_mat, ncol_block_mat, &
578 dimen_ia_for_block_size=dimen_ia(1), &
579 do_im_time=do_im_time, fm_mat_q_gemm=fm_mat_q_gemm, fm_mat_q=fm_mat_q, qs_env=qs_env)
581 DEALLOCATE (bib_c_2d, my_ia_end, my_ia_size, my_ia_start)
584 ALLOCATE (fm_mat_s_gw(nspins))
585 IF (my_do_gw .AND. .NOT. do_im_time)
THEN
587 CALL create_integ_mat(bib_c_2d_gw, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
588 dimen_ri_red, [dimen_nm_gw, dimen_nm_gw], color_rpa_group, &
589 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
590 [my_nm_gw_size, my_nm_gw_size], [my_nm_gw_start, my_nm_gw_start], [my_nm_gw_end, my_nm_gw_end], &
591 my_group_l_size, my_group_l_start, my_group_l_end, &
592 para_env_rpa, fm_mat_s_gw, nrow_block_mat, ncol_block_mat, &
593 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context, &
594 fm_mat_q=fm_mat_r_gw)
595 DEALLOCATE (bib_c_2d_gw)
601 CALL create_integ_mat(bib_c_2d_bse_ij, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
602 dimen_ri_red, [dimen_homo_square], color_rpa_group, &
603 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
604 [my_ij_comb_bse_size], [my_ij_comb_bse_start], [my_ij_comb_bse_end], &
605 my_group_l_size, my_group_l_start, my_group_l_end, &
606 para_env_rpa, fm_mat_s_ij_bse, nrow_block_mat, ncol_block_mat, &
607 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context)
609 CALL create_integ_mat(bib_c_2d_bse_ab, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
610 dimen_ri_red, [dimen_virt_square], color_rpa_group, &
611 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
612 [my_ab_comb_bse_size], [my_ab_comb_bse_start], [my_ab_comb_bse_end], &
613 my_group_l_size, my_group_l_start, my_group_l_end, &
614 para_env_rpa, fm_mat_s_ab_bse, nrow_block_mat, ncol_block_mat, &
615 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context)
619 do_kpoints_from_gamma = qs_env%mp2_env%ri_rpa_im_time%do_kpoints_from_Gamma
620 IF (do_kpoints_from_gamma)
THEN
626 CALL rpa_num_int(qs_env, erpa, mp2_env, para_env, para_env_rpa, para_env_sub, unit_nr, &
627 homo, virtual, dimen_ri, dimen_ri_red, dimen_ia, dimen_nm_gw, &
628 eigenval_kp, num_integ_points, num_integ_group, color_rpa_group, &
629 fm_matrix_pq, fm_mat_s, fm_mat_q_gemm, fm_mat_q, fm_mat_s_gw, fm_mat_r_gw(1), &
630 fm_mat_s_ij_bse(1), fm_mat_s_ab_bse(1), &
631 my_do_gw, do_bse, gw_corr_lev_occ, gw_corr_lev_virt, &
634 do_im_time, mo_coeff, &
635 fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, &
636 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, mat_munu, mat_p_global, &
637 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
638 starts_array_mc, ends_array_mc, &
639 starts_array_mc_block, ends_array_mc_block, &
640 matrix_s, do_kpoints_from_gamma, kpoints, gd_array, color_sub, &
641 do_ri_sos_laplace_mp2=do_ri_sos_laplace_mp2, calc_forces=calc_forces)
647 IF (.NOT. do_im_time)
THEN
653 IF (my_do_gw .AND. .NOT. do_im_time)
THEN
663 CALL timestop(handle)
684 SUBROUTINE calculate_bib_c_2d(BIb_C_2D, BIb_C, para_env_sub, dimen_ia, homo, virtual, &
686 my_ia_size, my_ia_start, my_ia_end, my_group_L_size)
688 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :), &
689 INTENT(OUT) :: bib_c_2d
690 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
693 INTEGER,
INTENT(IN) :: dimen_ia, homo, virtual
695 INTEGER :: my_ia_size, my_ia_start, my_ia_end, &
698 INTEGER,
PARAMETER :: occ_chunk = 128
700 INTEGER :: ia_global, iib, itmp(2), jjb, my_b_size, my_b_virtual_start, occ_high, occ_low, &
701 proc_receive, proc_send, proc_shift, rec_b_size, rec_b_virtual_end, rec_b_virtual_start
702 REAL(kind=
dp),
ALLOCATABLE,
DIMENSION(:),
TARGET :: bib_c_rec_1d
703 REAL(kind=
dp),
DIMENSION(:, :, :),
POINTER :: bib_c_rec
705 itmp =
get_limit(dimen_ia, para_env_sub%num_pe, para_env_sub%mepos)
706 my_ia_start = itmp(1)
708 my_ia_size = my_ia_end - my_ia_start + 1
710 CALL get_group_dist(gd_b_virtual, para_env_sub%mepos, sizes=my_b_size, starts=my_b_virtual_start)
713 ALLOCATE (bib_c_2d(my_group_l_size, my_ia_size))
719 DO jjb = 1, my_b_size
720 ia_global = (iib - 1)*virtual + my_b_virtual_start + jjb - 1
721 IF (ia_global >= my_ia_start .AND. ia_global <= my_ia_end)
THEN
722 bib_c_2d(1:my_group_l_size, ia_global - my_ia_start + 1) = bib_c(1:my_group_l_size, jjb, iib)
727 IF (para_env_sub%num_pe > 1)
THEN
728 ALLOCATE (bib_c_rec_1d(int(my_group_l_size,
int_8)*
maxsize(gd_b_virtual)*min(homo, occ_chunk)))
729 DO proc_shift = 1, para_env_sub%num_pe - 1
730 proc_send =
modulo(para_env_sub%mepos + proc_shift, para_env_sub%num_pe)
731 proc_receive =
modulo(para_env_sub%mepos - proc_shift, para_env_sub%num_pe)
733 CALL get_group_dist(gd_b_virtual, proc_receive, rec_b_virtual_start, rec_b_virtual_end, rec_b_size)
736 DO occ_low = 1, homo, occ_chunk
737 occ_high = min(homo, occ_low + occ_chunk - 1)
738 bib_c_rec(1:my_group_l_size, 1:rec_b_size, 1:occ_high - occ_low + 1) => &
739 bib_c_rec_1d(1:int(my_group_l_size,
int_8)*rec_b_size*(occ_high - occ_low + 1))
740 CALL para_env_sub%sendrecv(bib_c(:, :, occ_low:occ_high), proc_send, &
741 bib_c_rec(:, :, 1:occ_high - occ_low + 1), proc_receive)
745 DO iib = occ_low, occ_high
746 DO jjb = 1, rec_b_size
747 ia_global = (iib - 1)*virtual + rec_b_virtual_start + jjb - 1
748 IF (ia_global >= my_ia_start .AND. ia_global <= my_ia_end)
THEN
749 bib_c_2d(1:my_group_l_size, ia_global - my_ia_start + 1) = bib_c_rec(1:my_group_l_size, jjb, iib - occ_low + 1)
756 DEALLOCATE (bib_c_rec_1d)
759 END SUBROUTINE calculate_bib_c_2d
793 SUBROUTINE create_integ_mat(BIb_C_2D, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
794 dimen_RI, dimen_ia, color_rpa_group, &
795 ext_row_block_size, ext_col_block_size, unit_nr, &
796 my_ia_size, my_ia_start, my_ia_end, &
797 my_group_L_size, my_group_L_start, my_group_L_end, &
798 para_env_RPA, fm_mat_S, nrow_block_mat, ncol_block_mat, &
799 blacs_env_ext, blacs_env_ext_S, dimen_ia_for_block_size, &
800 do_im_time, fm_mat_Q_gemm, fm_mat_Q, qs_env)
803 INTENT(INOUT) :: bib_c_2d
805 INTEGER,
INTENT(IN) :: color_sub, ngroup, integ_group_size, &
807 INTEGER,
DIMENSION(:),
INTENT(IN) :: dimen_ia
808 INTEGER,
INTENT(IN) :: color_rpa_group, ext_row_block_size, &
809 ext_col_block_size, unit_nr
810 INTEGER,
DIMENSION(:),
INTENT(IN) :: my_ia_size, my_ia_start, my_ia_end
811 INTEGER,
INTENT(IN) :: my_group_l_size, my_group_l_start, &
814 TYPE(
cp_fm_type),
DIMENSION(:),
INTENT(INOUT) :: fm_mat_s
815 INTEGER,
INTENT(INOUT) :: nrow_block_mat, ncol_block_mat
817 INTEGER,
INTENT(IN),
OPTIONAL :: dimen_ia_for_block_size
818 LOGICAL,
INTENT(IN),
OPTIONAL :: do_im_time
819 TYPE(
cp_fm_type),
DIMENSION(:),
OPTIONAL :: fm_mat_q_gemm, fm_mat_q
823 CHARACTER(LEN=*),
PARAMETER :: routinen =
'create_integ_mat'
825 INTEGER :: col_row_proc_ratio, grid_2d(2), handle, &
826 iproc, iproc_col, iproc_row, ispin, &
828 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: group_grid_2_mepos
829 LOGICAL :: my_blacs_ext, my_blacs_s_ext, &
835 CALL timeset(routinen, handle)
837 cpassert(
PRESENT(blacs_env_ext) .OR.
PRESENT(dimen_ia_for_block_size))
839 my_blacs_ext = .false.
840 IF (
PRESENT(blacs_env_ext)) my_blacs_ext = .true.
842 my_blacs_s_ext = .false.
843 IF (
PRESENT(blacs_env_ext_s)) my_blacs_s_ext = .true.
845 my_do_im_time = .false.
846 IF (
PRESENT(do_im_time)) my_do_im_time = do_im_time
850 IF (my_blacs_s_ext)
THEN
851 blacs_env => blacs_env_ext_s
853 IF (para_env_rpa%num_pe > 1)
THEN
854 col_row_proc_ratio = max(1, dimen_ia_for_block_size/dimen_ri)
856 iproc_col = min(max(int(sqrt(real(para_env_rpa%num_pe*col_row_proc_ratio, kind=
dp))), 1), para_env_rpa%num_pe) + 1
857 DO iproc = 1, para_env_rpa%num_pe
858 iproc_col = iproc_col - 1
859 IF (mod(para_env_rpa%num_pe, iproc_col) == 0)
EXIT
862 iproc_row = para_env_rpa%num_pe/iproc_col
863 grid_2d(1) = iproc_row
864 grid_2d(2) = iproc_col
870 IF (unit_nr > 0 .AND. .NOT. my_do_im_time)
THEN
871 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
872 "MATRIX_INFO| Number row processes:", grid_2d(1)
873 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
874 "MATRIX_INFO| Number column processes:", grid_2d(2)
878 IF (ext_row_block_size > 0)
THEN
879 nrow_block_mat = ext_row_block_size
881 nrow_block_mat = max(1, dimen_ri/grid_2d(1)/2)
885 IF (ext_col_block_size > 0)
THEN
886 ncol_block_mat = ext_col_block_size
888 ncol_block_mat = max(1, dimen_ia_for_block_size/grid_2d(2)/2)
891 IF (unit_nr > 0 .AND. .NOT. my_do_im_time)
THEN
892 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
893 "MATRIX_INFO| Row block size:", nrow_block_mat
894 WRITE (unit=unit_nr, fmt=
"(T3,A,T75,i6)") &
895 "MATRIX_INFO| Column block size:", ncol_block_mat
899 IF (.NOT. my_do_im_time)
THEN
900 DO ispin = 1,
SIZE(bib_c_2d)
902 IF (my_blacs_ext)
THEN
904 ncol_global=dimen_ia(ispin), para_env=para_env_rpa)
907 ncol_global=dimen_ia(ispin), para_env=para_env_rpa, &
908 nrow_block=nrow_block_mat, ncol_block=ncol_block_mat, force_block=.true.)
912 CALL create_group_dist(gd_ia, my_ia_start(ispin), my_ia_end(ispin), my_ia_size(ispin), para_env_rpa)
913 CALL create_group_dist(gd_l, my_group_l_start, my_group_l_end, my_group_l_size, para_env_rpa)
917 mepos_in_rpa_group = mod(color_sub, integ_group_size)
918 ALLOCATE (group_grid_2_mepos(0:integ_group_size - 1, 0:para_env_sub%num_pe - 1))
919 group_grid_2_mepos = 0
920 group_grid_2_mepos(mepos_in_rpa_group, para_env_sub%mepos) = para_env_rpa%mepos
921 CALL para_env_rpa%sum(group_grid_2_mepos)
923 CALL array2fm(bib_c_2d(ispin)%array, fm_struct, my_group_l_start, my_group_l_end, &
924 my_ia_start(ispin), my_ia_end(ispin), gd_l, gd_ia, &
925 group_grid_2_mepos, ngroup, para_env_sub%num_pe, fm_mat_s(ispin), &
926 integ_group_size, color_rpa_group)
928 DEALLOCATE (group_grid_2_mepos)
936 IF (para_env_rpa%num_pe /= para_env%num_pe)
THEN
939 comm_exchange = fm_mat_s(ispin)%matrix_struct%context%interconnect(para_env)
940 CALL comm_exchange%sum(fm_mat_s(ispin)%local_data)
941 CALL comm_exchange%free()
947 IF (
PRESENT(fm_mat_q_gemm) .AND. .NOT. my_do_im_time)
THEN
951 ncol_global=dimen_ri, para_env=para_env_rpa, &
952 nrow_block=nrow_block_mat, ncol_block=ncol_block_mat, force_block=.true.)
953 DO ispin = 1,
SIZE(fm_mat_q_gemm)
954 CALL cp_fm_create(fm_mat_q_gemm(ispin), fm_struct, name=
"fm_mat_Q_gemm")
959 IF (
PRESENT(fm_mat_q))
THEN
960 NULLIFY (blacs_env_q)
961 IF (my_blacs_ext)
THEN
962 blacs_env_q => blacs_env_ext
963 ELSE IF (para_env_rpa%num_pe == para_env%num_pe .AND.
PRESENT(qs_env))
THEN
964 CALL get_qs_env(qs_env, blacs_env=blacs_env_q)
970 ncol_global=dimen_ri, para_env=para_env_rpa)
971 DO ispin = 1,
SIZE(fm_mat_q)
972 CALL cp_fm_create(fm_mat_q(ispin), fm_struct, name=
"fm_mat_Q")
977 IF (.NOT. (my_blacs_ext .OR. (para_env_rpa%num_pe == para_env%num_pe .AND.
PRESENT(qs_env)))) &
982 IF (.NOT. my_blacs_s_ext)
THEN
988 CALL timestop(handle)
990 END SUBROUTINE create_integ_mat
1049 SUBROUTINE rpa_num_int(qs_env, Erpa, mp2_env, para_env, para_env_RPA, para_env_sub, unit_nr, &
1050 homo, virtual, dimen_RI, dimen_RI_red, dimen_ia, dimen_nm_gw, &
1051 Eigenval, num_integ_points, num_integ_group, color_rpa_group, &
1052 fm_matrix_PQ, fm_mat_S, fm_mat_Q_gemm, fm_mat_Q, fm_mat_S_gw, fm_mat_R_gw, &
1053 fm_mat_S_ij_bse, fm_mat_S_ab_bse, &
1054 my_do_gw, do_bse, gw_corr_lev_occ, gw_corr_lev_virt, &
1056 do_minimax_quad, do_im_time, mo_coeff, &
1057 fm_matrix_L_kpoints, fm_matrix_Minv_L_kpoints, &
1058 fm_matrix_Minv, fm_matrix_Minv_Vtrunc_Minv, mat_munu, mat_P_global, &
1059 t_3c_M, t_3c_O, t_3c_O_compressed, t_3c_O_ind, &
1060 starts_array_mc, ends_array_mc, &
1061 starts_array_mc_block, ends_array_mc_block, &
1062 matrix_s, do_kpoints_from_Gamma, kpoints, gd_array, color_sub, &
1063 do_ri_sos_laplace_mp2, calc_forces)
1065 TYPE(qs_environment_type),
POINTER :: qs_env
1066 REAL(kind=dp),
INTENT(OUT) :: erpa
1067 TYPE(mp2_type) :: mp2_env
1068 TYPE(mp_para_env_type),
POINTER :: para_env, para_env_rpa, para_env_sub
1069 INTEGER,
INTENT(IN) :: unit_nr
1070 INTEGER,
DIMENSION(:),
INTENT(IN) :: homo, virtual
1071 INTEGER,
INTENT(IN) :: dimen_ri, dimen_ri_red
1072 INTEGER,
DIMENSION(:),
INTENT(IN) :: dimen_ia
1073 INTEGER,
INTENT(IN) :: dimen_nm_gw
1074 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
1075 INTENT(INOUT) :: eigenval
1076 INTEGER,
INTENT(IN) :: num_integ_points, num_integ_group, &
1078 TYPE(cp_fm_type),
INTENT(IN) :: fm_matrix_pq
1079 TYPE(cp_fm_type),
DIMENSION(:),
INTENT(INOUT) :: fm_mat_s
1080 TYPE(cp_fm_type),
DIMENSION(:),
INTENT(IN) :: fm_mat_q_gemm, fm_mat_q, fm_mat_s_gw
1081 TYPE(cp_fm_type),
INTENT(IN) :: fm_mat_r_gw, fm_mat_s_ij_bse, &
1083 LOGICAL,
INTENT(IN) :: my_do_gw, do_bse
1084 INTEGER,
DIMENSION(:),
INTENT(IN) :: gw_corr_lev_occ, gw_corr_lev_virt
1085 INTEGER,
INTENT(IN) :: bse_lev_virt
1086 LOGICAL,
INTENT(IN) :: do_minimax_quad, do_im_time
1087 TYPE(cp_fm_type),
DIMENSION(:),
INTENT(IN) :: mo_coeff
1088 TYPE(cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_matrix_l_kpoints, &
1089 fm_matrix_minv_l_kpoints, &
1091 fm_matrix_minv_vtrunc_minv
1092 TYPE(dbcsr_p_type),
INTENT(IN) :: mat_munu
1093 TYPE(dbcsr_p_type),
INTENT(INOUT) :: mat_p_global
1094 TYPE(dbt_type),
INTENT(INOUT) :: t_3c_m
1095 TYPE(dbt_type),
ALLOCATABLE,
DIMENSION(:, :), &
1096 INTENT(INOUT) :: t_3c_o
1097 TYPE(hfx_compression_type),
ALLOCATABLE, &
1098 DIMENSION(:, :, :),
INTENT(INOUT) :: t_3c_o_compressed
1099 TYPE(block_ind_type),
ALLOCATABLE, &
1100 DIMENSION(:, :, :),
INTENT(INOUT) :: t_3c_o_ind
1101 INTEGER,
ALLOCATABLE,
DIMENSION(:),
INTENT(IN) :: starts_array_mc, ends_array_mc, &
1102 starts_array_mc_block, &
1104 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_s
1105 LOGICAL :: do_kpoints_from_gamma
1106 TYPE(kpoint_type),
POINTER :: kpoints
1107 TYPE(group_dist_d1_type),
INTENT(IN) :: gd_array
1108 INTEGER,
INTENT(IN) :: color_sub
1109 LOGICAL,
INTENT(IN) :: do_ri_sos_laplace_mp2, calc_forces
1111 CHARACTER(LEN=*),
PARAMETER :: routinen =
'rpa_num_int'
1113 COMPLEX(KIND=dp),
ALLOCATABLE, &
1114 DIMENSION(:, :, :, :) :: vec_sigma_c_gw
1115 INTEGER :: count_ev_sc_gw, cut_memory, group_size_p, gw_corr_lev_tot, handle, handle3, i, &
1116 ikp_local, ispin, iter_evgw, iter_sc_gw0, j, jquad, min_bsize, mm_style, nkp, &
1117 nkp_self_energy, nmo, nspins, num_3c_repl, num_cells_dm, num_fit_points, pspin, qspin, &
1119 INTEGER(int_8) :: dbcsr_nflop
1120 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: index_to_cell_3c
1121 INTEGER,
ALLOCATABLE,
DIMENSION(:, :, :) :: cell_to_index_3c
1122 INTEGER,
DIMENSION(:),
POINTER :: col_blk_size, prim_blk_sizes, &
1124 LOGICAL :: do_apply_ic_corr_to_gw, do_gw_im_time, do_ic_model, do_kpoints_cubic_rpa, &
1125 do_periodic, do_print, do_ri_sigma_x, exit_ev_gw, first_cycle, &
1126 first_cycle_periodic_correction, my_open_shell, print_ic_values
1127 LOGICAL,
ALLOCATABLE,
DIMENSION(:, :, :, :, :) :: has_mat_p_blocks
1128 REAL(kind=dp) :: a_scaling, alpha, dbcsr_time, e_exchange, e_exchange_corr, eps_filter, &
1129 eps_filter_im_time, ext_scaling, fermi_level_offset, fermi_level_offset_input, &
1130 my_flop_rate, omega, omega_max_fit, omega_old, tau, tau_old
1131 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:) :: delta_corr, e_fermi, tau_tj, tau_wj, tj, &
1132 trace_qomega, vec_omega_fit_gw, wj, &
1134 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :) :: vec_w_gw, weights_cos_tf_t_to_w, &
1135 weights_cos_tf_w_to_t, &
1136 weights_sin_tf_t_to_w
1137 REAL(kind=dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: eigenval_last, eigenval_scf, &
1139 TYPE(cp_cfm_type) :: cfm_mat_q
1140 TYPE(cp_fm_type) :: fm_mat_q_static_bse_gemm, fm_mat_ri_global_work, fm_mat_s_ia_bse, &
1141 fm_mat_work, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, fm_scaled_dm_occ_tau, &
1142 fm_scaled_dm_virt_tau
1143 TYPE(cp_fm_type),
ALLOCATABLE,
DIMENSION(:) :: fm_mat_s_gw_work, fm_mat_w, &
1144 fm_mo_coeff_occ, fm_mo_coeff_virt
1145 TYPE(cp_fm_type),
ALLOCATABLE,
DIMENSION(:, :) :: fm_mat_l_kpoints, fm_mat_minv_l_kpoints
1146 TYPE(dbcsr_p_type) :: mat_dm, mat_l, mat_m_p_munu_occ, &
1147 mat_m_p_munu_virt, mat_minvvminv
1148 TYPE(dbcsr_p_type),
ALLOCATABLE, &
1149 DIMENSION(:, :, :) :: mat_p_omega
1150 TYPE(dbcsr_p_type),
DIMENSION(:),
POINTER :: matrix_berry_im_mo_mo, &
1151 matrix_berry_re_mo_mo
1152 TYPE(dbcsr_p_type),
DIMENSION(:, :),
POINTER :: mat_p_omega_kp
1153 TYPE(dbcsr_type),
POINTER :: mat_w, mat_work
1154 TYPE(dbt_type) :: t_3c_overl_int_ao_mo
1155 TYPE(dbt_type),
ALLOCATABLE,
DIMENSION(:) :: t_3c_overl_int_gw_ao, &
1156 t_3c_overl_int_gw_ri, &
1157 t_3c_overl_nnp_ic, &
1158 t_3c_overl_nnp_ic_reflected
1159 TYPE(dgemm_counter_type) :: dgemm_counter
1160 TYPE(hfx_compression_type),
ALLOCATABLE, &
1161 DIMENSION(:) :: t_3c_o_mo_compressed
1162 TYPE(im_time_force_type) :: force_data
1163 TYPE(rpa_exchange_work_type) :: exchange_work
1165 TYPE(rpa_sigma_type) :: rpa_sigma
1166 TYPE(two_dim_int_array),
ALLOCATABLE,
DIMENSION(:) :: t_3c_o_mo_ind
1168 CALL timeset(routinen, handle)
1171 nmo = homo(1) + virtual(1)
1173 my_open_shell = (nspins == 2)
1175 do_gw_im_time = my_do_gw .AND. do_im_time
1176 do_ri_sigma_x = mp2_env%ri_g0w0%do_ri_Sigma_x
1177 do_ic_model = mp2_env%ri_g0w0%do_ic_model
1178 print_ic_values = mp2_env%ri_g0w0%print_ic_values
1179 do_periodic = mp2_env%ri_g0w0%do_periodic
1180 do_kpoints_cubic_rpa = mp2_env%ri_rpa_im_time%do_im_time_kpoints
1183 mm_style = wfc_mm_style_gemm
1184 IF (.NOT. do_ri_sos_laplace_mp2) mm_style = mp2_env%ri_rpa%mm_style
1187 ext_scaling = 0.2_dp
1188 omega_max_fit = mp2_env%ri_g0w0%omega_max_fit
1189 fermi_level_offset_input = mp2_env%ri_g0w0%fermi_level_offset
1190 iter_evgw = mp2_env%ri_g0w0%iter_evGW
1191 iter_sc_gw0 = mp2_env%ri_g0w0%iter_sc_GW0
1192 IF ((.NOT. do_im_time))
THEN
1193 IF (iter_sc_gw0 .NE. 1 .AND. iter_evgw .NE. 1) cpabort(
"Mixed scGW0/evGW not implemented.")
1195 IF (iter_sc_gw0 .NE. 1) iter_evgw = iter_sc_gw0
1198 ext_scaling = 0.0_dp
1203 IF (do_kpoints_cubic_rpa .AND. do_ri_sos_laplace_mp2)
THEN
1204 cpabort(
"RI-SOS-Laplace-MP2 with k-point-sampling is not implemented.")
1207 do_apply_ic_corr_to_gw = .false.
1208 IF (mp2_env%ri_g0w0%ic_corr_list(1)%array(1) > 0.0_dp) do_apply_ic_corr_to_gw = .true.
1210 IF (do_im_time)
THEN
1211 cpassert(do_minimax_quad .OR. do_ri_sos_laplace_mp2)
1213 group_size_p = mp2_env%ri_rpa_im_time%group_size_P
1214 cut_memory = mp2_env%ri_rpa_im_time%cut_memory
1215 eps_filter = mp2_env%ri_rpa_im_time%eps_filter
1216 eps_filter_im_time = mp2_env%ri_rpa_im_time%eps_filter* &
1217 mp2_env%ri_rpa_im_time%eps_filter_factor
1219 min_bsize = mp2_env%ri_rpa_im_time%min_bsize
1221 CALL alloc_im_time(qs_env, para_env, dimen_ri, dimen_ri_red, &
1222 num_integ_points, nspins, fm_mat_q(1), fm_mo_coeff_occ, fm_mo_coeff_virt, &
1223 fm_matrix_minv_l_kpoints, fm_matrix_l_kpoints, mat_p_global, &
1224 t_3c_o, matrix_s, kpoints, eps_filter_im_time, &
1225 cut_memory, nkp, num_cells_dm, num_3c_repl, &
1226 size_p, ikp_local, &
1230 do_ic_model, do_kpoints_cubic_rpa, &
1231 do_kpoints_from_gamma, do_ri_sigma_x, my_open_shell, &
1232 has_mat_p_blocks, wkp_w, &
1233 cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1234 fm_mat_ri_global_work, fm_mat_work, fm_mo_coeff_occ_scaled, &
1235 fm_mo_coeff_virt_scaled, mat_dm, mat_l, mat_m_p_munu_occ, mat_m_p_munu_virt, &
1236 mat_minvvminv, mat_p_omega, mat_p_omega_kp, mat_work, mo_coeff, &
1237 fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, homo, nmo)
1239 IF (calc_forces)
CALL init_im_time_forces(force_data, fm_matrix_pq, t_3c_m, unit_nr, mp2_env, qs_env)
1243 CALL dbcsr_get_info(mat_p_global%matrix, &
1244 row_blk_size=ri_blk_sizes)
1246 CALL dbcsr_get_info(matrix_s(1)%matrix, &
1247 row_blk_size=prim_blk_sizes)
1249 gw_corr_lev_tot = gw_corr_lev_occ(1) + gw_corr_lev_virt(1)
1251 IF (.NOT. do_kpoints_cubic_rpa)
THEN
1252 CALL allocate_matrices_gw_im_time(gw_corr_lev_occ, gw_corr_lev_virt, homo, nmo, &
1253 num_integ_points, unit_nr, &
1254 ri_blk_sizes, do_ic_model, &
1255 para_env, fm_mat_w, fm_mat_q(1), &
1257 t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, &
1258 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1259 starts_array_mc, ends_array_mc, &
1260 t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, &
1261 matrix_s, mat_w, t_3c_o, &
1262 t_3c_o_compressed, t_3c_o_ind, &
1269 IF (do_ic_model)
THEN
1271 cpassert(do_gw_im_time)
1272 cpassert(.NOT. do_periodic)
1273 IF (cut_memory .NE. 1) cpabort(
"For IC, use MEMORY_CUT 1 in the LOW_SCALING section.")
1276 ALLOCATE (e_fermi(nspins))
1277 IF (do_minimax_quad .OR. do_ri_sos_laplace_mp2)
THEN
1278 do_print = .NOT. do_ic_model
1279 CALL get_minimax_grid(para_env, unit_nr, homo, eigenval, num_integ_points, do_im_time, &
1280 do_ri_sos_laplace_mp2, do_print, &
1281 tau_tj, tau_wj, qs_env, do_gw_im_time, &
1282 do_kpoints_cubic_rpa, e_fermi(1), tj, wj, &
1283 weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, &
1284 qs_env%mp2_env%ri_g0w0%regularization_minimax)
1287 IF (qs_env%mp2_env%ri_rpa_im_time%keep_quad)
THEN
1288 CALL keep_initial_quad(tj, wj, tau_tj, tau_wj, weights_cos_tf_t_to_w, &
1289 weights_cos_tf_w_to_t, do_ri_sos_laplace_mp2, do_im_time, &
1290 num_integ_points, unit_nr, qs_env)
1293 IF (calc_forces) cpabort(
"Forces with Clenshaw-Curtis grid not implemented.")
1294 CALL get_clenshaw_grid(para_env, para_env_rpa, unit_nr, homo, virtual, eigenval, num_integ_points, &
1295 num_integ_group, color_rpa_group, fm_mat_s, my_do_gw, &
1296 ext_scaling, a_scaling, tj, wj)
1300 IF (.NOT. do_ri_sos_laplace_mp2)
THEN
1301 ALLOCATE (trace_qomega(dimen_ri_red))
1304 IF (do_ri_sos_laplace_mp2 .AND. .NOT. do_im_time)
THEN
1306 ELSE IF (my_open_shell .OR. do_ri_sos_laplace_mp2)
THEN
1312 CALL allocate_matrices_gw(vec_sigma_c_gw, color_rpa_group, dimen_nm_gw, &
1313 gw_corr_lev_occ, gw_corr_lev_virt, homo, &
1314 nmo, num_integ_group, num_integ_points, unit_nr, &
1315 gw_corr_lev_tot, num_fit_points, omega_max_fit, &
1316 do_minimax_quad, do_periodic, do_ri_sigma_x,.NOT. do_im_time, &
1317 first_cycle_periodic_correction, &
1318 a_scaling, eigenval, tj, vec_omega_fit_gw, vec_sigma_x_gw, &
1319 delta_corr, eigenval_last, eigenval_scf, vec_w_gw, &
1320 fm_mat_s_gw, fm_mat_s_gw_work, &
1321 para_env, mp2_env, kpoints, nkp, nkp_self_energy, &
1322 do_kpoints_cubic_rpa, do_kpoints_from_gamma)
1326 CALL cp_fm_create(fm_mat_q_static_bse_gemm, fm_mat_q_gemm(1)%matrix_struct)
1327 CALL cp_fm_to_fm(fm_mat_q_gemm(1), fm_mat_q_static_bse_gemm)
1328 CALL cp_fm_set_all(fm_mat_q_static_bse_gemm, 0.0_dp)
1334 IF (calc_forces .AND. .NOT. do_im_time)
CALL rpa_grad_create(
rpa_grad, fm_mat_q(1), &
1335 fm_mat_s, homo, virtual, mp2_env, eigenval(:, 1, :), &
1336 unit_nr, do_ri_sos_laplace_mp2)
1337 IF (.NOT. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2) &
1338 CALL exchange_work%create(qs_env, para_env_sub, mat_munu, dimen_ri_red, &
1339 fm_mat_s, fm_mat_q(1), fm_mat_q_gemm(1), homo, virtual)
1341 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) e_exchange = 0.0_dp
1342 first_cycle = .true.
1344 CALL dgemm_counter_init(dgemm_counter, unit_nr, mp2_env%ri_rpa%print_dgemm_info)
1346 DO count_ev_sc_gw = 1, iter_evgw
1350 IF (do_ic_model) cycle
1355 vec_sigma_c_gw = z_zero
1356 first_cycle = .true.
1360 IF (do_im_time)
THEN
1362 IF (.NOT. do_kpoints_cubic_rpa)
THEN
1363 DO ispin = 1, nspins
1364 e_fermi(ispin) = (eigenval(homo(ispin), 1, ispin) + eigenval(homo(ispin) + 1, 1, ispin))*0.5_dp
1371 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(/T3,A,T66,i15)") &
1372 "MEMORY_INFO| Memory cut:", cut_memory
1373 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3,A,T66,ES15.2)") &
1374 "SPARSITY_INFO| Eps filter for M virt/occ tensors:", eps_filter
1375 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3,A,T66,ES15.2)") &
1376 "SPARSITY_INFO| Eps filter for P matrix:", eps_filter_im_time
1377 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3,A,T66,i15)") &
1378 "SPARSITY_INFO| Minimum tensor block size:", min_bsize
1381 CALL zero_mat_p_omega(mat_p_omega(:, :, 1))
1384 CALL compute_mat_p_omega(mat_p_omega(:, :, 1), fm_scaled_dm_occ_tau, &
1385 fm_scaled_dm_virt_tau, fm_mo_coeff_occ(1), fm_mo_coeff_virt(1), &
1386 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1387 mat_p_global, matrix_s, 1, &
1388 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
1389 starts_array_mc, ends_array_mc, &
1390 starts_array_mc_block, ends_array_mc_block, &
1391 weights_cos_tf_t_to_w, tj, tau_tj, e_fermi(1), eps_filter, alpha, &
1392 eps_filter_im_time, eigenval(:, 1, 1), nmo, &
1393 num_integ_points, cut_memory, &
1394 unit_nr, mp2_env, para_env, &
1395 qs_env, do_kpoints_from_gamma, &
1396 index_to_cell_3c, cell_to_index_3c, &
1397 has_mat_p_blocks, do_ri_sos_laplace_mp2, &
1398 dbcsr_time, dbcsr_nflop)
1401 IF (my_open_shell)
THEN
1402 CALL zero_mat_p_omega(mat_p_omega(:, :, 2))
1403 CALL compute_mat_p_omega(mat_p_omega(:, :, 2), fm_scaled_dm_occ_tau, &
1404 fm_scaled_dm_virt_tau, fm_mo_coeff_occ(2), &
1405 fm_mo_coeff_virt(2), &
1406 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1407 mat_p_global, matrix_s, 2, &
1408 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
1409 starts_array_mc, ends_array_mc, &
1410 starts_array_mc_block, ends_array_mc_block, &
1411 weights_cos_tf_t_to_w, tj, tau_tj, e_fermi(2), eps_filter, alpha, &
1412 eps_filter_im_time, eigenval(:, 1, 2), nmo, &
1413 num_integ_points, cut_memory, &
1414 unit_nr, mp2_env, para_env, &
1415 qs_env, do_kpoints_from_gamma, &
1416 index_to_cell_3c, cell_to_index_3c, &
1417 has_mat_p_blocks, do_ri_sos_laplace_mp2, &
1418 dbcsr_time, dbcsr_nflop)
1421 IF (.NOT. do_ri_sos_laplace_mp2)
THEN
1422 DO j = 1,
SIZE(mat_p_omega, 2)
1423 DO i = 1,
SIZE(mat_p_omega, 1)
1424 CALL dbcsr_add(mat_p_omega(i, j, 1)%matrix, mat_p_omega(i, j, 2)%matrix, 1.0_dp, 1.0_dp)
1425 IF (.NOT. calc_forces)
CALL dbcsr_clear(mat_p_omega(i, j, 2)%matrix)
1433 IF (mp2_env%ri_rpa%sigma_param /= sigma_none)
THEN
1434 CALL rpa_sigma_create(rpa_sigma, mp2_env%ri_rpa%sigma_param, fm_mat_q(1), unit_nr, para_env)
1437 DO jquad = 1, num_integ_points
1438 IF (
modulo(jquad, num_integ_group) /= color_rpa_group) cycle
1440 CALL timeset(routinen//
"_RPA_matrix_operations", handle3)
1442 IF (do_ri_sos_laplace_mp2)
THEN
1443 omega = tau_tj(jquad)
1445 IF (do_minimax_quad)
THEN
1448 omega = a_scaling/tan(tj(jquad))
1452 IF (do_im_time)
THEN
1455 IF (.NOT. (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma))
THEN
1457 DO ispin = 1,
SIZE(mat_p_omega, 3)
1458 CALL contract_p_omega_with_mat_l(mat_p_omega(jquad, 1, ispin)%matrix, mat_l%matrix, mat_work, &
1459 eps_filter_im_time, fm_mat_work, dimen_ri, dimen_ri_red, &
1460 fm_mat_minv_l_kpoints(1, 1), fm_mat_q(ispin))
1465 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3, A, 1X, I3, 1X, A, 1X, I3)") &
1466 "INTEG_INFO| Started with Integration point", jquad,
"of", num_integ_points
1468 IF (first_cycle .AND. count_ev_sc_gw > 1)
THEN
1469 IF (iter_sc_gw0 == 1)
THEN
1470 DO ispin = 1, nspins
1471 CALL remove_scaling_factor_rpa(fm_mat_s(ispin), virtual(ispin), &
1472 eigenval_last(:, 1, ispin), homo(ispin), omega_old)
1475 DO ispin = 1, nspins
1476 CALL remove_scaling_factor_rpa(fm_mat_s(ispin), virtual(ispin), &
1477 eigenval_scf(:, 1, ispin), homo(ispin), omega_old)
1482 IF (iter_sc_gw0 > 1)
THEN
1483 DO ispin = 1, nspins
1484 CALL calc_mat_q(fm_mat_s(ispin), do_ri_sos_laplace_mp2, first_cycle, virtual(ispin), &
1485 eigenval_scf(:, 1, ispin), homo(ispin), omega, omega_old, jquad, mm_style, &
1486 dimen_ri_red, dimen_ia(ispin), alpha, fm_mat_q(ispin), &
1487 fm_mat_q_gemm(ispin), do_bse, fm_mat_q_static_bse_gemm, dgemm_counter, &
1488 num_integ_points, count_ev_sc_gw)
1492 IF (.NOT. do_ri_sos_laplace_mp2)
THEN
1493 DO ispin = 2, nspins
1494 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_q(1), beta=1.0_dp, matrix_b=fm_mat_q(ispin))
1498 DO ispin = 1, nspins
1499 CALL calc_mat_q(fm_mat_s(ispin), do_ri_sos_laplace_mp2, first_cycle, virtual(ispin), &
1500 eigenval(:, 1, ispin), homo(ispin), omega, omega_old, jquad, mm_style, &
1501 dimen_ri_red, dimen_ia(ispin), alpha, fm_mat_q(ispin), &
1502 fm_mat_q_gemm(ispin), do_bse, fm_mat_q_static_bse_gemm, dgemm_counter, &
1503 num_integ_points, count_ev_sc_gw)
1507 IF (.NOT. do_ri_sos_laplace_mp2)
THEN
1508 DO ispin = 2, nspins
1509 CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=fm_mat_q(1), beta=1.0_dp, matrix_b=fm_mat_q(ispin))
1518 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none)
THEN
1519 e_exchange_corr = 0.0_dp
1520 CALL exchange_work%compute(fm_mat_q(1), eigenval(:, 1, :), fm_mat_s, omega, e_exchange_corr, mp2_env)
1523 e_exchange = e_exchange + e_exchange_corr*wj(jquad)
1527 IF (mp2_env%ri_rpa%sigma_param /= sigma_none)
THEN
1528 CALL rpa_sigma_matrix_spectral(rpa_sigma, fm_mat_q(1), wj(jquad), para_env_rpa)
1531 IF (do_ri_sos_laplace_mp2)
THEN
1533 CALL sos_mp2_postprocessing(fm_mat_q, erpa, tau_wj(jquad))
1535 IF (calc_forces .AND. .NOT. do_im_time)
CALL rpa_grad_matrix_operations(mp2_env,
rpa_grad, do_ri_sos_laplace_mp2, &
1536 fm_mat_q, fm_mat_q_gemm, dgemm_counter, fm_mat_s, omega, homo, virtual, &
1537 eigenval(:, 1, :), tau_wj(jquad), unit_nr)
1539 IF (calc_forces .AND. .NOT. do_im_time)
CALL rpa_grad_copy_q(fm_mat_q(1),
rpa_grad)
1541 CALL q_trace_and_add_unit_matrix(dimen_ri_red, trace_qomega, fm_mat_q(1))
1543 IF (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma)
THEN
1544 CALL invert_eps_compute_w_and_erpa_kp(dimen_ri, num_integ_points, jquad, nkp, count_ev_sc_gw, para_env, &
1545 erpa, tau_tj, tj, wj, weights_cos_tf_w_to_t, &
1546 wkp_w, do_gw_im_time, do_ri_sigma_x, do_kpoints_from_gamma, &
1547 cfm_mat_q, ikp_local, &
1548 mat_p_omega(:, :, 1), mat_p_omega_kp, qs_env, eps_filter_im_time, unit_nr, &
1549 kpoints, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1550 fm_mat_w, fm_mat_ri_global_work, mat_minvvminv, &
1551 fm_matrix_minv, fm_matrix_minv_vtrunc_minv)
1553 CALL compute_erpa_by_freq_int(dimen_ri_red, trace_qomega, fm_mat_q(1), para_env_rpa, erpa, wj(jquad))
1556 IF (calc_forces .AND. .NOT. do_im_time)
CALL rpa_grad_matrix_operations(mp2_env,
rpa_grad, do_ri_sos_laplace_mp2, &
1557 fm_mat_q, fm_mat_q_gemm, dgemm_counter, fm_mat_s, omega, homo, virtual, &
1558 eigenval(:, 1, :), wj(jquad), unit_nr)
1562 first_cycle = .false.
1565 CALL timestop(handle3)
1569 CALL get_fermi_level_offset(fermi_level_offset, fermi_level_offset_input, eigenval(:, 1, :), homo)
1572 IF (do_im_time)
THEN
1574 IF (.NOT. (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma))
THEN
1575 CALL compute_w_cubic_gw(fm_mat_w, fm_mat_q(1), fm_mat_work, dimen_ri, fm_mat_minv_l_kpoints, num_integ_points, &
1576 tj, tau_tj, weights_cos_tf_w_to_t, jquad, omega)
1579 CALL compute_gw_self_energy(vec_sigma_c_gw, dimen_nm_gw, dimen_ri_red, gw_corr_lev_occ, &
1580 gw_corr_lev_virt, homo, jquad, nmo, num_fit_points, &
1581 do_im_time, do_periodic, first_cycle_periodic_correction, &
1582 fermi_level_offset, &
1583 omega, eigenval(:, 1, :), delta_corr, vec_omega_fit_gw, vec_w_gw, wj, &
1584 fm_mat_q(1), fm_mat_r_gw, fm_mat_s_gw, &
1585 fm_mat_s_gw_work, mo_coeff(1), para_env, &
1586 para_env_rpa, matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, &
1587 kpoints, qs_env, mp2_env)
1591 IF (unit_nr > 0)
CALL m_flush(unit_nr)
1592 CALL para_env%sync()
1596 IF (mp2_env%ri_rpa%sigma_param /= sigma_none)
THEN
1597 CALL finalize_rpa_sigma(rpa_sigma, unit_nr, mp2_env%ri_rpa%e_sigma_corr, para_env, do_minimax_quad)
1598 IF (do_minimax_quad) mp2_env%ri_rpa%e_sigma_corr = mp2_env%ri_rpa%e_sigma_corr/2.0_dp
1599 CALL para_env%sum(mp2_env%ri_rpa%e_sigma_corr)
1602 CALL para_env%sum(erpa)
1604 IF (.NOT. (do_ri_sos_laplace_mp2))
THEN
1605 erpa = erpa/(pi*2.0_dp)
1606 IF (do_minimax_quad) erpa = erpa/2.0_dp
1609 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none)
THEN
1610 CALL para_env%sum(e_exchange)
1611 e_exchange = e_exchange/(pi*2.0_dp)
1612 IF (do_minimax_quad) e_exchange = e_exchange/2.0_dp
1613 mp2_env%ri_rpa%ener_exchange = e_exchange
1616 IF (calc_forces .AND. do_ri_sos_laplace_mp2 .AND. do_im_time)
THEN
1617 IF (my_open_shell)
THEN
1620 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1621 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1622 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1623 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1624 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1625 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1626 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1627 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1630 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1631 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1632 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1633 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1634 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1635 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1636 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1637 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1642 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1643 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1644 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1645 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1646 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1647 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1648 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1649 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1651 CALL calc_post_loop_forces(force_data, unit_nr, qs_env)
1654 IF (calc_forces .AND. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2)
THEN
1655 DO ispin = 1, nspins
1656 CALL calc_rpa_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1657 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1658 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1659 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1660 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1661 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1662 e_fermi(ispin), weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, tj, &
1663 wj, tau_tj, cut_memory, ispin, my_open_shell, unit_nr, dbcsr_time, &
1664 dbcsr_nflop, mp2_env, qs_env)
1666 CALL calc_post_loop_forces(force_data, unit_nr, qs_env)
1669 IF (do_im_time)
THEN
1671 my_flop_rate = real(dbcsr_nflop, dp)/(1.0e09_dp*dbcsr_time)
1672 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(/T3,A,T73,ES8.2)") &
1673 "PERFORMANCE| DBCSR total number of flops:", real(dbcsr_nflop*para_env%num_pe, dp)
1674 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3,A,T66,F15.2)") &
1675 "PERFORMANCE| DBCSR total execution time:", dbcsr_time
1676 IF (unit_nr > 0)
WRITE (unit=unit_nr, fmt=
"(T3,A,T66,F15.2)") &
1677 "PERFORMANCE| DBCSR flop rate (Gflops / MPI rank):", my_flop_rate
1681 CALL dgemm_counter_write(dgemm_counter, para_env)
1690 CALL compute_qp_energies(vec_sigma_c_gw, count_ev_sc_gw, gw_corr_lev_occ, &
1691 gw_corr_lev_tot, gw_corr_lev_virt, homo, &
1692 nmo, num_fit_points, num_integ_points, &
1693 unit_nr, do_apply_ic_corr_to_gw, do_im_time, &
1694 do_periodic, do_ri_sigma_x, first_cycle_periodic_correction, &
1695 e_fermi, eps_filter, fermi_level_offset, &
1696 delta_corr, eigenval, &
1697 eigenval_last, eigenval_scf, iter_sc_gw0, exit_ev_gw, tau_tj, tj, &
1698 vec_omega_fit_gw, vec_sigma_x_gw, mp2_env%ri_g0w0%ic_corr_list, &
1699 weights_cos_tf_t_to_w, weights_sin_tf_t_to_w, &
1700 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, fm_mo_coeff_occ, &
1701 fm_mo_coeff_virt, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, &
1702 mo_coeff(1), fm_mat_w, para_env, para_env_rpa, mat_dm, mat_minvvminv, &
1703 t_3c_o, t_3c_m, t_3c_overl_int_ao_mo, t_3c_o_compressed, t_3c_o_mo_compressed, &
1704 t_3c_o_ind, t_3c_o_mo_ind, &
1705 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1706 matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, mat_w, matrix_s, &
1707 kpoints, mp2_env, qs_env, nkp_self_energy, do_kpoints_cubic_rpa, &
1708 starts_array_mc, ends_array_mc)
1711 IF (exit_ev_gw)
EXIT
1717 IF (do_ic_model)
THEN
1719 IF (my_open_shell)
THEN
1721 CALL calculate_ic_correction(eigenval(:, 1, 1), mat_minvvminv%matrix, &
1722 t_3c_overl_nnp_ic(1), t_3c_overl_nnp_ic_reflected(1), &
1724 gw_corr_lev_occ(1), gw_corr_lev_virt(1), homo(1), unit_nr, &
1725 print_ic_values, para_env, do_alpha=.true.)
1727 CALL calculate_ic_correction(eigenval(:, 1, 2), mat_minvvminv%matrix, &
1728 t_3c_overl_nnp_ic(2), t_3c_overl_nnp_ic_reflected(2), &
1730 gw_corr_lev_occ(2), gw_corr_lev_virt(2), homo(2), unit_nr, &
1731 print_ic_values, para_env, do_beta=.true.)
1735 CALL calculate_ic_correction(eigenval(:, 1, 1), mat_minvvminv%matrix, &
1736 t_3c_overl_nnp_ic(1), t_3c_overl_nnp_ic_reflected(1), &
1738 gw_corr_lev_occ(1), gw_corr_lev_virt(1), homo(1), unit_nr, &
1739 print_ic_values, para_env)
1749 IF (mp2_env%ri_g0w0%iter_evGW > 1)
THEN
1750 IF (unit_nr > 0)
THEN
1751 CALL cp_warn(__location__, &
1752 "BSE@evGW applies W0, i.e. screening with DFT energies to the BSE!")
1756 CALL cp_fm_create(fm_mat_s_ia_bse, fm_mat_s(1)%matrix_struct)
1757 CALL cp_fm_to_fm(fm_mat_s(1), fm_mat_s_ia_bse)
1759 IF (iter_sc_gw0 == 1)
THEN
1760 CALL remove_scaling_factor_rpa(fm_mat_s_ia_bse, virtual(1), &
1761 eigenval_last(:, 1, 1), homo(1), omega)
1763 CALL remove_scaling_factor_rpa(fm_mat_s_ia_bse, virtual(1), &
1764 eigenval_scf(:, 1, 1), homo(1), omega)
1767 CALL start_bse_calculation(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, &
1768 fm_mat_q_static_bse_gemm, &
1769 eigenval, eigenval_scf, &
1770 homo, virtual, dimen_ri, dimen_ri_red, bse_lev_virt, &
1771 gd_array, color_sub, mp2_env, qs_env, mo_coeff, unit_nr)
1773 CALL cp_fm_release(fm_mat_s_ia_bse)
1777 CALL deallocate_matrices_gw(fm_mat_s_gw_work, vec_w_gw, vec_sigma_c_gw, vec_omega_fit_gw, &
1778 mp2_env%ri_g0w0%vec_Sigma_x_minus_vxc_gw, &
1779 eigenval_last, eigenval_scf, do_periodic, matrix_berry_re_mo_mo, matrix_berry_im_mo_mo, &
1780 kpoints, vec_sigma_x_gw,.NOT. do_im_time)
1783 IF (do_im_time)
THEN
1785 CALL dealloc_im_time(fm_mo_coeff_occ, fm_mo_coeff_virt, &
1786 fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, index_to_cell_3c, &
1787 cell_to_index_3c, do_ic_model, &
1788 do_kpoints_cubic_rpa, do_kpoints_from_gamma, do_ri_sigma_x, &
1790 wkp_w, cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1791 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, fm_mat_ri_global_work, fm_mat_work, &
1792 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, mat_dm, mat_l, &
1793 mat_minvvminv, mat_p_omega, mat_p_omega_kp, &
1794 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, mat_work, qs_env)
1797 CALL deallocate_matrices_gw_im_time(weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, do_ic_model, &
1798 do_kpoints_cubic_rpa, fm_mat_w, &
1799 t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, &
1800 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1801 t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, &
1807 IF (.NOT. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2)
CALL exchange_work%release()
1809 IF (.NOT. do_ri_sos_laplace_mp2)
THEN
1812 DEALLOCATE (trace_qomega)
1815 IF (do_im_time .OR. do_ri_sos_laplace_mp2)
THEN
1820 IF (do_im_time .AND. calc_forces)
THEN
1821 CALL im_time_force_release(force_data)
1824 IF (calc_forces .AND. .NOT. do_im_time)
CALL rpa_grad_finalize(
rpa_grad, mp2_env, para_env_sub, para_env, &
1825 qs_env, gd_array, color_sub, do_ri_sos_laplace_mp2, &
1828 CALL timestop(handle)
1830 END SUBROUTINE rpa_num_int
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public wilhelm2016b
integer, save, public delben2015
integer, save, public wilhelm2016a
integer, save, public wilhelm2018
integer, save, public delben2013
integer, save, public ren2013
integer, save, public wilhelm2017
integer, save, public ren2011
integer, save, public gruneis2009
integer, save, public freeman1977
integer, save, public bates2013
Main routines for GW + Bethe-Salpeter for computing electronic excitations.
subroutine, public start_bse_calculation(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, fm_mat_q_static_bse_gemm, eigenval, eigenval_scf, homo, virtual, dimen_ri, dimen_ri_red, bse_lev_virt, gd_array, color_sub, mp2_env, qs_env, mo_coeff, unit_nr)
Main subroutine managing BSE calculations.
methods related to the blacs parallel environment
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Represents a complex full matrix distributed on many processors.
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
subroutine, public dbcsr_clear(matrix)
...
subroutine, public dbcsr_add(matrix_a, matrix_b, alpha_scalar, beta_scalar)
...
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
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
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
This is the start of a dbt_api, all publically needed functions are exported here....
Counters to determine the performance of parallel DGEMMs.
elemental subroutine, public dgemm_counter_init(dgemm_counter, unit_nr, print_info)
Initialize a dgemm_counter.
subroutine, public dgemm_counter_write(dgemm_counter, para_env)
calculate and print flop rates
Types to describe group distributions.
elemental integer function, public maxsize(this)
...
Types and set/get functions for HFX.
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
Types and basic routines needed for a kpoint calculation.
Machine interface based on Fortran 2003 and POSIX.
subroutine, public m_memory(mem)
Returns the total amount of memory [bytes] in use, if known, zero otherwise.
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
complex(kind=dp), parameter, public z_zero
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
Routines to calculate the minimax coefficients in order to approximate 1/x as a sum over exponential ...
subroutine, public check_exp_minimax_range(k, rc, ierr)
Check that a minimax approximation is available for given input k, Rc. ierr == 0: everything ok ierr ...
Routines to calculate frequency and time grids (integration points and weights) for correlation metho...
subroutine, public get_minimax_grid(para_env, unit_nr, homo, eigenval, num_integ_points, do_im_time, do_ri_sos_laplace_mp2, do_print, tau_tj, tau_wj, qs_env, do_gw_im_time, do_kpoints_cubic_rpa, e_fermi, tj, wj, weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, regularization)
...
subroutine, public get_clenshaw_grid(para_env, para_env_rpa, unit_nr, homo, virtual, eigenval, num_integ_points, num_integ_group, color_rpa_group, fm_mat_s, my_do_gw, ext_scaling, a_scaling, tj, wj)
...
Routines to calculate MP2 energy with laplace approach.
subroutine, public sos_mp2_postprocessing(fm_mat_q, erpa, tau_wjquad)
...
Routines for calculating RI-MP2 gradients.
subroutine, public array2fm(mat2d, fm_struct, my_start_row, my_end_row, my_start_col, my_end_col, gd_row, gd_col, group_grid_2_mepos, ngroup_row, ngroup_col, fm_mat, integ_group_size, color_group, do_release_mat)
redistribute local part of array to fm
Types needed for MP2 calculations.
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.
Auxiliary routines needed for RPA-exchange given blacs_env to another.
subroutine, public rpa_exchange_needed_mem(mp2_env, homo, virtual, dimen_ri, para_env, mem_per_rank, mem_per_repl)
...
Routines to calculate RI-RPA and SOS-MP2 gradients.
subroutine, public rpa_grad_finalize(rpa_grad, mp2_env, para_env_sub, para_env, qs_env, gd_array, color_sub, do_ri_sos_laplace_mp2, homo, virtual)
...
subroutine, public rpa_grad_copy_q(fm_mat_q, rpa_grad)
...
pure subroutine, public rpa_grad_needed_mem(homo, virtual, dimen_ri, mem_per_rank, mem_per_repl, do_ri_sos_laplace_mp2)
Calculates the necessary minimum memory for the Gradient code ion MiB.
subroutine, public rpa_grad_create(rpa_grad, fm_mat_q, fm_mat_s, homo, virtual, mp2_env, eigenval, unit_nr, do_ri_sos_laplace_mp2)
Creates the arrays of a rpa_grad_type.
subroutine, public rpa_grad_matrix_operations(mp2_env, rpa_grad, do_ri_sos_laplace_mp2, fm_mat_q, fm_mat_q_gemm, dgemm_counter, fm_mat_s, omega, homo, virtual, eigenval, weight, unit_nr)
...
Routines to calculate image charge corrections.
subroutine, public calculate_ic_correction(eigenval, mat_sinvvsinv, t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, gw_corr_lev_tot, gw_corr_lev_occ, gw_corr_lev_virt, homo, unit_nr, print_ic_values, para_env, do_alpha, do_beta)
...
Routines treating GW and RPA calculations with kpoints.
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 get_bandstruc_and_k_dependent_mos(qs_env, eigenval_kp)
...
Routines for GW, continuous development [Jan Wilhelm].
subroutine, public deallocate_matrices_gw(fm_mat_s_gw_work, vec_w_gw, vec_sigma_c_gw, vec_omega_fit_gw, vec_sigma_x_minus_vxc_gw, eigenval_last, eigenval_scf, do_periodic, matrix_berry_re_mo_mo, matrix_berry_im_mo_mo, kpoints, vec_sigma_x_gw, my_do_gw)
...
subroutine, public allocate_matrices_gw(vec_sigma_c_gw, color_rpa_group, dimen_nm_gw, gw_corr_lev_occ, gw_corr_lev_virt, homo, nmo, num_integ_group, num_integ_points, unit_nr, gw_corr_lev_tot, num_fit_points, omega_max_fit, do_minimax_quad, do_periodic, do_ri_sigma_x, my_do_gw, first_cycle_periodic_correction, a_scaling, eigenval, tj, vec_omega_fit_gw, vec_sigma_x_gw, delta_corr, eigenval_last, eigenval_scf, vec_w_gw, fm_mat_s_gw, fm_mat_s_gw_work, para_env, mp2_env, kpoints, nkp, nkp_self_energy, do_kpoints_cubic_rpa, do_kpoints_from_gamma)
...
subroutine, public compute_gw_self_energy(vec_sigma_c_gw, dimen_nm_gw, dimen_ri, gw_corr_lev_occ, gw_corr_lev_virt, homo, jquad, nmo, num_fit_points, do_im_time, do_periodic, first_cycle_periodic_correction, fermi_level_offset, omega, eigenval, delta_corr, vec_omega_fit_gw, vec_w_gw, wj, fm_mat_q, fm_mat_r_gw, fm_mat_s_gw, fm_mat_s_gw_work, mo_coeff, para_env, para_env_rpa, matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, kpoints, qs_env, mp2_env)
...
subroutine, public get_fermi_level_offset(fermi_level_offset, fermi_level_offset_input, eigenval, homo)
...
subroutine, public allocate_matrices_gw_im_time(gw_corr_lev_occ, gw_corr_lev_virt, homo, nmo, num_integ_points, unit_nr, ri_blk_sizes, do_ic_model, para_env, fm_mat_w, fm_mat_q, mo_coeff, t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, starts_array_mc, ends_array_mc, t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, matrix_s, mat_w, t_3c_overl_int, t_3c_o_compressed, t_3c_o_ind, qs_env)
...
subroutine, public deallocate_matrices_gw_im_time(weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, do_ic_model, do_kpoints_cubic_rpa, fm_mat_w, t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, mat_w, qs_env)
...
subroutine, public compute_qp_energies(vec_sigma_c_gw, count_ev_sc_gw, gw_corr_lev_occ, gw_corr_lev_tot, gw_corr_lev_virt, homo, nmo, num_fit_points, num_integ_points, unit_nr, do_apply_ic_corr_to_gw, do_im_time, do_periodic, do_ri_sigma_x, first_cycle_periodic_correction, e_fermi, eps_filter, fermi_level_offset, delta_corr, eigenval, eigenval_last, eigenval_scf, iter_sc_gw0, exit_ev_gw, tau_tj, tj, vec_omega_fit_gw, vec_sigma_x_gw, ic_corr_list, weights_cos_tf_t_to_w, weights_sin_tf_t_to_w, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, fm_mo_coeff_occ, fm_mo_coeff_virt, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, mo_coeff, fm_mat_w, para_env, para_env_rpa, mat_dm, mat_minvvminv, t_3c_o, t_3c_m, t_3c_overl_int_ao_mo, t_3c_o_compressed, t_3c_o_mo_compressed, t_3c_o_ind, t_3c_o_mo_ind, t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, mat_w, matrix_s, kpoints, mp2_env, qs_env, nkp_self_energy, do_kpoints_cubic_rpa, starts_array_mc, ends_array_mc)
...
subroutine, public compute_w_cubic_gw(fm_mat_w, fm_mat_q, fm_mat_work, dimen_ri, fm_mat_l, num_integ_points, tj, tau_tj, weights_cos_tf_w_to_t, jquad, omega)
...
Routines needed for cubic-scaling RPA and SOS-Laplace-MP2 forces.
subroutine, public calc_laplace_loop_forces(force_data, mat_p_omega, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, num_integ_points, nmo, eigenval, tau_tj, tau_wj, cut_memory, pspin, qspin, open_shell, unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
Updates the cubic-scaling SOS-Laplace-MP2 contribution to the forces at each quadrature point.
subroutine, public calc_post_loop_forces(force_data, unit_nr, qs_env)
All the forces that can be calculated after the loop on the Laplace quaradture, using terms collected...
subroutine, public calc_rpa_loop_forces(force_data, mat_p_omega, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, num_integ_points, nmo, eigenval, e_fermi, weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, tj, wj, tau_tj, cut_memory, ispin, open_shell, unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
Updates the cubic-scaling RPA contribution to the forces at each quadrature point....
subroutine, public keep_initial_quad(tj, wj, tau_tj, tau_wj, weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, do_laplace, do_im_time, num_integ_points, unit_nr, qs_env)
Overwrites the "optimal" Laplace quadrature with that of the first step.
subroutine, public init_im_time_forces(force_data, fm_matrix_pq, t_3c_m, unit_nr, mp2_env, qs_env)
Initializes and pre-calculates all needed tensors for the forces.
Types needed for cubic-scaling RPA and SOS-Laplace-MP2 forces.
subroutine, public im_time_force_release(force_data)
Cleans everything up.
Routines for low-scaling RPA/GW with imaginary time.
subroutine, public zero_mat_p_omega(mat_p_omega)
...
subroutine, public compute_mat_p_omega(mat_p_omega, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, mat_p_global, matrix_s, ispin, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, weights_cos_tf_t_to_w, tj, tau_tj, e_fermi, eps_filter, alpha, eps_filter_im_time, eigenval, nmo, num_integ_points, cut_memory, unit_nr, mp2_env, para_env, qs_env, do_kpoints_from_gamma, index_to_cell_3c, cell_to_index_3c, has_mat_p_blocks, do_ri_sos_laplace_mp2, dbcsr_time, dbcsr_nflop)
...
Routines to calculate RI-RPA energy.
subroutine, public rpa_ri_compute_en(qs_env, erpa, mp2_env, bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, para_env, para_env_sub, color_sub, gd_array, gd_b_virtual, gd_b_all, gd_b_occ_bse, gd_b_virt_bse, mo_coeff, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, kpoints, eigenval, nmo, homo, dimen_ri, dimen_ri_red, gw_corr_lev_occ, gw_corr_lev_virt, bse_lev_virt, unit_nr, do_ri_sos_laplace_mp2, my_do_gw, do_im_time, do_bse, matrix_s, mat_munu, mat_p_global, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, calc_forces)
...
Routines to calculate RI-RPA energy and Sigma correction to the RPA energies using the cubic spline b...
subroutine, public rpa_sigma_create(rpa_sigma, sigma_param, fm_mat_q, unit_nr, para_env)
... Collect the Q(w) (fm_mat_Q) matrix to create rpa_sigma a derived type variable....
subroutine, public finalize_rpa_sigma(rpa_sigma, unit_nr, e_sigma_corr, para_env, do_minimax_quad)
... Save the calculated value of E_c correction to the global variable and memory clean.
subroutine, public rpa_sigma_matrix_spectral(rpa_sigma, fm_mat_q, wj, para_env_rpa)
... Diagonalize and store the eigenvalues of fm_mat_Q in rpa_sigmasigma_eigenvalue.
Utility functions for RPA calculations.
subroutine, public compute_erpa_by_freq_int(dimen_ri, trace_qomega, fm_mat_q, para_env_rpa, erpa, wjquad)
...
subroutine, public alloc_im_time(qs_env, para_env, dimen_ri, dimen_ri_red, num_integ_points, nspins, fm_mat_q, fm_mo_coeff_occ, fm_mo_coeff_virt, fm_matrix_minv_l_kpoints, fm_matrix_l_kpoints, mat_p_global, t_3c_o, matrix_s, kpoints, eps_filter_im_time, cut_memory, nkp, num_cells_dm, num_3c_repl, size_p, ikp_local, index_to_cell_3c, cell_to_index_3c, col_blk_size, do_ic_model, do_kpoints_cubic_rpa, do_kpoints_from_gamma, do_ri_sigma_x, my_open_shell, has_mat_p_blocks, wkp_w, cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, fm_mat_ri_global_work, fm_mat_work, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, mat_dm, mat_l, mat_m_p_munu_occ, mat_m_p_munu_virt, mat_minvvminv, mat_p_omega, mat_p_omega_kp, mat_work, mo_coeff, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, homo, nmo)
...
subroutine, public dealloc_im_time(fm_mo_coeff_occ, fm_mo_coeff_virt, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, index_to_cell_3c, cell_to_index_3c, do_ic_model, do_kpoints_cubic_rpa, do_kpoints_from_gamma, do_ri_sigma_x, has_mat_p_blocks, wkp_w, cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, fm_mat_ri_global_work, fm_mat_work, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, mat_dm, mat_l, mat_minvvminv, mat_p_omega, mat_p_omega_kp, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, mat_work, qs_env)
...
subroutine, public q_trace_and_add_unit_matrix(dimen_ri, trace_qomega, fm_mat_q)
...
subroutine, public calc_mat_q(fm_mat_s, do_ri_sos_laplace_mp2, first_cycle, virtual, eigenval, homo, omega, omega_old, jquad, mm_style, dimen_ri, dimen_ia, alpha, fm_mat_q, fm_mat_q_gemm, do_bse, fm_mat_q_static_bse_gemm, dgemm_counter, num_integ_points, count_ev_sc_gw)
...
subroutine, public contract_p_omega_with_mat_l(mat_p_omega, mat_l, mat_work, eps_filter_im_time, fm_mat_work, dimen_ri, dimen_ri_red, fm_mat_l, fm_mat_q)
...
subroutine, public remove_scaling_factor_rpa(fm_mat_s, virtual, eigenval_last, homo, omega_old)
...
All kind of helpful little routines.
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
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
Contains information about kpoints.
stores all the informations relevant to an mpi environment