(git:12618bd)
Loading...
Searching...
No Matches
rpa_main.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines to calculate RI-RPA energy
10!> \par History
11!> 06.2012 created [Mauro Del Ben]
12!> 04.2015 GW routines added [Jan Wilhelm]
13!> 10.2015 Cubic-scaling RPA routines added [Jan Wilhelm]
14!> 10.2018 Cubic-scaling SOS-MP2 added [Frederick Stein]
15!> 03.2019 Refactoring [Frederick Stein]
16! **************************************************************************************************
18 USE bibliography, ONLY: &
25 USE cp_cfm_types, ONLY: cp_cfm_type
26 USE cp_dbcsr_api, ONLY: dbcsr_add,&
35 USE cp_fm_types, ONLY: cp_fm_create,&
40 USE dbt_api, ONLY: dbt_type
47 maxsize,&
49 USE hfx_types, ONLY: block_ind_type,&
56 USE kinds, ONLY: dp,&
57 int_8
58 USE kpoint_types, ONLY: kpoint_type
59 USE machine, ONLY: m_flush,&
61 USE mathconstants, ONLY: pi,&
62 z_zero
63 USE message_passing, ONLY: mp_comm_type,&
67 USE mp2_grids, ONLY: get_clenshaw_grid,&
70 USE mp2_ri_grad_util, ONLY: array2fm
71 USE mp2_types, ONLY: mp2_type,&
79 USE rpa_grad, ONLY: rpa_grad_copy_q,&
85 USE rpa_gw, ONLY: allocate_matrices_gw,&
111 calc_mat_q,&
116 USE util, ONLY: get_limit
117#include "./base/base_uses.f90"
118
119 IMPLICIT NONE
120
121 PRIVATE
122
123 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rpa_main'
124
125 PUBLIC :: rpa_ri_compute_en
126
127CONTAINS
128
129! **************************************************************************************************
130!> \brief ...
131!> \param qs_env ...
132!> \param Erpa ...
133!> \param mp2_env ...
134!> \param BIb_C ...
135!> \param BIb_C_gw ...
136!> \param BIb_C_bse_ij ...
137!> \param BIb_C_bse_ab ...
138!> \param para_env ...
139!> \param para_env_sub ...
140!> \param color_sub ...
141!> \param gd_array ...
142!> \param gd_B_virtual ...
143!> \param gd_B_all ...
144!> \param gd_B_occ_bse ...
145!> \param gd_B_virt_bse ...
146!> \param mo_coeff ...
147!> \param fm_matrix_PQ ...
148!> \param fm_matrix_L_kpoints ...
149!> \param fm_matrix_Minv_L_kpoints ...
150!> \param fm_matrix_Minv ...
151!> \param fm_matrix_Minv_Vtrunc_Minv ...
152!> \param kpoints ...
153!> \param Eigenval ...
154!> \param nmo ...
155!> \param homo ...
156!> \param dimen_RI ...
157!> \param dimen_RI_red ...
158!> \param gw_corr_lev_occ ...
159!> \param gw_corr_lev_virt ...
160!> \param bse_lev_virt ...
161!> \param unit_nr ...
162!> \param do_ri_sos_laplace_mp2 ...
163!> \param my_do_gw ...
164!> \param do_im_time ...
165!> \param do_bse ...
166!> \param matrix_s ...
167!> \param mat_munu ...
168!> \param mat_P_global ...
169!> \param t_3c_M ...
170!> \param t_3c_O ...
171!> \param t_3c_O_compressed ...
172!> \param t_3c_O_ind ...
173!> \param starts_array_mc ...
174!> \param ends_array_mc ...
175!> \param starts_array_mc_block ...
176!> \param ends_array_mc_block ...
177!> \param calc_forces ...
178! **************************************************************************************************
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, &
185 bse_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)
190
191 TYPE(qs_environment_type), POINTER :: qs_env
192 REAL(kind=dp), INTENT(OUT) :: erpa
193 TYPE(mp2_type), INTENT(INOUT) :: mp2_env
194 TYPE(three_dim_real_array), DIMENSION(:), &
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
198 TYPE(mp_para_env_type), POINTER :: para_env, para_env_sub
199 INTEGER, INTENT(INOUT) :: color_sub
200 TYPE(group_dist_d1_type), INTENT(INOUT) :: gd_array
201 TYPE(group_dist_d1_type), DIMENSION(:), &
202 INTENT(INOUT) :: gd_b_virtual
203 TYPE(group_dist_d1_type), INTENT(INOUT) :: gd_b_all, gd_b_occ_bse, gd_b_virt_bse
204 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
205 TYPE(cp_fm_type), INTENT(IN) :: fm_matrix_pq
206 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_l_kpoints, &
207 fm_matrix_minv_l_kpoints, &
208 fm_matrix_minv, &
209 fm_matrix_minv_vtrunc_minv
210 TYPE(kpoint_type), POINTER :: kpoints
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, &
219 do_im_time, do_bse
220 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
221 TYPE(dbcsr_p_type), INTENT(IN) :: mat_munu
222 TYPE(dbcsr_p_type), INTENT(INOUT) :: mat_p_global
223 TYPE(dbt_type) :: t_3c_m
224 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :) :: t_3c_o
225 TYPE(hfx_compression_type), ALLOCATABLE, &
226 DIMENSION(:, :, :), INTENT(INOUT) :: t_3c_o_compressed
227 TYPE(block_ind_type), ALLOCATABLE, &
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, &
231 ends_array_mc_block
232 LOGICAL, INTENT(IN) :: calc_forces
233
234 CHARACTER(LEN=*), PARAMETER :: routinen = 'rpa_ri_compute_en'
235
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, &
245 my_ia_start, virtual
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, &
252 fm_mat_s_gw
253 TYPE(cp_fm_type), DIMENSION(1) :: fm_mat_r_gw, fm_mat_s_ab_bse, &
254 fm_mat_s_ij_bse
255 TYPE(mp_para_env_type), POINTER :: para_env_rpa
256 TYPE(two_dim_real_array), ALLOCATABLE, &
257 DIMENSION(:) :: bib_c_2d, bib_c_2d_gw
258 TYPE(two_dim_real_array), DIMENSION(1) :: bib_c_2d_bse_ab, bib_c_2d_bse_ij
259
260 CALL timeset(routinen, handle)
261
262 CALL cite_reference(delben2013)
263 CALL cite_reference(delben2015)
264
265 IF (mp2_env%ri_rpa%exchange_correction == rpa_exchange_axk) THEN
266 CALL cite_reference(bates2013)
267 ELSE IF (mp2_env%ri_rpa%exchange_correction == rpa_exchange_sosex) THEN
268 CALL cite_reference(freeman1977)
269 CALL cite_reference(gruneis2009)
270 END IF
271 IF (mp2_env%ri_rpa%do_rse) THEN
272 CALL cite_reference(ren2011)
273 CALL cite_reference(ren2013)
274 END IF
275
276 IF (my_do_gw) THEN
277 CALL cite_reference(wilhelm2016a)
278 CALL cite_reference(wilhelm2017)
279 CALL cite_reference(wilhelm2018)
280 END IF
281
282 IF (do_im_time) THEN
283 CALL cite_reference(wilhelm2016b)
284 END IF
285
286 nspins = SIZE(homo)
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(:)
291
292 ALLOCATE (eigenval_kp(nmo, 1, nspins))
293 eigenval_kp(:, 1, :) = eigenval(:, :)
294
295 IF (do_im_time) mp2_env%ri_rpa%minimax_quad = .true.
296 do_minimax_quad = mp2_env%ri_rpa%minimax_quad
297
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
301
302 ! check the range for the minimax approximation
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
305 emin = huge(dp)
306 emax = 0.0_dp
307 DO ispin = 1, nspins
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))))
311 END IF
312 END DO
313 e_range = emax/emin
314 END IF
315 IF (e_range < 2.0_dp) e_range = 2.0_dp
316 ierr = 0
317 CALL check_exp_minimax_range(num_integ_points, e_range, ierr)
318 IF (ierr /= 0) THEN
319 jjb = num_integ_points - 1
320 DO iib = 1, jjb
321 num_integ_points = num_integ_points - 1
322 ierr = 0
323 CALL check_exp_minimax_range(num_integ_points, e_range, ierr)
324 IF (ierr == 0) EXIT
325 END DO
326 END IF
327 cpassert(num_integ_points >= 1)
328 ELSE
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
333 IF (unit_nr > 0) &
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
338 END IF
339 ELSE
340 IF (do_minimax_quad .AND. num_integ_points > 20) THEN
341 IF (unit_nr > 0) &
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
346 END IF
347 END IF
348 END IF
349 allowed_memory = mp2_env%mp2_memory
350
351 CALL get_group_dist(gd_array, color_sub, my_group_l_start, my_group_l_end, my_group_l_size)
352
353 ngroup = para_env%num_pe/para_env_sub%num_pe
354
355 ! for imaginary time or periodic GW or BSE, we use all processors for a single frequency/time point
356 IF (do_im_time .OR. mp2_env%ri_g0w0%do_periodic .OR. do_bse) THEN
357
358 integ_group_size = ngroup
359 best_num_integ_point = num_integ_points
360
361 ELSE
362
363 ! Calculate available memory and create integral group according to that
364 ! mem_for_iaK is the memory needed for storing the 3 centre integrals
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)
367
368 CALL m_memory(mem)
369 mem_real = (mem + 1024*1024 - 1)/(1024*1024)
370 CALL para_env%min(mem_real)
371
372 mem_per_rank = 0.0_dp
373
374 ! B_ia_P
375 mem_per_repl = mem_for_iak
376 ! Q (regular and for dgemm)
377 mem_per_repl = mem_per_repl + 2.0_dp*mem_for_qk
378
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)
380 CALL rpa_exchange_needed_mem(mp2_env, homo, virtual, dimen_ri_red, para_env, mem_per_rank, mem_per_repl)
381
382 mem_min = mem_per_repl/para_env%num_pe + mem_per_rank
383
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'
387 END IF
388
389 ! Use only the allowed amount of memory
390 mem_real = min(mem_real, allowed_memory)
391 ! For the memory estimate, we require the amount of required memory per replication group and the available memory
392 mem_real = mem_real - mem_per_rank
393
394 mem_per_group = mem_real*para_env_sub%num_pe
395
396 ! here we try to find the best rpa/laplace group size
397 skip_integ_group_opt = .false.
398
399 ! Check the input number of integration groups
400 IF (input_num_integ_groups > 0) THEN
401 IF (num_integ_points < input_num_integ_groups) 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 - 1)/input_num_integ_groups
405 skip_integ_group_opt = .true.
406 ELSE
407 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Total number of groups not multiple of NUM_INTEG_GROUPS'
408 END IF
409 ELSE
410 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Too many integration groups for the given number of quadrature points'
411 END IF
412 END IF
413
414 IF (.NOT. skip_integ_group_opt) THEN
415 best_integ_group_size = ngroup
416 best_num_integ_point = num_integ_points
417
418 min_integ_group_size = max(1, ngroup/num_integ_points)
419
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
423
424 ! check that the ngroup is a multiple of integ_group_size
425 IF (mod(ngroup, integ_group_size) /= 0) cycle
426
427 ! check for memory
428 avail_mem = integ_group_size*mem_per_group
429 IF (avail_mem < mem_per_repl) cycle
430
431 ! check that the integration groups have the same size
432 num_integ_group = ngroup/integ_group_size
433
434 best_num_integ_point = (num_integ_points + num_integ_group - 1)/num_integ_group
435 best_integ_group_size = integ_group_size
436
437 EXIT
438
439 END DO
440 END IF
441
442 integ_group_size = best_integ_group_size
443
444 END IF
445
446 IF (unit_nr > 0 .AND. .NOT. do_im_time) THEN
447 IF (do_ri_sos_laplace_mp2) THEN
448 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
449 "RI_INFO| Group size for laplace numerical integration:", integ_group_size*para_env_sub%num_pe
450 WRITE (unit=unit_nr, fmt="(T3,A)") &
451 "INTEG_INFO| MINIMAX approximation"
452 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
453 "INTEG_INFO| Number of integration points:", num_integ_points
454 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
455 "INTEG_INFO| Max. number of integration points per Laplace group:", best_num_integ_point
456 ELSE
457 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
458 "RI_INFO| Group size for frequency integration:", integ_group_size*para_env_sub%num_pe
459 IF (do_minimax_quad) THEN
460 WRITE (unit=unit_nr, fmt="(T3,A)") &
461 "INTEG_INFO| MINIMAX quadrature"
462 ELSE
463 WRITE (unit=unit_nr, fmt="(T3,A)") &
464 "INTEG_INFO| Clenshaw-Curtius quadrature"
465 END IF
466 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
467 "INTEG_INFO| Number of integration points:", num_integ_points
468 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
469 "INTEG_INFO| Max. number of integration points per RPA group:", best_num_integ_point
470 END IF
471 CALL m_flush(unit_nr)
472 END IF
473
474 num_integ_group = ngroup/integ_group_size
475
476 pos_integ_group = mod(color_sub, integ_group_size)
477 color_rpa_group = color_sub/integ_group_size
478
479 CALL timeset(routinen//"_reorder", handle2)
480
481 ! not necessary for imaginary time
482
483 ALLOCATE (bib_c_2d(nspins))
484
485 IF (.NOT. do_im_time) THEN
486
487 ! reorder the local data in such a way to help the next stage of matrix creation
488 ! now the data inside the group are divided into a ia x K matrix
489 DO ispin = 1, nspins
490 CALL calculate_bib_c_2d(bib_c_2d(ispin)%array, bib_c(ispin)%array, para_env_sub, dimen_ia(ispin), &
491 homo(ispin), virtual(ispin), gd_b_virtual(ispin), &
492 my_ia_size(ispin), my_ia_start(ispin), my_ia_end(ispin), my_group_l_size)
493
494 DEALLOCATE (bib_c(ispin)%array)
495 CALL release_group_dist(gd_b_virtual(ispin))
496
497 END DO
498
499 ! in the GW case, BIb_C_2D_gw is an nm x K matrix, with n: number of corr GW levels, m=nmo
500 IF (my_do_gw) THEN
501 ALLOCATE (bib_c_2d_gw(nspins))
502
503 CALL timeset(routinen//"_reorder_gw", handle3)
504
505 dimen_nm_gw = nmo*(gw_corr_lev_occ(1) + gw_corr_lev_virt(1))
506
507 ! The same for open shell
508 DO ispin = 1, nspins
509 CALL calculate_bib_c_2d(bib_c_2d_gw(ispin)%array, bib_c_gw(ispin)%array, para_env_sub, dimen_nm_gw, &
510 gw_corr_lev_occ(ispin) + gw_corr_lev_virt(ispin), nmo, gd_b_all, &
511 my_nm_gw_size, my_nm_gw_start, my_nm_gw_end, my_group_l_size)
512 DEALLOCATE (bib_c_gw(ispin)%array)
513 END DO
514
515 CALL release_group_dist(gd_b_all)
516
517 CALL timestop(handle3)
518
519 END IF
520 END IF
521
522 IF (do_bse) THEN
523
524 CALL timeset(routinen//"_reorder_bse1", handle3)
525
526 dimen_homo_square = homo(1)**2
527 ! We do not implement an explicit bse_lev_occ different to homo here, because the small number of occupied levels
528 ! does not critically influence the memory
529 CALL calculate_bib_c_2d(bib_c_2d_bse_ij(1)%array, bib_c_bse_ij, para_env_sub, dimen_homo_square, &
530 homo(1), homo(1), gd_b_occ_bse, &
531 my_ij_comb_bse_size, my_ij_comb_bse_start, my_ij_comb_bse_end, my_group_l_size)
532
533 DEALLOCATE (bib_c_bse_ij)
534 CALL release_group_dist(gd_b_occ_bse)
535
536 CALL timestop(handle3)
537
538 CALL timeset(routinen//"_reorder_bse2", handle3)
539
540 dimen_virt_square = bse_lev_virt**2
541
542 CALL calculate_bib_c_2d(bib_c_2d_bse_ab(1)%array, bib_c_bse_ab, para_env_sub, dimen_virt_square, &
543 bse_lev_virt, bse_lev_virt, gd_b_virt_bse, &
544 my_ab_comb_bse_size, my_ab_comb_bse_start, my_ab_comb_bse_end, my_group_l_size)
545
546 DEALLOCATE (bib_c_bse_ab)
547 CALL release_group_dist(gd_b_virt_bse)
548
549 CALL timestop(handle3)
550
551 END IF
552
553 CALL timestop(handle2)
554
555 IF (num_integ_group > 1) THEN
556 ALLOCATE (para_env_rpa)
557 CALL para_env_rpa%from_split(para_env, color_rpa_group)
558 ELSE
559 para_env_rpa => para_env
560 END IF
561
562 ! now create the matrices needed for the calculation, Q, S and G
563 ! Q and G will have omega dependence
564
565 IF (do_im_time) THEN
566 ALLOCATE (fm_mat_q(nspins), fm_mat_q_gemm(1), fm_mat_s(1))
567 ELSE
568 ALLOCATE (fm_mat_q(nspins), fm_mat_q_gemm(nspins), fm_mat_s(nspins))
569 END IF
570
571 CALL create_integ_mat(bib_c_2d, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
572 dimen_ri_red, dimen_ia, color_rpa_group, &
573 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
574 my_ia_size, my_ia_start, my_ia_end, &
575 my_group_l_size, my_group_l_start, my_group_l_end, &
576 para_env_rpa, fm_mat_s, nrow_block_mat, ncol_block_mat, &
577 dimen_ia_for_block_size=dimen_ia(1), &
578 do_im_time=do_im_time, fm_mat_q_gemm=fm_mat_q_gemm, fm_mat_q=fm_mat_q, qs_env=qs_env)
579
580 DEALLOCATE (bib_c_2d, my_ia_end, my_ia_size, my_ia_start)
581
582 ! for GW, we need other matrix fm_mat_S, always allocate the container to prevent crying compilers
583 ALLOCATE (fm_mat_s_gw(nspins))
584 IF (my_do_gw .AND. .NOT. do_im_time) THEN
585
586 CALL create_integ_mat(bib_c_2d_gw, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
587 dimen_ri_red, [dimen_nm_gw, dimen_nm_gw], color_rpa_group, &
588 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
589 [my_nm_gw_size, my_nm_gw_size], [my_nm_gw_start, my_nm_gw_start], [my_nm_gw_end, my_nm_gw_end], &
590 my_group_l_size, my_group_l_start, my_group_l_end, &
591 para_env_rpa, fm_mat_s_gw, nrow_block_mat, ncol_block_mat, &
592 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context, &
593 fm_mat_q=fm_mat_r_gw)
594 DEALLOCATE (bib_c_2d_gw)
595
596 END IF
597
598 ! for Bethe-Salpeter, we need other matrix fm_mat_S
599 IF (do_bse) THEN
600 CALL create_integ_mat(bib_c_2d_bse_ij, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
601 dimen_ri_red, [dimen_homo_square], color_rpa_group, &
602 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
603 [my_ij_comb_bse_size], [my_ij_comb_bse_start], [my_ij_comb_bse_end], &
604 my_group_l_size, my_group_l_start, my_group_l_end, &
605 para_env_rpa, fm_mat_s_ij_bse, nrow_block_mat, ncol_block_mat, &
606 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context)
607
608 CALL create_integ_mat(bib_c_2d_bse_ab, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
609 dimen_ri_red, [dimen_virt_square], color_rpa_group, &
610 mp2_env%block_size_row, mp2_env%block_size_col, unit_nr, &
611 [my_ab_comb_bse_size], [my_ab_comb_bse_start], [my_ab_comb_bse_end], &
612 my_group_l_size, my_group_l_start, my_group_l_end, &
613 para_env_rpa, fm_mat_s_ab_bse, nrow_block_mat, ncol_block_mat, &
614 fm_mat_q(1)%matrix_struct%context, fm_mat_q(1)%matrix_struct%context)
615
616 END IF
617
618 do_kpoints_from_gamma = qs_env%mp2_env%ri_rpa_im_time%do_kpoints_from_Gamma
619 IF (do_kpoints_from_gamma) THEN
620 CALL get_bandstruc_and_k_dependent_mos(qs_env, eigenval_kp)
621 END IF
622
623 ! Now start the RPA calculation
624 ! fm_mo_coeff_occ, fm_mo_coeff_virt will be deallocated here
625 CALL rpa_num_int(qs_env, erpa, mp2_env, para_env, para_env_rpa, para_env_sub, unit_nr, &
626 homo, virtual, dimen_ri, dimen_ri_red, dimen_ia, dimen_nm_gw, &
627 eigenval_kp, num_integ_points, num_integ_group, color_rpa_group, &
628 fm_matrix_pq, fm_mat_s, fm_mat_q_gemm, fm_mat_q, fm_mat_s_gw, fm_mat_r_gw(1), &
629 fm_mat_s_ij_bse(1), fm_mat_s_ab_bse(1), &
630 my_do_gw, do_bse, gw_corr_lev_occ, gw_corr_lev_virt, &
631 bse_lev_virt, &
632 do_minimax_quad, &
633 do_im_time, mo_coeff, &
634 fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, &
635 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, mat_munu, mat_p_global, &
636 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
637 starts_array_mc, ends_array_mc, &
638 starts_array_mc_block, ends_array_mc_block, &
639 matrix_s, do_kpoints_from_gamma, kpoints, gd_array, color_sub, &
640 do_ri_sos_laplace_mp2=do_ri_sos_laplace_mp2, calc_forces=calc_forces)
641
642 CALL release_group_dist(gd_array)
643
644 IF (num_integ_group > 1) CALL mp_para_env_release(para_env_rpa)
645
646 IF (.NOT. do_im_time) THEN
647 CALL cp_fm_release(fm_mat_q_gemm)
648 CALL cp_fm_release(fm_mat_s)
649 END IF
650 CALL cp_fm_release(fm_mat_q)
651
652 IF (my_do_gw .AND. .NOT. do_im_time) THEN
653 CALL cp_fm_release(fm_mat_s_gw)
654 CALL cp_fm_release(fm_mat_r_gw(1))
655 END IF
656
657 IF (do_bse) THEN
658 CALL cp_fm_release(fm_mat_s_ij_bse(1))
659 CALL cp_fm_release(fm_mat_s_ab_bse(1))
660 END IF
661
662 CALL timestop(handle)
663
664 END SUBROUTINE rpa_ri_compute_en
665
666! **************************************************************************************************
667!> \brief reorder the local data in such a way to help the next stage of matrix creation;
668!> now the data inside the group are divided into a ia x K matrix (BIb_C_2D);
669!> Subroutine created to avoid massive double coding
670!> \param BIb_C_2D ...
671!> \param BIb_C ...
672!> \param para_env_sub ...
673!> \param dimen_ia ...
674!> \param homo ...
675!> \param virtual ...
676!> \param gd_B_virtual ...
677!> \param my_ia_size ...
678!> \param my_ia_start ...
679!> \param my_ia_end ...
680!> \param my_group_L_size ...
681!> \author Jan Wilhelm, 03/2015
682! **************************************************************************************************
683 SUBROUTINE calculate_bib_c_2d(BIb_C_2D, BIb_C, para_env_sub, dimen_ia, homo, virtual, &
684 gd_B_virtual, &
685 my_ia_size, my_ia_start, my_ia_end, my_group_L_size)
686
687 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
688 INTENT(OUT) :: bib_c_2d
689 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :), &
690 INTENT(IN) :: bib_c
691 TYPE(mp_para_env_type), INTENT(IN) :: para_env_sub
692 INTEGER, INTENT(IN) :: dimen_ia, homo, virtual
693 TYPE(group_dist_d1_type), INTENT(INOUT) :: gd_b_virtual
694 INTEGER :: my_ia_size, my_ia_start, my_ia_end, &
695 my_group_l_size
696
697 INTEGER, PARAMETER :: occ_chunk = 128
698
699 INTEGER :: ia_global, iib, itmp(2), jjb, my_b_size, my_b_virtual_start, occ_high, occ_low, &
700 proc_receive, proc_send, proc_shift, rec_b_size, rec_b_virtual_end, rec_b_virtual_start
701 REAL(kind=dp), ALLOCATABLE, DIMENSION(:), TARGET :: bib_c_rec_1d
702 REAL(kind=dp), DIMENSION(:, :, :), POINTER :: bib_c_rec
703
704 itmp = get_limit(dimen_ia, para_env_sub%num_pe, para_env_sub%mepos)
705 my_ia_start = itmp(1)
706 my_ia_end = itmp(2)
707 my_ia_size = my_ia_end - my_ia_start + 1
708
709 CALL get_group_dist(gd_b_virtual, para_env_sub%mepos, sizes=my_b_size, starts=my_b_virtual_start)
710
711 ! reorder data
712 ALLOCATE (bib_c_2d(my_group_l_size, my_ia_size))
713
714!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,ia_global) &
715!$OMP SHARED(homo,my_B_size,virtual,my_B_virtual_start,my_ia_start,my_ia_end,BIb_C,BIb_C_2D,&
716!$OMP my_group_L_size)
717 DO iib = 1, homo
718 DO jjb = 1, my_b_size
719 ia_global = (iib - 1)*virtual + my_b_virtual_start + jjb - 1
720 IF (ia_global >= my_ia_start .AND. ia_global <= my_ia_end) THEN
721 bib_c_2d(1:my_group_l_size, ia_global - my_ia_start + 1) = bib_c(1:my_group_l_size, jjb, iib)
722 END IF
723 END DO
724 END DO
725
726 IF (para_env_sub%num_pe > 1) THEN
727 ALLOCATE (bib_c_rec_1d(int(my_group_l_size, int_8)*maxsize(gd_b_virtual)*min(homo, occ_chunk)))
728 DO proc_shift = 1, para_env_sub%num_pe - 1
729 proc_send = modulo(para_env_sub%mepos + proc_shift, para_env_sub%num_pe)
730 proc_receive = modulo(para_env_sub%mepos - proc_shift, para_env_sub%num_pe)
731
732 CALL get_group_dist(gd_b_virtual, proc_receive, rec_b_virtual_start, rec_b_virtual_end, rec_b_size)
733
734 ! do this in chunks to avoid high memory overhead
735 DO occ_low = 1, homo, occ_chunk
736 occ_high = min(homo, occ_low + occ_chunk - 1)
737 bib_c_rec(1:my_group_l_size, 1:rec_b_size, 1:occ_high - occ_low + 1) => &
738 bib_c_rec_1d(1:int(my_group_l_size, int_8)*rec_b_size*(occ_high - occ_low + 1))
739 CALL para_env_sub%sendrecv(bib_c(:, :, occ_low:occ_high), proc_send, &
740 bib_c_rec(:, :, 1:occ_high - occ_low + 1), proc_receive)
741!$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,ia_global) &
742!$OMP SHARED(occ_low,occ_high,rec_B_size,virtual,rec_B_virtual_start,my_ia_start,my_ia_end,BIb_C_rec,BIb_C_2D,&
743!$OMP my_group_L_size)
744 DO iib = occ_low, occ_high
745 DO jjb = 1, rec_b_size
746 ia_global = (iib - 1)*virtual + rec_b_virtual_start + jjb - 1
747 IF (ia_global >= my_ia_start .AND. ia_global <= my_ia_end) THEN
748 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)
749 END IF
750 END DO
751 END DO
752 END DO
753
754 END DO
755 DEALLOCATE (bib_c_rec_1d)
756 END IF
757
758 END SUBROUTINE calculate_bib_c_2d
759
760! **************************************************************************************************
761!> \brief ...
762!> \param BIb_C_2D ...
763!> \param para_env ...
764!> \param para_env_sub ...
765!> \param color_sub ...
766!> \param ngroup ...
767!> \param integ_group_size ...
768!> \param dimen_RI ...
769!> \param dimen_ia ...
770!> \param color_rpa_group ...
771!> \param ext_row_block_size ...
772!> \param ext_col_block_size ...
773!> \param unit_nr ...
774!> \param my_ia_size ...
775!> \param my_ia_start ...
776!> \param my_ia_end ...
777!> \param my_group_L_size ...
778!> \param my_group_L_start ...
779!> \param my_group_L_end ...
780!> \param para_env_RPA ...
781!> \param fm_mat_S ...
782!> \param nrow_block_mat ...
783!> \param ncol_block_mat ...
784!> \param blacs_env_ext ...
785!> \param blacs_env_ext_S ...
786!> \param dimen_ia_for_block_size ...
787!> \param do_im_time ...
788!> \param fm_mat_Q_gemm ...
789!> \param fm_mat_Q ...
790!> \param qs_env ...
791! **************************************************************************************************
792 SUBROUTINE create_integ_mat(BIb_C_2D, para_env, para_env_sub, color_sub, ngroup, integ_group_size, &
793 dimen_RI, dimen_ia, color_rpa_group, &
794 ext_row_block_size, ext_col_block_size, unit_nr, &
795 my_ia_size, my_ia_start, my_ia_end, &
796 my_group_L_size, my_group_L_start, my_group_L_end, &
797 para_env_RPA, fm_mat_S, nrow_block_mat, ncol_block_mat, &
798 blacs_env_ext, blacs_env_ext_S, dimen_ia_for_block_size, &
799 do_im_time, fm_mat_Q_gemm, fm_mat_Q, qs_env)
800
801 TYPE(two_dim_real_array), DIMENSION(:), &
802 INTENT(INOUT) :: bib_c_2d
803 TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_sub
804 INTEGER, INTENT(IN) :: color_sub, ngroup, integ_group_size, &
805 dimen_ri
806 INTEGER, DIMENSION(:), INTENT(IN) :: dimen_ia
807 INTEGER, INTENT(IN) :: color_rpa_group, ext_row_block_size, &
808 ext_col_block_size, unit_nr
809 INTEGER, DIMENSION(:), INTENT(IN) :: my_ia_size, my_ia_start, my_ia_end
810 INTEGER, INTENT(IN) :: my_group_l_size, my_group_l_start, &
811 my_group_l_end
812 TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env_rpa
813 TYPE(cp_fm_type), DIMENSION(:), INTENT(INOUT) :: fm_mat_s
814 INTEGER, INTENT(INOUT) :: nrow_block_mat, ncol_block_mat
815 TYPE(cp_blacs_env_type), OPTIONAL, POINTER :: blacs_env_ext, blacs_env_ext_s
816 INTEGER, INTENT(IN), OPTIONAL :: dimen_ia_for_block_size
817 LOGICAL, INTENT(IN), OPTIONAL :: do_im_time
818 TYPE(cp_fm_type), DIMENSION(:), OPTIONAL :: fm_mat_q_gemm, fm_mat_q
819 TYPE(qs_environment_type), INTENT(IN), OPTIONAL, &
820 POINTER :: qs_env
821
822 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_integ_mat'
823
824 INTEGER :: col_row_proc_ratio, grid_2d(2), handle, &
825 iproc, iproc_col, iproc_row, ispin, &
826 mepos_in_rpa_group
827 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: group_grid_2_mepos
828 LOGICAL :: my_blacs_ext, my_blacs_s_ext, &
829 my_do_im_time
830 TYPE(cp_blacs_env_type), POINTER :: blacs_env, blacs_env_q
831 TYPE(cp_fm_struct_type), POINTER :: fm_struct
832 TYPE(group_dist_d1_type) :: gd_ia, gd_l
833
834 CALL timeset(routinen, handle)
835
836 cpassert(PRESENT(blacs_env_ext) .OR. PRESENT(dimen_ia_for_block_size))
837
838 my_blacs_ext = .false.
839 IF (PRESENT(blacs_env_ext)) my_blacs_ext = .true.
840
841 my_blacs_s_ext = .false.
842 IF (PRESENT(blacs_env_ext_s)) my_blacs_s_ext = .true.
843
844 my_do_im_time = .false.
845 IF (PRESENT(do_im_time)) my_do_im_time = do_im_time
846
847 NULLIFY (blacs_env)
848 ! create the RPA blacs env
849 IF (my_blacs_s_ext) THEN
850 blacs_env => blacs_env_ext_s
851 ELSE
852 IF (para_env_rpa%num_pe > 1) THEN
853 col_row_proc_ratio = max(1, dimen_ia_for_block_size/dimen_ri)
854
855 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
856 DO iproc = 1, para_env_rpa%num_pe
857 iproc_col = iproc_col - 1
858 IF (mod(para_env_rpa%num_pe, iproc_col) == 0) EXIT
859 END DO
860
861 iproc_row = para_env_rpa%num_pe/iproc_col
862 grid_2d(1) = iproc_row
863 grid_2d(2) = iproc_col
864 ELSE
865 grid_2d = 1
866 END IF
867 CALL cp_blacs_env_create(blacs_env=blacs_env, para_env=para_env_rpa, grid_2d=grid_2d)
868
869 IF (unit_nr > 0 .AND. .NOT. my_do_im_time) THEN
870 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
871 "MATRIX_INFO| Number row processes:", grid_2d(1)
872 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
873 "MATRIX_INFO| Number column processes:", grid_2d(2)
874 END IF
875
876 ! define the block_size for the row
877 IF (ext_row_block_size > 0) THEN
878 nrow_block_mat = ext_row_block_size
879 ELSE
880 nrow_block_mat = max(1, dimen_ri/grid_2d(1)/2)
881 END IF
882
883 ! define the block_size for the column
884 IF (ext_col_block_size > 0) THEN
885 ncol_block_mat = ext_col_block_size
886 ELSE
887 ncol_block_mat = max(1, dimen_ia_for_block_size/grid_2d(2)/2)
888 END IF
889
890 IF (unit_nr > 0 .AND. .NOT. my_do_im_time) THEN
891 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
892 "MATRIX_INFO| Row block size:", nrow_block_mat
893 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
894 "MATRIX_INFO| Column block size:", ncol_block_mat
895 END IF
896 END IF
897
898 IF (.NOT. my_do_im_time) THEN
899 DO ispin = 1, SIZE(bib_c_2d)
900 NULLIFY (fm_struct)
901 IF (my_blacs_ext) THEN
902 CALL cp_fm_struct_create(fm_struct, context=blacs_env, nrow_global=dimen_ri, &
903 ncol_global=dimen_ia(ispin), para_env=para_env_rpa)
904 ELSE
905 CALL cp_fm_struct_create(fm_struct, context=blacs_env, nrow_global=dimen_ri, &
906 ncol_global=dimen_ia(ispin), para_env=para_env_rpa, &
907 nrow_block=nrow_block_mat, ncol_block=ncol_block_mat, force_block=.true.)
908
909 END IF ! external blacs_env
910
911 CALL create_group_dist(gd_ia, my_ia_start(ispin), my_ia_end(ispin), my_ia_size(ispin), para_env_rpa)
912 CALL create_group_dist(gd_l, my_group_l_start, my_group_l_end, my_group_l_size, para_env_rpa)
913
914 ! create the info array
915
916 mepos_in_rpa_group = mod(color_sub, integ_group_size)
917 ALLOCATE (group_grid_2_mepos(0:integ_group_size - 1, 0:para_env_sub%num_pe - 1))
918 group_grid_2_mepos = 0
919 group_grid_2_mepos(mepos_in_rpa_group, para_env_sub%mepos) = para_env_rpa%mepos
920 CALL para_env_rpa%sum(group_grid_2_mepos)
921
922 CALL array2fm(bib_c_2d(ispin)%array, fm_struct, my_group_l_start, my_group_l_end, &
923 my_ia_start(ispin), my_ia_end(ispin), gd_l, gd_ia, &
924 group_grid_2_mepos, ngroup, para_env_sub%num_pe, fm_mat_s(ispin), &
925 integ_group_size, color_rpa_group)
926
927 DEALLOCATE (group_grid_2_mepos)
928 CALL cp_fm_struct_release(fm_struct)
929
930 ! deallocate the info array
931 CALL release_group_dist(gd_l)
932 CALL release_group_dist(gd_ia)
933
934 ! sum the local data across processes belonging to different RPA group.
935 IF (para_env_rpa%num_pe /= para_env%num_pe) THEN
936 block
937 TYPE(mp_comm_type) :: comm_exchange
938 comm_exchange = fm_mat_s(ispin)%matrix_struct%context%interconnect(para_env)
939 CALL comm_exchange%sum(fm_mat_s(ispin)%local_data)
940 CALL comm_exchange%free()
941 END block
942 END IF
943 END DO
944 END IF
945
946 IF (PRESENT(fm_mat_q_gemm) .AND. .NOT. my_do_im_time) THEN
947 ! create the Q matrix dimen_RIxdimen_RI where the result of the mat-mat-mult will be stored
948 NULLIFY (fm_struct)
949 CALL cp_fm_struct_create(fm_struct, context=blacs_env, nrow_global=dimen_ri, &
950 ncol_global=dimen_ri, para_env=para_env_rpa, &
951 nrow_block=nrow_block_mat, ncol_block=ncol_block_mat, force_block=.true.)
952 DO ispin = 1, SIZE(fm_mat_q_gemm)
953 CALL cp_fm_create(fm_mat_q_gemm(ispin), fm_struct, name="fm_mat_Q_gemm")
954 END DO
955 CALL cp_fm_struct_release(fm_struct)
956 END IF
957
958 IF (PRESENT(fm_mat_q)) THEN
959 NULLIFY (blacs_env_q)
960 IF (my_blacs_ext) THEN
961 blacs_env_q => blacs_env_ext
962 ELSE IF (para_env_rpa%num_pe == para_env%num_pe .AND. PRESENT(qs_env)) THEN
963 CALL get_qs_env(qs_env, blacs_env=blacs_env_q)
964 ELSE
965 CALL cp_blacs_env_create(blacs_env=blacs_env_q, para_env=para_env_rpa)
966 END IF
967 NULLIFY (fm_struct)
968 CALL cp_fm_struct_create(fm_struct, context=blacs_env_q, nrow_global=dimen_ri, &
969 ncol_global=dimen_ri, para_env=para_env_rpa)
970 DO ispin = 1, SIZE(fm_mat_q)
971 CALL cp_fm_create(fm_mat_q(ispin), fm_struct, name="fm_mat_Q", set_zero=.true.)
972 END DO
973
974 CALL cp_fm_struct_release(fm_struct)
975
976 IF (.NOT. (my_blacs_ext .OR. (para_env_rpa%num_pe == para_env%num_pe .AND. PRESENT(qs_env)))) &
977 CALL cp_blacs_env_release(blacs_env_q)
978 END IF
979
980 ! release blacs_env
981 IF (.NOT. my_blacs_s_ext) THEN
982 CALL cp_blacs_env_release(blacs_env)
983 ELSE
984 NULLIFY (blacs_env)
985 END IF
986
987 CALL timestop(handle)
988
989 END SUBROUTINE create_integ_mat
990
991! **************************************************************************************************
992!> \brief ...
993!> \param qs_env ...
994!> \param Erpa ...
995!> \param mp2_env ...
996!> \param para_env ...
997!> \param para_env_RPA ...
998!> \param para_env_sub ...
999!> \param unit_nr ...
1000!> \param homo ...
1001!> \param virtual ...
1002!> \param dimen_RI ...
1003!> \param dimen_RI_red ...
1004!> \param dimen_ia ...
1005!> \param dimen_nm_gw ...
1006!> \param Eigenval ...
1007!> \param num_integ_points ...
1008!> \param num_integ_group ...
1009!> \param color_rpa_group ...
1010!> \param fm_matrix_PQ ...
1011!> \param fm_mat_S ...
1012!> \param fm_mat_Q_gemm ...
1013!> \param fm_mat_Q ...
1014!> \param fm_mat_S_gw ...
1015!> \param fm_mat_R_gw ...
1016!> \param fm_mat_S_ij_bse ...
1017!> \param fm_mat_S_ab_bse ...
1018!> \param my_do_gw ...
1019!> \param do_bse ...
1020!> \param gw_corr_lev_occ ...
1021!> \param gw_corr_lev_virt ...
1022!> \param bse_lev_virt ...
1023!> \param do_minimax_quad ...
1024!> \param do_im_time ...
1025!> \param mo_coeff ...
1026!> \param fm_matrix_L_kpoints ...
1027!> \param fm_matrix_Minv_L_kpoints ...
1028!> \param fm_matrix_Minv ...
1029!> \param fm_matrix_Minv_Vtrunc_Minv ...
1030!> \param mat_munu ...
1031!> \param mat_P_global ...
1032!> \param t_3c_M ...
1033!> \param t_3c_O ...
1034!> \param t_3c_O_compressed ...
1035!> \param t_3c_O_ind ...
1036!> \param starts_array_mc ...
1037!> \param ends_array_mc ...
1038!> \param starts_array_mc_block ...
1039!> \param ends_array_mc_block ...
1040!> \param matrix_s ...
1041!> \param do_kpoints_from_Gamma ...
1042!> \param kpoints ...
1043!> \param gd_array ...
1044!> \param color_sub ...
1045!> \param do_ri_sos_laplace_mp2 ...
1046!> \param calc_forces ...
1047! **************************************************************************************************
1048 SUBROUTINE rpa_num_int(qs_env, Erpa, mp2_env, para_env, para_env_RPA, para_env_sub, unit_nr, &
1049 homo, virtual, dimen_RI, dimen_RI_red, dimen_ia, dimen_nm_gw, &
1050 Eigenval, num_integ_points, num_integ_group, color_rpa_group, &
1051 fm_matrix_PQ, fm_mat_S, fm_mat_Q_gemm, fm_mat_Q, fm_mat_S_gw, fm_mat_R_gw, &
1052 fm_mat_S_ij_bse, fm_mat_S_ab_bse, &
1053 my_do_gw, do_bse, gw_corr_lev_occ, gw_corr_lev_virt, &
1054 bse_lev_virt, &
1055 do_minimax_quad, do_im_time, mo_coeff, &
1056 fm_matrix_L_kpoints, fm_matrix_Minv_L_kpoints, &
1057 fm_matrix_Minv, fm_matrix_Minv_Vtrunc_Minv, mat_munu, mat_P_global, &
1058 t_3c_M, t_3c_O, t_3c_O_compressed, t_3c_O_ind, &
1059 starts_array_mc, ends_array_mc, &
1060 starts_array_mc_block, ends_array_mc_block, &
1061 matrix_s, do_kpoints_from_Gamma, kpoints, gd_array, color_sub, &
1062 do_ri_sos_laplace_mp2, calc_forces)
1063
1064 TYPE(qs_environment_type), POINTER :: qs_env
1065 REAL(kind=dp), INTENT(OUT) :: erpa
1066 TYPE(mp2_type) :: mp2_env
1067 TYPE(mp_para_env_type), POINTER :: para_env, para_env_rpa, para_env_sub
1068 INTEGER, INTENT(IN) :: unit_nr
1069 INTEGER, DIMENSION(:), INTENT(IN) :: homo, virtual
1070 INTEGER, INTENT(IN) :: dimen_ri, dimen_ri_red
1071 INTEGER, DIMENSION(:), INTENT(IN) :: dimen_ia
1072 INTEGER, INTENT(IN) :: dimen_nm_gw
1073 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :), &
1074 INTENT(INOUT) :: eigenval
1075 INTEGER, INTENT(IN) :: num_integ_points, num_integ_group, &
1076 color_rpa_group
1077 TYPE(cp_fm_type), INTENT(IN) :: fm_matrix_pq
1078 TYPE(cp_fm_type), DIMENSION(:), INTENT(INOUT) :: fm_mat_s
1079 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_mat_q_gemm, fm_mat_q, fm_mat_s_gw
1080 TYPE(cp_fm_type), INTENT(IN) :: fm_mat_r_gw, fm_mat_s_ij_bse, &
1081 fm_mat_s_ab_bse
1082 LOGICAL, INTENT(IN) :: my_do_gw, do_bse
1083 INTEGER, DIMENSION(:), INTENT(IN) :: gw_corr_lev_occ, gw_corr_lev_virt
1084 INTEGER, INTENT(IN) :: bse_lev_virt
1085 LOGICAL, INTENT(IN) :: do_minimax_quad, do_im_time
1086 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
1087 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_l_kpoints, &
1088 fm_matrix_minv_l_kpoints, &
1089 fm_matrix_minv, &
1090 fm_matrix_minv_vtrunc_minv
1091 TYPE(dbcsr_p_type), INTENT(IN) :: mat_munu
1092 TYPE(dbcsr_p_type), INTENT(INOUT) :: mat_p_global
1093 TYPE(dbt_type), INTENT(INOUT) :: t_3c_m
1094 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :), &
1095 INTENT(INOUT) :: t_3c_o
1096 TYPE(hfx_compression_type), ALLOCATABLE, &
1097 DIMENSION(:, :, :), INTENT(INOUT) :: t_3c_o_compressed
1098 TYPE(block_ind_type), ALLOCATABLE, &
1099 DIMENSION(:, :, :), INTENT(INOUT) :: t_3c_o_ind
1100 INTEGER, ALLOCATABLE, DIMENSION(:), INTENT(IN) :: starts_array_mc, ends_array_mc, &
1101 starts_array_mc_block, &
1102 ends_array_mc_block
1103 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
1104 LOGICAL :: do_kpoints_from_gamma
1105 TYPE(kpoint_type), POINTER :: kpoints
1106 TYPE(group_dist_d1_type), INTENT(IN) :: gd_array
1107 INTEGER, INTENT(IN) :: color_sub
1108 LOGICAL, INTENT(IN) :: do_ri_sos_laplace_mp2, calc_forces
1109
1110 CHARACTER(LEN=*), PARAMETER :: routinen = 'rpa_num_int'
1111
1112 COMPLEX(KIND=dp), ALLOCATABLE, &
1113 DIMENSION(:, :, :, :) :: vec_sigma_c_gw
1114 INTEGER :: count_ev_sc_gw, cut_memory, group_size_p, gw_corr_lev_tot, handle, handle3, i, &
1115 ikp_local, ispin, iter_evgw, iter_sc_gw0, j, jquad, min_bsize, mm_style, nkp, &
1116 nkp_self_energy, nmo, nspins, num_3c_repl, num_cells_dm, num_fit_points, pspin, qspin, &
1117 size_p
1118 INTEGER(int_8) :: dbcsr_nflop
1119 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: index_to_cell_3c
1120 INTEGER, ALLOCATABLE, DIMENSION(:, :, :) :: cell_to_index_3c
1121 INTEGER, DIMENSION(:), POINTER :: col_blk_size, prim_blk_sizes, &
1122 ri_blk_sizes
1123 LOGICAL :: do_apply_ic_corr_to_gw, do_gw_im_time, do_ic_model, do_kpoints_cubic_rpa, &
1124 do_periodic, do_print, do_ri_sigma_x, exit_ev_gw, first_cycle, &
1125 first_cycle_periodic_correction, my_open_shell, print_ic_values
1126 LOGICAL, ALLOCATABLE, DIMENSION(:, :, :, :, :) :: has_mat_p_blocks
1127 REAL(kind=dp) :: a_scaling, alpha, dbcsr_time, e_exchange, e_exchange_corr, eps_filter, &
1128 eps_filter_im_time, ext_scaling, fermi_level_offset, fermi_level_offset_input, &
1129 my_flop_rate, omega, omega_max_fit, omega_old, tau, tau_old
1130 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: delta_corr, e_fermi, tau_tj, tau_wj, tj, &
1131 trace_qomega, vec_omega_fit_gw, wj, &
1132 wkp_w
1133 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: vec_w_gw, weights_cos_tf_t_to_w, &
1134 weights_cos_tf_w_to_t, &
1135 weights_sin_tf_t_to_w
1136 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: eigenval_last, eigenval_scf, &
1137 vec_sigma_x_gw
1138 TYPE(cp_cfm_type) :: cfm_mat_q
1139 TYPE(cp_fm_type) :: fm_mat_q_static_bse_gemm, fm_mat_ri_global_work, fm_mat_s_ia_bse, &
1140 fm_mat_work, fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, fm_scaled_dm_occ_tau, &
1141 fm_scaled_dm_virt_tau
1142 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: fm_mat_s_gw_work, fm_mat_w, &
1143 fm_mo_coeff_occ, fm_mo_coeff_virt
1144 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_mat_l_kpoints, fm_mat_minv_l_kpoints
1145 TYPE(dbcsr_p_type) :: mat_dm, mat_l, mat_m_p_munu_occ, &
1146 mat_m_p_munu_virt, mat_minvvminv
1147 TYPE(dbcsr_p_type), ALLOCATABLE, &
1148 DIMENSION(:, :, :) :: mat_p_omega
1149 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_berry_im_mo_mo, &
1150 matrix_berry_re_mo_mo
1151 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_p_omega_kp
1152 TYPE(dbcsr_type), POINTER :: mat_w, mat_work
1153 TYPE(dbt_type) :: t_3c_overl_int_ao_mo
1154 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:) :: t_3c_overl_int_gw_ao, &
1155 t_3c_overl_int_gw_ri, &
1156 t_3c_overl_nnp_ic, &
1157 t_3c_overl_nnp_ic_reflected
1158 TYPE(dgemm_counter_type) :: dgemm_counter
1159 TYPE(hfx_compression_type), ALLOCATABLE, &
1160 DIMENSION(:) :: t_3c_o_mo_compressed
1161 TYPE(im_time_force_type) :: force_data
1162 TYPE(rpa_exchange_work_type) :: exchange_work
1163 TYPE(rpa_grad_type) :: rpa_grad
1164 TYPE(rpa_sigma_type) :: rpa_sigma
1165 TYPE(two_dim_int_array), ALLOCATABLE, DIMENSION(:) :: t_3c_o_mo_ind
1166
1167 CALL timeset(routinen, handle)
1168
1169 nspins = SIZE(homo)
1170 nmo = homo(1) + virtual(1)
1171
1172 my_open_shell = (nspins == 2)
1173
1174 do_gw_im_time = my_do_gw .AND. do_im_time
1175 do_ri_sigma_x = mp2_env%ri_g0w0%do_ri_Sigma_x
1176 do_ic_model = mp2_env%ri_g0w0%do_ic_model
1177 print_ic_values = mp2_env%ri_g0w0%print_ic_values
1178 do_periodic = mp2_env%ri_g0w0%do_periodic
1179 do_kpoints_cubic_rpa = mp2_env%ri_rpa_im_time%do_im_time_kpoints
1180
1181 ! For SOS-MP2 only gemm is implemented
1182 mm_style = wfc_mm_style_gemm
1183 IF (.NOT. do_ri_sos_laplace_mp2) mm_style = mp2_env%ri_rpa%mm_style
1184
1185 IF (my_do_gw) THEN
1186 ext_scaling = 0.2_dp
1187 omega_max_fit = mp2_env%ri_g0w0%omega_max_fit
1188 fermi_level_offset_input = mp2_env%ri_g0w0%fermi_level_offset
1189 iter_evgw = mp2_env%ri_g0w0%iter_evGW
1190 iter_sc_gw0 = mp2_env%ri_g0w0%iter_sc_GW0
1191 IF ((.NOT. do_im_time)) THEN
1192 IF (iter_sc_gw0 /= 1 .AND. iter_evgw /= 1) cpabort("Mixed scGW0/evGW not implemented.")
1193 ! in case of scGW0 with the N^4 algorithm, we use the evGW code but use the DFT eigenvalues for W
1194 IF (iter_sc_gw0 /= 1) iter_evgw = iter_sc_gw0
1195 END IF
1196 ELSE
1197 ext_scaling = 0.0_dp
1198 iter_evgw = 1
1199 iter_sc_gw0 = 1
1200 END IF
1201
1202 IF (do_kpoints_cubic_rpa .AND. do_ri_sos_laplace_mp2) THEN
1203 cpabort("RI-SOS-Laplace-MP2 with k-point-sampling is not implemented.")
1204 END IF
1205
1206 do_apply_ic_corr_to_gw = .false.
1207 IF (mp2_env%ri_g0w0%ic_corr_list(1)%array(1) > 0.0_dp) do_apply_ic_corr_to_gw = .true.
1208
1209 IF (do_im_time) THEN
1210 cpassert(do_minimax_quad .OR. do_ri_sos_laplace_mp2)
1211
1212 group_size_p = mp2_env%ri_rpa_im_time%group_size_P
1213 cut_memory = mp2_env%ri_rpa_im_time%cut_memory
1214 eps_filter = mp2_env%ri_rpa_im_time%eps_filter
1215 eps_filter_im_time = mp2_env%ri_rpa_im_time%eps_filter* &
1216 mp2_env%ri_rpa_im_time%eps_filter_factor
1217
1218 min_bsize = mp2_env%ri_rpa_im_time%min_bsize
1219
1220 CALL alloc_im_time(qs_env, para_env, dimen_ri, dimen_ri_red, &
1221 num_integ_points, nspins, fm_mat_q(1), fm_mo_coeff_occ, fm_mo_coeff_virt, &
1222 fm_matrix_minv_l_kpoints, fm_matrix_l_kpoints, mat_p_global, &
1223 t_3c_o, matrix_s, kpoints, eps_filter_im_time, &
1224 cut_memory, nkp, num_cells_dm, num_3c_repl, &
1225 size_p, ikp_local, &
1226 index_to_cell_3c, &
1227 cell_to_index_3c, &
1228 col_blk_size, &
1229 do_ic_model, do_kpoints_cubic_rpa, &
1230 do_kpoints_from_gamma, do_ri_sigma_x, my_open_shell, &
1231 has_mat_p_blocks, wkp_w, &
1232 cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1233 fm_mat_ri_global_work, fm_mat_work, fm_mo_coeff_occ_scaled, &
1234 fm_mo_coeff_virt_scaled, mat_dm, mat_l, mat_m_p_munu_occ, mat_m_p_munu_virt, &
1235 mat_minvvminv, mat_p_omega, mat_p_omega_kp, mat_work, mo_coeff, &
1236 fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, homo, nmo)
1237
1238 IF (calc_forces) CALL init_im_time_forces(force_data, fm_matrix_pq, t_3c_m, unit_nr, mp2_env, qs_env)
1239
1240 IF (my_do_gw) THEN
1241
1242 CALL dbcsr_get_info(mat_p_global%matrix, &
1243 row_blk_size=ri_blk_sizes)
1244
1245 CALL dbcsr_get_info(matrix_s(1)%matrix, &
1246 row_blk_size=prim_blk_sizes)
1247
1248 gw_corr_lev_tot = gw_corr_lev_occ(1) + gw_corr_lev_virt(1)
1249
1250 IF (.NOT. do_kpoints_cubic_rpa) THEN
1251 CALL allocate_matrices_gw_im_time(gw_corr_lev_occ, gw_corr_lev_virt, homo, nmo, &
1252 num_integ_points, unit_nr, &
1253 ri_blk_sizes, do_ic_model, &
1254 para_env, fm_mat_w, fm_mat_q(1), &
1255 mo_coeff, &
1256 t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, &
1257 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1258 starts_array_mc, ends_array_mc, &
1259 t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, &
1260 matrix_s, mat_w, t_3c_o, &
1261 t_3c_o_compressed, t_3c_o_ind, &
1262 qs_env)
1263
1264 END IF
1265 END IF
1266
1267 END IF
1268 IF (do_ic_model) THEN
1269 ! image charge model only implemented for cubic scaling GW
1270 cpassert(do_gw_im_time)
1271 cpassert(.NOT. do_periodic)
1272 IF (cut_memory /= 1) cpabort("For IC, use MEMORY_CUT 1 in the LOW_SCALING section.")
1273 END IF
1274
1275 ALLOCATE (e_fermi(nspins))
1276 IF (do_minimax_quad .OR. do_ri_sos_laplace_mp2) THEN
1277 do_print = .NOT. do_ic_model
1278 CALL get_minimax_grid(para_env, unit_nr, homo, eigenval, num_integ_points, do_im_time, &
1279 do_ri_sos_laplace_mp2, do_print, &
1280 tau_tj, tau_wj, qs_env, do_gw_im_time, &
1281 do_kpoints_cubic_rpa, e_fermi(1), tj, wj, &
1282 weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, &
1283 qs_env%mp2_env%ri_g0w0%regularization_minimax)
1284
1285 !For sos_laplace_mp2 and low-scaling RPA, potentially need to store/retrieve the initial weights
1286 IF (qs_env%mp2_env%ri_rpa_im_time%keep_quad) THEN
1287 CALL keep_initial_quad(tj, wj, tau_tj, tau_wj, weights_cos_tf_t_to_w, &
1288 weights_cos_tf_w_to_t, do_ri_sos_laplace_mp2, do_im_time, &
1289 num_integ_points, unit_nr, qs_env)
1290 END IF
1291 ELSE
1292 IF (calc_forces) cpabort("Forces with Clenshaw-Curtis grid not implemented.")
1293 CALL get_clenshaw_grid(para_env, para_env_rpa, unit_nr, homo, virtual, eigenval, num_integ_points, &
1294 num_integ_group, color_rpa_group, fm_mat_s, my_do_gw, &
1295 ext_scaling, a_scaling, tj, wj)
1296 END IF
1297
1298 ! This array is needed for RPA
1299 IF (.NOT. do_ri_sos_laplace_mp2) THEN
1300 ALLOCATE (trace_qomega(dimen_ri_red))
1301 END IF
1302
1303 IF (do_ri_sos_laplace_mp2 .AND. .NOT. do_im_time) THEN
1304 alpha = 1.0_dp
1305 ELSE IF (my_open_shell .OR. do_ri_sos_laplace_mp2) THEN
1306 alpha = 2.0_dp
1307 ELSE
1308 alpha = 4.0_dp
1309 END IF
1310 IF (my_do_gw) THEN
1311 CALL allocate_matrices_gw(vec_sigma_c_gw, color_rpa_group, dimen_nm_gw, &
1312 gw_corr_lev_occ, gw_corr_lev_virt, homo, &
1313 nmo, num_integ_group, num_integ_points, unit_nr, &
1314 gw_corr_lev_tot, num_fit_points, omega_max_fit, &
1315 do_minimax_quad, do_periodic, do_ri_sigma_x,.NOT. do_im_time, &
1316 first_cycle_periodic_correction, &
1317 a_scaling, eigenval, tj, vec_omega_fit_gw, vec_sigma_x_gw, &
1318 delta_corr, eigenval_last, eigenval_scf, vec_w_gw, &
1319 fm_mat_s_gw, fm_mat_s_gw_work, &
1320 para_env, mp2_env, kpoints, nkp, nkp_self_energy, &
1321 do_kpoints_cubic_rpa, do_kpoints_from_gamma)
1322
1323 IF (do_bse) THEN
1324
1325 CALL cp_fm_create(fm_mat_q_static_bse_gemm, fm_mat_q_gemm(1)%matrix_struct)
1326 CALL cp_fm_to_fm(fm_mat_q_gemm(1), fm_mat_q_static_bse_gemm)
1327 CALL cp_fm_set_all(fm_mat_q_static_bse_gemm, 0.0_dp)
1328
1329 END IF
1330
1331 END IF
1332
1333 IF (calc_forces .AND. .NOT. do_im_time) CALL rpa_grad_create(rpa_grad, fm_mat_q(1), &
1334 fm_mat_s, homo, virtual, mp2_env, eigenval(:, 1, :), &
1335 unit_nr, do_ri_sos_laplace_mp2)
1336 IF (.NOT. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2) &
1337 CALL exchange_work%create(qs_env, para_env_sub, mat_munu, dimen_ri_red, &
1338 fm_mat_s, fm_mat_q(1), fm_mat_q_gemm(1), homo, virtual)
1339 erpa = 0.0_dp
1340 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) e_exchange = 0.0_dp
1341 first_cycle = .true.
1342 omega_old = 0.0_dp
1343 CALL dgemm_counter_init(dgemm_counter, unit_nr, mp2_env%ri_rpa%print_dgemm_info)
1344
1345 DO count_ev_sc_gw = 1, iter_evgw
1346 dbcsr_time = 0.0_dp
1347 dbcsr_nflop = 0
1348
1349 IF (do_ic_model) cycle
1350
1351 ! reset some values, important when doing eigenvalue self-consistent GW
1352 IF (my_do_gw) THEN
1353 erpa = 0.0_dp
1354 vec_sigma_c_gw = z_zero
1355 first_cycle = .true.
1356 END IF
1357
1358 ! calculate Q_PQ(it)
1359 IF (do_im_time) THEN ! not using Imaginary time
1360
1361 IF (.NOT. do_kpoints_cubic_rpa) THEN
1362 DO ispin = 1, nspins
1363 e_fermi(ispin) = (eigenval(homo(ispin), 1, ispin) + eigenval(homo(ispin) + 1, 1, ispin))*0.5_dp
1364 END DO
1365 END IF
1366
1367 tau = 0.0_dp
1368 tau_old = 0.0_dp
1369
1370 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(/T3,A,T66,i15)") &
1371 "MEMORY_INFO| Memory cut:", cut_memory
1372 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T66,ES15.2)") &
1373 "SPARSITY_INFO| Eps filter for M virt/occ tensors:", eps_filter
1374 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T66,ES15.2)") &
1375 "SPARSITY_INFO| Eps filter for P matrix:", eps_filter_im_time
1376 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T66,i15)") &
1377 "SPARSITY_INFO| Minimum tensor block size:", min_bsize
1378
1379 ! for evGW, we have to ensure that mat_P_omega is zero
1380 CALL zero_mat_p_omega(mat_p_omega(:, :, 1))
1381
1382 ! compute the matrix Q(it) and Fourier transform it directly to mat_P_omega(iw)
1383 CALL compute_mat_p_omega(mat_p_omega(:, :, 1), fm_scaled_dm_occ_tau, &
1384 fm_scaled_dm_virt_tau, fm_mo_coeff_occ(1), fm_mo_coeff_virt(1), &
1385 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1386 mat_p_global, matrix_s, 1, &
1387 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
1388 starts_array_mc, ends_array_mc, &
1389 starts_array_mc_block, ends_array_mc_block, &
1390 weights_cos_tf_t_to_w, tj, tau_tj, e_fermi(1), eps_filter, alpha, &
1391 eps_filter_im_time, eigenval(:, 1, 1), nmo, &
1392 num_integ_points, cut_memory, &
1393 unit_nr, mp2_env, para_env, &
1394 qs_env, do_kpoints_from_gamma, &
1395 index_to_cell_3c, cell_to_index_3c, &
1396 has_mat_p_blocks, do_ri_sos_laplace_mp2, &
1397 dbcsr_time, dbcsr_nflop)
1398
1399 ! the same for open shell, use fm_mo_coeff_occ_beta and fm_mo_coeff_virt_beta
1400 IF (my_open_shell) THEN
1401 CALL zero_mat_p_omega(mat_p_omega(:, :, 2))
1402 CALL compute_mat_p_omega(mat_p_omega(:, :, 2), fm_scaled_dm_occ_tau, &
1403 fm_scaled_dm_virt_tau, fm_mo_coeff_occ(2), &
1404 fm_mo_coeff_virt(2), &
1405 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1406 mat_p_global, matrix_s, 2, &
1407 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
1408 starts_array_mc, ends_array_mc, &
1409 starts_array_mc_block, ends_array_mc_block, &
1410 weights_cos_tf_t_to_w, tj, tau_tj, e_fermi(2), eps_filter, alpha, &
1411 eps_filter_im_time, eigenval(:, 1, 2), nmo, &
1412 num_integ_points, cut_memory, &
1413 unit_nr, mp2_env, para_env, &
1414 qs_env, do_kpoints_from_gamma, &
1415 index_to_cell_3c, cell_to_index_3c, &
1416 has_mat_p_blocks, do_ri_sos_laplace_mp2, &
1417 dbcsr_time, dbcsr_nflop)
1418
1419 !For RPA, we sum up the P matrices. If no force needed, can clean-up the beta spin one
1420 IF (.NOT. do_ri_sos_laplace_mp2) THEN
1421 DO j = 1, SIZE(mat_p_omega, 2)
1422 DO i = 1, SIZE(mat_p_omega, 1)
1423 CALL dbcsr_add(mat_p_omega(i, j, 1)%matrix, mat_p_omega(i, j, 2)%matrix, 1.0_dp, 1.0_dp)
1424 IF (.NOT. calc_forces) CALL dbcsr_clear(mat_p_omega(i, j, 2)%matrix)
1425 END DO
1426 END DO
1427 END IF
1428 END IF ! my_open_shell
1429
1430 END IF ! do im time
1431
1432 IF (mp2_env%ri_rpa%sigma_param /= sigma_none) THEN
1433 CALL rpa_sigma_create(rpa_sigma, mp2_env%ri_rpa%sigma_param, fm_mat_q(1), unit_nr, para_env)
1434 END IF
1435
1436 DO jquad = 1, num_integ_points
1437 IF (modulo(jquad, num_integ_group) /= color_rpa_group) cycle
1438
1439 CALL timeset(routinen//"_RPA_matrix_operations", handle3)
1440
1441 IF (do_ri_sos_laplace_mp2) THEN
1442 omega = tau_tj(jquad)
1443 ELSE
1444 IF (do_minimax_quad) THEN
1445 omega = tj(jquad)
1446 ELSE
1447 omega = a_scaling/tan(tj(jquad))
1448 END IF
1449 END IF ! do_ri_sos_laplace_mp2
1450
1451 IF (do_im_time) THEN
1452 ! in case we do imag time, we already calculated fm_mat_Q by a Fourier transform from im. time
1453
1454 IF (.NOT. (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma)) THEN
1455
1456 DO ispin = 1, SIZE(mat_p_omega, 3)
1457 CALL contract_p_omega_with_mat_l(mat_p_omega(jquad, 1, ispin)%matrix, mat_l%matrix, mat_work, &
1458 eps_filter_im_time, fm_mat_work, dimen_ri, dimen_ri_red, &
1459 fm_mat_minv_l_kpoints(1, 1), fm_mat_q(ispin))
1460 END DO
1461 END IF
1462
1463 ELSE
1464 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3, A, 1X, I3, 1X, A, 1X, I3)") &
1465 "INTEG_INFO| Started with Integration point", jquad, "of", num_integ_points
1466
1467 IF (first_cycle .AND. count_ev_sc_gw > 1) THEN
1468 IF (iter_sc_gw0 == 1) THEN
1469 DO ispin = 1, nspins
1470 CALL remove_scaling_factor_rpa(fm_mat_s(ispin), virtual(ispin), &
1471 eigenval_last(:, 1, ispin), homo(ispin), omega_old)
1472 END DO
1473 ELSE
1474 DO ispin = 1, nspins
1475 CALL remove_scaling_factor_rpa(fm_mat_s(ispin), virtual(ispin), &
1476 eigenval_scf(:, 1, ispin), homo(ispin), omega_old)
1477 END DO
1478 END IF
1479 END IF
1480
1481 IF (iter_sc_gw0 > 1) THEN
1482 DO ispin = 1, nspins
1483 CALL calc_mat_q(fm_mat_s(ispin), do_ri_sos_laplace_mp2, first_cycle, virtual(ispin), &
1484 eigenval_scf(:, 1, ispin), homo(ispin), omega, omega_old, jquad, mm_style, &
1485 dimen_ri_red, dimen_ia(ispin), alpha, fm_mat_q(ispin), &
1486 fm_mat_q_gemm(ispin), do_bse, fm_mat_q_static_bse_gemm, dgemm_counter, &
1487 num_integ_points, count_ev_sc_gw)
1488 END DO
1489
1490 ! For SOS-MP2 we need both matrices separately
1491 IF (.NOT. do_ri_sos_laplace_mp2) THEN
1492 DO ispin = 2, nspins
1493 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))
1494 END DO
1495 END IF
1496 ELSE
1497 DO ispin = 1, nspins
1498 CALL calc_mat_q(fm_mat_s(ispin), do_ri_sos_laplace_mp2, first_cycle, virtual(ispin), &
1499 eigenval(:, 1, ispin), homo(ispin), omega, omega_old, jquad, mm_style, &
1500 dimen_ri_red, dimen_ia(ispin), alpha, fm_mat_q(ispin), &
1501 fm_mat_q_gemm(ispin), do_bse, fm_mat_q_static_bse_gemm, dgemm_counter, &
1502 num_integ_points, count_ev_sc_gw)
1503 END DO
1504
1505 ! For SOS-MP2 we need both matrices separately
1506 IF (.NOT. do_ri_sos_laplace_mp2) THEN
1507 DO ispin = 2, nspins
1508 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))
1509 END DO
1510 END IF
1511
1512 END IF
1513
1514 END IF ! im time
1515
1516 ! Calculate RPA exchange energy correction
1517 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) THEN
1518 e_exchange_corr = 0.0_dp
1519 CALL exchange_work%compute(fm_mat_q(1), eigenval(:, 1, :), fm_mat_s, omega, e_exchange_corr, mp2_env)
1520
1521 ! Evaluate the final exchange energy correction
1522 e_exchange = e_exchange + e_exchange_corr*wj(jquad)
1523 END IF
1524
1525 ! for developing Sigma functional closed and open shell are taken cared for
1526 IF (mp2_env%ri_rpa%sigma_param /= sigma_none) THEN
1527 CALL rpa_sigma_matrix_spectral(rpa_sigma, fm_mat_q(1), wj(jquad), para_env_rpa)
1528 END IF
1529
1530 IF (do_ri_sos_laplace_mp2) THEN
1531
1532 CALL sos_mp2_postprocessing(fm_mat_q, erpa, tau_wj(jquad))
1533
1534 IF (calc_forces .AND. .NOT. do_im_time) CALL rpa_grad_matrix_operations(mp2_env, rpa_grad, do_ri_sos_laplace_mp2, &
1535 fm_mat_q, fm_mat_q_gemm, dgemm_counter, fm_mat_s, omega, homo, virtual, &
1536 eigenval(:, 1, :), tau_wj(jquad), unit_nr)
1537 ELSE
1538 IF (calc_forces .AND. .NOT. do_im_time) CALL rpa_grad_copy_q(fm_mat_q(1), rpa_grad)
1539
1540 CALL q_trace_and_add_unit_matrix(dimen_ri_red, trace_qomega, fm_mat_q(1))
1541
1542 IF (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma) THEN
1543 CALL invert_eps_compute_w_and_erpa_kp(dimen_ri, num_integ_points, jquad, nkp, count_ev_sc_gw, para_env, &
1544 erpa, tau_tj, tj, wj, weights_cos_tf_w_to_t, &
1545 wkp_w, do_gw_im_time, do_ri_sigma_x, do_kpoints_from_gamma, &
1546 cfm_mat_q, ikp_local, &
1547 mat_p_omega(:, :, 1), mat_p_omega_kp, qs_env, eps_filter_im_time, unit_nr, &
1548 kpoints, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1549 fm_mat_w, fm_mat_ri_global_work, mat_minvvminv, &
1550 fm_matrix_minv, fm_matrix_minv_vtrunc_minv)
1551 ELSE
1552 CALL compute_erpa_by_freq_int(dimen_ri_red, trace_qomega, fm_mat_q(1), para_env_rpa, erpa, wj(jquad))
1553 END IF
1554
1555 IF (calc_forces .AND. .NOT. do_im_time) CALL rpa_grad_matrix_operations(mp2_env, rpa_grad, do_ri_sos_laplace_mp2, &
1556 fm_mat_q, fm_mat_q_gemm, dgemm_counter, fm_mat_s, omega, homo, virtual, &
1557 eigenval(:, 1, :), wj(jquad), unit_nr)
1558 END IF ! do_ri_sos_laplace_mp2
1559
1560 ! save omega and reset the first_cycle flag
1561 first_cycle = .false.
1562 omega_old = omega
1563
1564 CALL timestop(handle3)
1565
1566 IF (my_do_gw) THEN
1567
1568 CALL get_fermi_level_offset(fermi_level_offset, fermi_level_offset_input, eigenval(:, 1, :), homo)
1569
1570 ! do_im_time = TRUE means low-scaling calculation
1571 IF (do_im_time) THEN
1572 ! only for molecules
1573 IF (.NOT. (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma)) THEN
1574 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, &
1575 tj, tau_tj, weights_cos_tf_w_to_t, jquad, omega)
1576 END IF
1577 ELSE
1578 CALL compute_gw_self_energy(vec_sigma_c_gw, dimen_nm_gw, dimen_ri_red, gw_corr_lev_occ, &
1579 gw_corr_lev_virt, homo, jquad, nmo, num_fit_points, &
1580 do_im_time, do_periodic, first_cycle_periodic_correction, &
1581 fermi_level_offset, &
1582 omega, eigenval(:, 1, :), delta_corr, vec_omega_fit_gw, vec_w_gw, wj, &
1583 fm_mat_q(1), fm_mat_r_gw, fm_mat_s_gw, &
1584 fm_mat_s_gw_work, mo_coeff(1), para_env, &
1585 para_env_rpa, matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, &
1586 kpoints, qs_env, mp2_env)
1587 END IF
1588 END IF
1589
1590 IF (unit_nr > 0) CALL m_flush(unit_nr)
1591 CALL para_env_rpa%sync() ! sync to see output
1592
1593 END DO ! jquad
1594
1595 IF (mp2_env%ri_rpa%sigma_param /= sigma_none) THEN
1596 CALL finalize_rpa_sigma(rpa_sigma, unit_nr, mp2_env%ri_rpa%e_sigma_corr, para_env, do_minimax_quad)
1597 IF (do_minimax_quad) mp2_env%ri_rpa%e_sigma_corr = mp2_env%ri_rpa%e_sigma_corr/2.0_dp
1598 CALL para_env%sum(mp2_env%ri_rpa%e_sigma_corr)
1599 END IF
1600
1601 CALL para_env%sum(erpa)
1602
1603 IF (.NOT. (do_ri_sos_laplace_mp2)) THEN
1604 erpa = erpa/(pi*2.0_dp)
1605 IF (do_minimax_quad) erpa = erpa/2.0_dp
1606 END IF
1607
1608 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) THEN
1609 CALL para_env%sum(e_exchange)
1610 e_exchange = e_exchange/(pi*2.0_dp)
1611 IF (do_minimax_quad) e_exchange = e_exchange/2.0_dp
1612 mp2_env%ri_rpa%ener_exchange = e_exchange
1613 END IF
1614
1615 IF (calc_forces .AND. do_ri_sos_laplace_mp2 .AND. do_im_time) THEN
1616 IF (my_open_shell) THEN
1617 pspin = 1
1618 qspin = 2
1619 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1620 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1621 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1622 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1623 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1624 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1625 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1626 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1627 pspin = 2
1628 qspin = 1
1629 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1630 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1631 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1632 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1633 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1634 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1635 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1636 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1637
1638 ELSE
1639 pspin = 1
1640 qspin = 1
1641 CALL calc_laplace_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1642 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1643 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1644 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1645 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1646 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1647 tau_tj, tau_wj, cut_memory, pspin, qspin, my_open_shell, &
1648 unit_nr, dbcsr_time, dbcsr_nflop, mp2_env, qs_env)
1649 END IF
1650 CALL calc_post_loop_forces(force_data, unit_nr, qs_env)
1651 END IF !laplace SOS-MP2
1652
1653 IF (calc_forces .AND. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2) THEN
1654 DO ispin = 1, nspins
1655 CALL calc_rpa_loop_forces(force_data, mat_p_omega(:, 1, :), t_3c_m, t_3c_o(1, 1), &
1656 t_3c_o_compressed(1, 1, :), t_3c_o_ind(1, 1, :), fm_scaled_dm_occ_tau, &
1657 fm_scaled_dm_virt_tau, fm_mo_coeff_occ, fm_mo_coeff_virt, &
1658 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, &
1659 starts_array_mc, ends_array_mc, starts_array_mc_block, &
1660 ends_array_mc_block, num_integ_points, nmo, eigenval(:, 1, :), &
1661 e_fermi(ispin), weights_cos_tf_t_to_w, weights_cos_tf_w_to_t, tj, &
1662 wj, tau_tj, cut_memory, ispin, my_open_shell, unit_nr, dbcsr_time, &
1663 dbcsr_nflop, mp2_env, qs_env)
1664 END DO
1665 CALL calc_post_loop_forces(force_data, unit_nr, qs_env)
1666 END IF
1667
1668 IF (do_im_time) THEN
1669
1670 my_flop_rate = real(dbcsr_nflop, dp)/(1.0e09_dp*dbcsr_time)
1671 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(/T3,A,T73,ES8.2)") &
1672 "PERFORMANCE| DBCSR total number of flops:", real(dbcsr_nflop*para_env%num_pe, dp)
1673 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T66,F15.2)") &
1674 "PERFORMANCE| DBCSR total execution time:", dbcsr_time
1675 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T66,F15.2)") &
1676 "PERFORMANCE| DBCSR flop rate (Gflops / MPI rank):", my_flop_rate
1677
1678 ELSE
1679
1680 CALL dgemm_counter_write(dgemm_counter, para_env)
1681
1682 END IF
1683
1684 ! GW: for low-scaling calculation: Compute self-energy Sigma(i*tau), Sigma(i*omega)
1685 ! for low-scaling and ordinary-scaling: analytic continuation from Sigma(iw) -> Sigma(w)
1686 ! and correction of quasiparticle energies e_n^GW
1687 IF (my_do_gw) THEN
1688
1689 CALL compute_qp_energies(vec_sigma_c_gw, count_ev_sc_gw, gw_corr_lev_occ, &
1690 gw_corr_lev_tot, gw_corr_lev_virt, homo, &
1691 nmo, num_fit_points, num_integ_points, &
1692 unit_nr, do_apply_ic_corr_to_gw, do_im_time, &
1693 do_periodic, do_ri_sigma_x, first_cycle_periodic_correction, &
1694 e_fermi, eps_filter, fermi_level_offset, &
1695 delta_corr, eigenval, &
1696 eigenval_last, eigenval_scf, iter_sc_gw0, exit_ev_gw, tau_tj, tj, &
1697 vec_omega_fit_gw, vec_sigma_x_gw, mp2_env%ri_g0w0%ic_corr_list, &
1698 weights_cos_tf_t_to_w, weights_sin_tf_t_to_w, &
1699 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, fm_mo_coeff_occ, &
1700 fm_mo_coeff_virt, fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, &
1701 mo_coeff(1), fm_mat_w, para_env, para_env_rpa, mat_dm, mat_minvvminv, &
1702 t_3c_o, t_3c_m, t_3c_overl_int_ao_mo, t_3c_o_compressed, t_3c_o_mo_compressed, &
1703 t_3c_o_ind, t_3c_o_mo_ind, &
1704 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1705 matrix_berry_im_mo_mo, matrix_berry_re_mo_mo, mat_w, matrix_s, &
1706 kpoints, mp2_env, qs_env, nkp_self_energy, do_kpoints_cubic_rpa, &
1707 starts_array_mc, ends_array_mc)
1708
1709 ! if HOMO-LUMO gap differs by less than mp2_env%ri_g0w0%eps_ev_sc_iter, exit ev sc GW loop
1710 IF (exit_ev_gw) EXIT
1711
1712 END IF ! my_do_gw if
1713
1714 END DO ! evGW loop
1715
1716 IF (do_ic_model) THEN
1717
1718 IF (my_open_shell) THEN
1719
1720 CALL calculate_ic_correction(eigenval(:, 1, 1), mat_minvvminv%matrix, &
1721 t_3c_overl_nnp_ic(1), t_3c_overl_nnp_ic_reflected(1), &
1722 gw_corr_lev_tot, &
1723 gw_corr_lev_occ(1), gw_corr_lev_virt(1), homo(1), unit_nr, &
1724 print_ic_values, para_env, do_alpha=.true.)
1725
1726 CALL calculate_ic_correction(eigenval(:, 1, 2), mat_minvvminv%matrix, &
1727 t_3c_overl_nnp_ic(2), t_3c_overl_nnp_ic_reflected(2), &
1728 gw_corr_lev_tot, &
1729 gw_corr_lev_occ(2), gw_corr_lev_virt(2), homo(2), unit_nr, &
1730 print_ic_values, para_env, do_beta=.true.)
1731
1732 ELSE
1733
1734 CALL calculate_ic_correction(eigenval(:, 1, 1), mat_minvvminv%matrix, &
1735 t_3c_overl_nnp_ic(1), t_3c_overl_nnp_ic_reflected(1), &
1736 gw_corr_lev_tot, &
1737 gw_corr_lev_occ(1), gw_corr_lev_virt(1), homo(1), unit_nr, &
1738 print_ic_values, para_env)
1739
1740 END IF
1741
1742 END IF
1743
1744 ! postprocessing after GW for Bethe-Salpeter
1745 IF (do_bse) THEN
1746 ! Check used GW flavor; in Case of evGW we use W0 for BSE
1747 ! Use environment variable, since local iter_evGW is overwritten if evGW0 is invoked
1748 IF (mp2_env%ri_g0w0%iter_evGW > 1) THEN
1749 IF (unit_nr > 0) THEN
1750 CALL cp_warn(__location__, &
1751 "BSE@evGW applies W0, i.e. screening with DFT energies to the BSE!")
1752 END IF
1753 END IF
1754 ! Create a copy of fm_mat_S for usage in BSE
1755 CALL cp_fm_create(fm_mat_s_ia_bse, fm_mat_s(1)%matrix_struct)
1756 CALL cp_fm_to_fm(fm_mat_s(1), fm_mat_s_ia_bse)
1757 ! Remove energy/frequency factor from 3c-Integral for BSE
1758 IF (iter_sc_gw0 == 1) THEN
1759 CALL remove_scaling_factor_rpa(fm_mat_s_ia_bse, virtual(1), &
1760 eigenval_last(:, 1, 1), homo(1), omega)
1761 ELSE
1762 CALL remove_scaling_factor_rpa(fm_mat_s_ia_bse, virtual(1), &
1763 eigenval_scf(:, 1, 1), homo(1), omega)
1764 END IF
1765 ! Main routine for all BSE postprocessing
1766 CALL start_bse_calculation(fm_mat_s_ia_bse, fm_mat_s_ij_bse, fm_mat_s_ab_bse, &
1767 fm_mat_q_static_bse_gemm, &
1768 eigenval, eigenval_scf, &
1769 homo, virtual, dimen_ri, dimen_ri_red, bse_lev_virt, &
1770 gd_array, color_sub, mp2_env, qs_env, mo_coeff, unit_nr)
1771 ! Release BSE-copy of fm_mat_S
1772 CALL cp_fm_release(fm_mat_s_ia_bse)
1773 END IF
1774
1775 IF (my_do_gw) THEN
1776 CALL deallocate_matrices_gw(fm_mat_s_gw_work, vec_w_gw, vec_sigma_c_gw, vec_omega_fit_gw, &
1777 mp2_env%ri_g0w0%vec_Sigma_x_minus_vxc_gw, &
1778 eigenval_last, eigenval_scf, do_periodic, matrix_berry_re_mo_mo, matrix_berry_im_mo_mo, &
1779 kpoints, vec_sigma_x_gw,.NOT. do_im_time)
1780 END IF
1781
1782 IF (do_im_time) THEN
1783
1784 CALL dealloc_im_time(fm_mo_coeff_occ, fm_mo_coeff_virt, &
1785 fm_scaled_dm_occ_tau, fm_scaled_dm_virt_tau, index_to_cell_3c, &
1786 cell_to_index_3c, do_ic_model, &
1787 do_kpoints_cubic_rpa, do_kpoints_from_gamma, do_ri_sigma_x, &
1788 has_mat_p_blocks, &
1789 wkp_w, cfm_mat_q, fm_mat_minv_l_kpoints, fm_mat_l_kpoints, &
1790 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, fm_mat_ri_global_work, fm_mat_work, &
1791 fm_mo_coeff_occ_scaled, fm_mo_coeff_virt_scaled, mat_dm, mat_l, &
1792 mat_minvvminv, mat_p_omega, mat_p_omega_kp, &
1793 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, mat_work, qs_env)
1794
1795 IF (my_do_gw) THEN
1796 CALL deallocate_matrices_gw_im_time(weights_cos_tf_w_to_t, weights_sin_tf_t_to_w, do_ic_model, &
1797 do_kpoints_cubic_rpa, fm_mat_w, &
1798 t_3c_overl_int_ao_mo, t_3c_o_mo_compressed, t_3c_o_mo_ind, &
1799 t_3c_overl_int_gw_ri, t_3c_overl_int_gw_ao, &
1800 t_3c_overl_nnp_ic, t_3c_overl_nnp_ic_reflected, &
1801 mat_w, qs_env)
1802 END IF
1803
1804 END IF
1805
1806 IF (.NOT. do_im_time .AND. .NOT. do_ri_sos_laplace_mp2) CALL exchange_work%release()
1807
1808 IF (.NOT. do_ri_sos_laplace_mp2) THEN
1809 DEALLOCATE (tj)
1810 DEALLOCATE (wj)
1811 DEALLOCATE (trace_qomega)
1812 END IF
1813
1814 IF (do_im_time .OR. do_ri_sos_laplace_mp2) THEN
1815 DEALLOCATE (tau_tj)
1816 DEALLOCATE (tau_wj)
1817 END IF
1818
1819 IF (do_im_time .AND. calc_forces) THEN
1820 CALL im_time_force_release(force_data)
1821 END IF
1822
1823 IF (calc_forces .AND. .NOT. do_im_time) CALL rpa_grad_finalize(rpa_grad, mp2_env, para_env_sub, para_env, &
1824 qs_env, gd_array, color_sub, do_ri_sos_laplace_mp2, &
1825 homo, virtual)
1826
1827 CALL timestop(handle)
1828
1829 END SUBROUTINE rpa_num_int
1830
1831END MODULE rpa_main
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.
Definition bse_main.F:14
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.
Definition bse_main.F:83
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
Definition cp_fm_types.F:15
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
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
This is the start of a dbt_api, all publically needed functions are exported here....
Definition dbt_api.F:17
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.
Definition hfx_types.F:15
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public rpa_exchange_none
integer, parameter, public sigma_none
integer, parameter, public rpa_exchange_sosex
integer, parameter, public wfc_mm_style_gemm
integer, parameter, public rpa_exchange_axk
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public int_8
Definition kinds.F:54
integer, parameter, public dp
Definition kinds.F:34
Types and basic routines needed for a kpoint calculation.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_memory(mem)
Returns the total amount of memory [bytes] in use, if known, zero otherwise.
Definition machine.F:452
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
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 ...
Definition minimax_exp.F:29
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...
Definition mp2_grids.F:14
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)
...
Definition mp2_grids.F:76
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)
...
Definition mp2_grids.F:300
Routines to calculate MP2 energy with laplace approach.
Definition mp2_laplace.F:13
subroutine, public sos_mp2_postprocessing(fm_mat_q, erpa, tau_wjquad)
...
Definition mp2_laplace.F:81
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.
Definition mp2_types.F:14
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, mimic, 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, xcint_weights, 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.
Definition rpa_grad.F:13
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)
...
Definition rpa_grad.F:2029
subroutine, public rpa_grad_copy_q(fm_mat_q, rpa_grad)
...
Definition rpa_grad.F:713
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.
Definition rpa_grad.F:122
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.
Definition rpa_grad.F:168
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)
...
Definition rpa_grad.F:738
Routines to calculate image charge corrections.
Definition rpa_gw_ic.F:13
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)
...
Definition rpa_gw_ic.F:58
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].
Definition rpa_gw.F:14
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)
...
Definition rpa_gw.F:615
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)
...
Definition rpa_gw.F:350
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)
...
Definition rpa_gw.F:809
subroutine, public get_fermi_level_offset(fermi_level_offset, fermi_level_offset_input, eigenval, homo)
...
Definition rpa_gw.F:913
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)
...
Definition rpa_gw.F:204
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)
...
Definition rpa_gw.F:690
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)
...
Definition rpa_gw.F:1234
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)
...
Definition rpa_gw.F:955
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.
Definition rpa_im_time.F:13
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.
Definition rpa_main.F:17
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)
...
Definition rpa_main.F:190
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.
Definition rpa_util.F:13
subroutine, public compute_erpa_by_freq_int(dimen_ri, trace_qomega, fm_mat_q, para_env_rpa, erpa, wjquad)
...
Definition rpa_util.F:929
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)
...
Definition rpa_util.F:147
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)
...
Definition rpa_util.F:1116
subroutine, public q_trace_and_add_unit_matrix(dimen_ri, trace_qomega, fm_mat_q)
...
Definition rpa_util.F:878
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)
...
Definition rpa_util.F:661
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)
...
Definition rpa_util.F:1266
subroutine, public remove_scaling_factor_rpa(fm_mat_s, virtual, eigenval_last, homo, omega_old)
...
Definition rpa_util.F:711
All kind of helpful little routines.
Definition util.F:14
pure integer function, dimension(2), public get_limit(m, n, me)
divide m entries into n parts, return size of part me
Definition util.F:333
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Represent a complex full matrix.
keeps the information about the structure of a full matrix
represent a full matrix
Contains information about kpoints.
stores all the informations relevant to an mpi environment