(git:374b731)
Loading...
Searching...
No Matches
mp2_integrals.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines to calculate and distribute 2c- and 3c- integrals for RI
10!> \par History
11!> 06.2012 created [Mauro Del Ben]
12!> 03.2019 separated from mp2_ri_gpw [Frederick Stein]
13! **************************************************************************************************
15 USE omp_lib, ONLY: omp_get_num_threads,&
16 omp_get_thread_num
20 USE bibliography, ONLY: delben2013,&
21 cite_reference
22 USE cell_types, ONLY: cell_type,&
33 USE cp_fm_types, ONLY: cp_fm_create,&
39 USE dbcsr_api, ONLY: &
40 dbcsr_copy, dbcsr_create, dbcsr_get_info, dbcsr_multiply, dbcsr_p_type, dbcsr_release, &
41 dbcsr_release_p, dbcsr_set, dbcsr_type, dbcsr_type_no_symmetry
42 USE dbt_api, ONLY: &
43 dbt_clear, dbt_contract, dbt_copy, dbt_create, dbt_destroy, dbt_distribution_destroy, &
44 dbt_distribution_new, dbt_distribution_type, dbt_filter, dbt_get_block, dbt_get_info, &
45 dbt_get_stored_coordinates, dbt_mp_environ_pgrid, dbt_pgrid_create, dbt_pgrid_destroy, &
46 dbt_pgrid_type, dbt_put_block, dbt_reserve_blocks, dbt_split_blocks, dbt_type
50 USE hfx_types, ONLY: alloc_containers,&
53 USE input_constants, ONLY: &
60 USE kinds, ONLY: default_string_length,&
61 dp,&
62 int_8
64 USE kpoint_types, ONLY: kpoint_type
67 USE machine, ONLY: m_flush
68 USE message_passing, ONLY: mp_cart_type,&
72 USE mp2_eri_gpw, ONLY: cleanup_gpw,&
79 USE pw_env_types, ONLY: pw_env_type
82 USE pw_types, ONLY: pw_c1d_gs_type,&
103 USE util, ONLY: get_limit
104#include "./base/base_uses.f90"
105
106 IMPLICIT NONE
107
108 PRIVATE
109
110 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'mp2_integrals'
111
113
114 TYPE intermediate_matrix_type
115 TYPE(dbcsr_type) :: matrix_ia_jnu, matrix_ia_jb
116 INTEGER :: max_row_col_local = 0
117 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: local_col_row_info
118 TYPE(cp_fm_type) :: fm_BIb_jb = cp_fm_type()
119 CHARACTER(LEN=default_string_length) :: descr = ""
120 END TYPE intermediate_matrix_type
121
122CONTAINS
123
124! **************************************************************************************************
125!> \brief with ri mp2 gpw
126!> \param BIb_C ...
127!> \param BIb_C_gw ...
128!> \param BIb_C_bse_ij ...
129!> \param BIb_C_bse_ab ...
130!> \param gd_array ...
131!> \param gd_B_virtual ...
132!> \param dimen_RI ...
133!> \param dimen_RI_red ...
134!> \param qs_env ...
135!> \param para_env ...
136!> \param para_env_sub ...
137!> \param color_sub ...
138!> \param cell ...
139!> \param particle_set ...
140!> \param atomic_kind_set ...
141!> \param qs_kind_set ...
142!> \param mo_coeff ...
143!> \param fm_matrix_PQ ...
144!> \param fm_matrix_L_kpoints ...
145!> \param fm_matrix_Minv_L_kpoints ...
146!> \param fm_matrix_Minv ...
147!> \param fm_matrix_Minv_Vtrunc_Minv ...
148!> \param nmo ...
149!> \param homo ...
150!> \param mat_munu ...
151!> \param sab_orb_sub ...
152!> \param mo_coeff_o ...
153!> \param mo_coeff_v ...
154!> \param mo_coeff_all ...
155!> \param mo_coeff_gw ...
156!> \param eps_filter ...
157!> \param unit_nr ...
158!> \param mp2_memory ...
159!> \param calc_PQ_cond_num ...
160!> \param calc_forces ...
161!> \param blacs_env_sub ...
162!> \param my_do_gw ...
163!> \param do_bse ...
164!> \param gd_B_all ...
165!> \param starts_array_mc ...
166!> \param ends_array_mc ...
167!> \param starts_array_mc_block ...
168!> \param ends_array_mc_block ...
169!> \param gw_corr_lev_occ ...
170!> \param gw_corr_lev_virt ...
171!> \param do_im_time ...
172!> \param do_kpoints_cubic_RPA ...
173!> \param kpoints ...
174!> \param t_3c_M ...
175!> \param t_3c_O ...
176!> \param t_3c_O_compressed ...
177!> \param t_3c_O_ind ...
178!> \param ri_metric ...
179!> \param gd_B_occ_bse ...
180!> \param gd_B_virt_bse ...
181!> \author Mauro Del Ben
182! **************************************************************************************************
183 SUBROUTINE mp2_ri_gpw_compute_in(BIb_C, BIb_C_gw, BIb_C_bse_ij, BIb_C_bse_ab, gd_array, gd_B_virtual, &
184 dimen_RI, dimen_RI_red, qs_env, para_env, para_env_sub, color_sub, &
185 cell, particle_set, atomic_kind_set, qs_kind_set, mo_coeff, &
186 fm_matrix_PQ, fm_matrix_L_kpoints, fm_matrix_Minv_L_kpoints, &
187 fm_matrix_Minv, fm_matrix_Minv_Vtrunc_Minv, &
188 nmo, homo, mat_munu, &
189 sab_orb_sub, mo_coeff_o, mo_coeff_v, mo_coeff_all, &
190 mo_coeff_gw, eps_filter, unit_nr, &
191 mp2_memory, calc_PQ_cond_num, calc_forces, blacs_env_sub, my_do_gw, do_bse, &
192 gd_B_all, starts_array_mc, ends_array_mc, &
193 starts_array_mc_block, ends_array_mc_block, &
194 gw_corr_lev_occ, gw_corr_lev_virt, &
195 do_im_time, do_kpoints_cubic_RPA, kpoints, &
196 t_3c_M, t_3c_O, t_3c_O_compressed, t_3c_O_ind, &
197 ri_metric, gd_B_occ_bse, gd_B_virt_bse)
198
199 TYPE(three_dim_real_array), ALLOCATABLE, &
200 DIMENSION(:), INTENT(OUT) :: bib_c, bib_c_gw
201 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :), &
202 INTENT(OUT) :: bib_c_bse_ij, bib_c_bse_ab
203 TYPE(group_dist_d1_type), INTENT(OUT) :: gd_array
204 TYPE(group_dist_d1_type), ALLOCATABLE, &
205 DIMENSION(:), INTENT(OUT) :: gd_b_virtual
206 INTEGER, INTENT(OUT) :: dimen_ri, dimen_ri_red
207 TYPE(qs_environment_type), POINTER :: qs_env
208 TYPE(mp_para_env_type), POINTER :: para_env, para_env_sub
209 INTEGER, INTENT(IN) :: color_sub
210 TYPE(cell_type), POINTER :: cell
211 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
212 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
213 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
214 TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
215 TYPE(cp_fm_type), INTENT(OUT) :: fm_matrix_pq
216 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_l_kpoints, &
217 fm_matrix_minv_l_kpoints, &
218 fm_matrix_minv, &
219 fm_matrix_minv_vtrunc_minv
220 INTEGER, INTENT(IN) :: nmo
221 INTEGER, DIMENSION(:), INTENT(IN) :: homo
222 TYPE(dbcsr_p_type), INTENT(INOUT) :: mat_munu
223 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
224 INTENT(IN), POINTER :: sab_orb_sub
225 TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN) :: mo_coeff_o, mo_coeff_v, mo_coeff_all, &
226 mo_coeff_gw
227 REAL(kind=dp), INTENT(IN) :: eps_filter
228 INTEGER, INTENT(IN) :: unit_nr
229 REAL(kind=dp), INTENT(IN) :: mp2_memory
230 LOGICAL, INTENT(IN) :: calc_pq_cond_num, calc_forces
231 TYPE(cp_blacs_env_type), POINTER :: blacs_env_sub
232 LOGICAL, INTENT(IN) :: my_do_gw, do_bse
233 TYPE(group_dist_d1_type), INTENT(OUT) :: gd_b_all
234 INTEGER, ALLOCATABLE, DIMENSION(:), INTENT(OUT) :: starts_array_mc, ends_array_mc, &
235 starts_array_mc_block, &
236 ends_array_mc_block
237 INTEGER, INTENT(IN) :: gw_corr_lev_occ, gw_corr_lev_virt
238 LOGICAL, INTENT(IN) :: do_im_time, do_kpoints_cubic_rpa
239 TYPE(kpoint_type), POINTER :: kpoints
240 TYPE(dbt_type), INTENT(OUT) :: t_3c_m
241 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :), &
242 INTENT(OUT) :: t_3c_o
243 TYPE(hfx_compression_type), ALLOCATABLE, &
244 DIMENSION(:, :, :), INTENT(INOUT) :: t_3c_o_compressed
245 TYPE(block_ind_type), ALLOCATABLE, &
246 DIMENSION(:, :, :) :: t_3c_o_ind
247 TYPE(libint_potential_type), INTENT(IN) :: ri_metric
248 TYPE(group_dist_d1_type), INTENT(OUT) :: gd_b_occ_bse, gd_b_virt_bse
249
250 CHARACTER(LEN=*), PARAMETER :: routinen = 'mp2_ri_gpw_compute_in'
251
252 INTEGER :: cm, cut_memory, cut_memory_int, eri_method, gw_corr_lev_total, handle, handle2, &
253 handle4, i, i_counter, i_mem, ibasis, ispin, itmp(2), j, jcell, kcell, lll, min_bsize, &
254 my_b_all_end, my_b_all_size, my_b_all_start, my_b_occ_bse_end, my_b_occ_bse_size, &
255 my_b_occ_bse_start, my_b_virt_bse_end, my_b_virt_bse_size, my_b_virt_bse_start, &
256 my_group_l_end, my_group_l_size, my_group_l_start, n_rep, natom, ngroup, nimg, nkind, &
257 nspins, potential_type, ri_metric_type
258 INTEGER(int_8) :: nze
259 INTEGER, ALLOCATABLE, DIMENSION(:) :: dist_ao_1, dist_ao_2, dist_ri, &
260 ends_array_mc_block_int, ends_array_mc_int, my_b_size, my_b_virtual_end, &
261 my_b_virtual_start, sizes_ao, sizes_ao_split, sizes_ri, sizes_ri_split, &
262 starts_array_mc_block_int, starts_array_mc_int, virtual
263 INTEGER, DIMENSION(2, 3) :: bounds
264 INTEGER, DIMENSION(3) :: bounds_3c, pcoord, pdims, pdims_t3c, &
265 periodic
266 LOGICAL :: do_gpw, do_kpoints_from_gamma, do_svd, &
267 memory_info
268 REAL(kind=dp) :: compression_factor, cutoff_old, eps_pgf_orb, eps_pgf_orb_old, mem_for_iak, &
269 memory_3c, occ, omega_pot, rc_ang, relative_cutoff_old
270 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: e_cutoff_old
271 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: my_lrows, my_vrows
272 TYPE(cp_eri_mme_param), POINTER :: eri_param
273 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mat_munu_local_l
274 TYPE(dbt_pgrid_type) :: pgrid_t3c_m, pgrid_t3c_overl
275 TYPE(dbt_type) :: t_3c_overl_int_template, t_3c_tmp
276 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :) :: t_3c_overl_int
277 TYPE(dft_control_type), POINTER :: dft_control
278 TYPE(distribution_3d_type) :: dist_3d
279 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_ao, basis_set_ri_aux
280 TYPE(gto_basis_set_type), POINTER :: orb_basis, ri_basis
281 TYPE(intermediate_matrix_type) :: intermed_mat_bse_ab, intermed_mat_bse_ij
282 TYPE(intermediate_matrix_type), ALLOCATABLE, &
283 DIMENSION(:) :: intermed_mat, intermed_mat_gw
284 TYPE(mp_cart_type) :: mp_comm_t3c_2
285 TYPE(neighbor_list_3c_type) :: nl_3c
286 TYPE(pw_c1d_gs_type) :: pot_g, rho_g
287 TYPE(pw_env_type), POINTER :: pw_env_sub
288 TYPE(pw_poisson_type), POINTER :: poisson_env
289 TYPE(pw_pool_type), POINTER :: auxbas_pw_pool
290 TYPE(pw_r3d_rs_type) :: psi_l, rho_r
291 TYPE(section_vals_type), POINTER :: qs_section
292 TYPE(task_list_type), POINTER :: task_list_sub
293
294 CALL timeset(routinen, handle)
295
296 CALL cite_reference(delben2013)
297
298 nspins = SIZE(homo)
299
300 ALLOCATE (virtual(nspins))
301 virtual(:) = nmo - homo(:)
302 gw_corr_lev_total = gw_corr_lev_virt + gw_corr_lev_occ
303
304 eri_method = qs_env%mp2_env%eri_method
305 eri_param => qs_env%mp2_env%eri_mme_param
306 do_svd = qs_env%mp2_env%do_svd
307 potential_type = qs_env%mp2_env%potential_parameter%potential_type
308 ri_metric_type = ri_metric%potential_type
309 omega_pot = qs_env%mp2_env%potential_parameter%omega
310
311 ! whether we need gpw integrals (plus pw stuff)
312 do_gpw = (eri_method == do_eri_gpw) .OR. &
313 ((potential_type == do_potential_long .OR. ri_metric_type == do_potential_long) &
314 .AND. qs_env%mp2_env%eri_method == do_eri_os) &
315 .OR. (ri_metric_type == do_potential_id .AND. qs_env%mp2_env%eri_method == do_eri_mme)
316
317 IF (do_svd .AND. calc_forces) THEN
318 cpabort("SVD not implemented for forces.!")
319 END IF
320
321 do_kpoints_from_gamma = qs_env%mp2_env%ri_rpa_im_time%do_kpoints_from_Gamma
322 IF (do_kpoints_cubic_rpa .OR. do_kpoints_from_gamma) THEN
323 CALL get_qs_env(qs_env=qs_env, &
324 kpoints=kpoints)
325 END IF
326 IF (do_kpoints_from_gamma) THEN
327 CALL compute_kpoints(qs_env, kpoints, unit_nr)
328 END IF
329
330 IF (do_bse) THEN
331 cpassert(my_do_gw)
332 cpassert(.NOT. do_im_time)
333 ! GPW integrals have to be implemented later
334 cpassert(.NOT. (eri_method == do_eri_gpw))
335 END IF
336
337 ngroup = para_env%num_pe/para_env_sub%num_pe
338
339 ! Preparations for MME method to compute ERIs
340 IF (qs_env%mp2_env%eri_method .EQ. do_eri_mme) THEN
341 ! cell might have changed, so we need to reset parameters
342 CALL cp_eri_mme_set_params(eri_param, cell, qs_kind_set, basis_type_1="ORB", basis_type_2="RI_AUX", para_env=para_env)
343 END IF
344
345 CALL get_cell(cell=cell, periodic=periodic)
346 ! for minimax Ewald summation, full periodicity is required
347 IF (eri_method == do_eri_mme) THEN
348 cpassert(periodic(1) == 1 .AND. periodic(2) == 1 .AND. periodic(3) == 1)
349 END IF
350
351 IF (do_svd .AND. (do_kpoints_from_gamma .OR. do_kpoints_cubic_rpa)) THEN
352 cpabort("SVD with kpoints not implemented yet!")
353 END IF
354
355 CALL get_2c_integrals(qs_env, eri_method, eri_param, para_env, para_env_sub, mp2_memory, &
356 my_lrows, my_vrows, fm_matrix_pq, ngroup, color_sub, dimen_ri, dimen_ri_red, &
357 kpoints, my_group_l_size, my_group_l_start, my_group_l_end, &
358 gd_array, calc_pq_cond_num .AND. .NOT. do_svd, do_svd, &
359 qs_env%mp2_env%potential_parameter, ri_metric, &
360 fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, &
361 do_im_time, do_kpoints_from_gamma .OR. do_kpoints_cubic_rpa, qs_env%mp2_env%mp2_gpw%eps_pgf_orb_S, &
362 qs_kind_set, sab_orb_sub, calc_forces, unit_nr)
363
364 IF (unit_nr > 0) THEN
365 associate(ri_metric => qs_env%mp2_env%ri_metric)
366 SELECT CASE (ri_metric%potential_type)
368 WRITE (unit_nr, fmt="(/T3,A,T74,A)") &
369 "RI_INFO| RI metric: ", "COULOMB"
370 CASE (do_potential_short)
371 WRITE (unit_nr, fmt="(T3,A,T71,A)") &
372 "RI_INFO| RI metric: ", "SHORTRANGE"
373 WRITE (unit_nr, '(T3,A,T61,F20.10)') &
374 "RI_INFO| Omega: ", ri_metric%omega
375 rc_ang = cp_unit_from_cp2k(ri_metric%cutoff_radius, "angstrom")
376 WRITE (unit_nr, '(T3,A,T61,F20.10)') &
377 "RI_INFO| Cutoff Radius [angstrom]: ", rc_ang
378 CASE (do_potential_long)
379 WRITE (unit_nr, fmt="(T3,A,T72,A)") &
380 "RI_INFO| RI metric: ", "LONGRANGE"
381 WRITE (unit_nr, '(T3,A,T61,F20.10)') &
382 "RI_INFO| Omega: ", ri_metric%omega
383 CASE (do_potential_id)
384 WRITE (unit_nr, fmt="(T3,A,T74,A)") &
385 "RI_INFO| RI metric: ", "OVERLAP"
387 WRITE (unit_nr, fmt="(T3,A,T64,A)") &
388 "RI_INFO| RI metric: ", "TRUNCATED COULOMB"
389 rc_ang = cp_unit_from_cp2k(ri_metric%cutoff_radius, "angstrom")
390 WRITE (unit_nr, '(T3,A,T61,F20.2)') &
391 "RI_INFO| Cutoff Radius [angstrom]: ", rc_ang
392 END SELECT
393 END associate
394 END IF
395
396 IF (calc_forces .AND. .NOT. do_im_time) THEN
397 ! we need (P|Q)^(-1/2) for future use, just save it
398 ! in a fully (home made) distributed way
399 itmp = get_limit(dimen_ri, para_env_sub%num_pe, para_env_sub%mepos)
400 lll = itmp(2) - itmp(1) + 1
401 ALLOCATE (qs_env%mp2_env%ri_grad%PQ_half(lll, my_group_l_size))
402 qs_env%mp2_env%ri_grad%PQ_half(:, :) = my_lrows(itmp(1):itmp(2), 1:my_group_l_size)
403 IF (.NOT. compare_potential_types(qs_env%mp2_env%ri_metric, qs_env%mp2_env%potential_parameter)) THEN
404 ALLOCATE (qs_env%mp2_env%ri_grad%operator_half(lll, my_group_l_size))
405 qs_env%mp2_env%ri_grad%operator_half(:, :) = my_vrows(itmp(1):itmp(2), 1:my_group_l_size)
406 DEALLOCATE (my_vrows)
407 END IF
408 END IF
409
410 IF (unit_nr > 0) THEN
411 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
412 "RI_INFO| Number of auxiliary basis functions:", dimen_ri, &
413 "GENERAL_INFO| Number of basis functions:", nmo, &
414 "GENERAL_INFO| Number of occupied orbitals:", homo(1), &
415 "GENERAL_INFO| Number of virtual orbitals:", virtual(1)
416 IF (do_svd) THEN
417 WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
418 "RI_INFO| Reduced auxiliary basis set size:", dimen_ri_red
419 END IF
420
421 mem_for_iak = dimen_ri*real(sum(homo*virtual), kind=dp)*8.0_dp/(1024_dp**2)
422
423 IF (.NOT. do_im_time) THEN
424 WRITE (unit_nr, '(T3,A,T66,F11.2,A4)') 'RI_INFO| Total memory for (ia|K) integrals:', &
425 mem_for_iak, ' MiB'
426 IF (my_do_gw .AND. .NOT. do_im_time) THEN
427 mem_for_iak = dimen_ri*real(nmo, kind=dp)*gw_corr_lev_total*8.0_dp/(1024_dp**2)
428
429 WRITE (unit_nr, '(T3,A,T66,F11.2,A4)') 'RI_INFO| Total memory for G0W0-(nm|K) integrals:', &
430 mem_for_iak, ' MiB'
431 END IF
432 END IF
433 CALL m_flush(unit_nr)
434 END IF
435
436 CALL para_env%sync() ! sync to see memory output
437
438 ! in case we do imaginary time, we need the overlap tensor (alpha beta P) or trunc. Coulomb tensor
439 IF (.NOT. do_im_time) THEN
440
441 ALLOCATE (gd_b_virtual(nspins), intermed_mat(nspins))
442 ALLOCATE (my_b_virtual_start(nspins), my_b_virtual_end(nspins), my_b_size(nspins))
443 DO ispin = 1, nspins
444
445 CALL create_intermediate_matrices(intermed_mat(ispin), mo_coeff_o(ispin)%matrix, virtual(ispin), homo(ispin), &
446 trim(adjustl(cp_to_string(ispin))), blacs_env_sub, para_env_sub)
447
448 CALL create_group_dist(gd_b_virtual(ispin), para_env_sub%num_pe, virtual(ispin))
449 CALL get_group_dist(gd_b_virtual(ispin), para_env_sub%mepos, my_b_virtual_start(ispin), my_b_virtual_end(ispin), &
450 my_b_size(ispin))
451
452 END DO
453
454 ! in the case of G0W0, we need (K|nm), n,m may be occ or virt (m restricted to corrected levels)
455 IF (my_do_gw) THEN
456
457 ALLOCATE (intermed_mat_gw(nspins))
458 DO ispin = 1, nspins
459 CALL create_intermediate_matrices(intermed_mat_gw(ispin), mo_coeff_gw(ispin)%matrix, &
460 nmo, gw_corr_lev_total, &
461 "gw_"//trim(adjustl(cp_to_string(ispin))), &
462 blacs_env_sub, para_env_sub)
463
464 END DO
465
466 CALL create_group_dist(gd_b_all, para_env_sub%num_pe, nmo)
467 CALL get_group_dist(gd_b_all, para_env_sub%mepos, my_b_all_start, my_b_all_end, my_b_all_size)
468
469 IF (do_bse) THEN
470
471 ! virt x virt matrices
472 CALL create_intermediate_matrices(intermed_mat_bse_ab, mo_coeff_v(1)%matrix, virtual(1), virtual(1), &
473 "bse_ab", blacs_env_sub, para_env_sub)
474
475 CALL create_group_dist(gd_b_virt_bse, para_env_sub%num_pe, virtual(1))
476 CALL get_group_dist(gd_b_virt_bse, para_env_sub%mepos, my_b_virt_bse_start, my_b_virt_bse_end, my_b_virt_bse_size)
477
478 ! occ x occ matrices
479 CALL create_intermediate_matrices(intermed_mat_bse_ij, mo_coeff_o(1)%matrix, homo(1), homo(1), &
480 "bse_ij", blacs_env_sub, para_env_sub)
481
482 CALL create_group_dist(gd_b_occ_bse, para_env_sub%num_pe, homo(1))
483 CALL get_group_dist(gd_b_occ_bse, para_env_sub%mepos, my_b_occ_bse_start, my_b_occ_bse_end, my_b_occ_bse_size)
484
485 END IF
486 END IF
487
488 ! array that will store the (ia|K) integrals
489 ALLOCATE (bib_c(nspins))
490 DO ispin = 1, nspins
491 ALLOCATE (bib_c(ispin)%array(my_group_l_size, my_b_size(ispin), homo(ispin)))
492 bib_c(ispin)%array = 0.0_dp
493 END DO
494
495 ! in the case of GW, we also need (nm|K)
496 IF (my_do_gw) THEN
497
498 ALLOCATE (bib_c_gw(nspins))
499 DO ispin = 1, nspins
500 ALLOCATE (bib_c_gw(ispin)%array(my_group_l_size, my_b_all_size, gw_corr_lev_total))
501 bib_c_gw(ispin)%array = 0.0_dp
502 END DO
503
504 END IF
505
506 IF (do_bse) THEN
507
508 ALLOCATE (bib_c_bse_ij(my_group_l_size, my_b_occ_bse_size, homo(1)))
509 bib_c_bse_ij = 0.0_dp
510
511 ALLOCATE (bib_c_bse_ab(my_group_l_size, my_b_virt_bse_size, virtual(1)))
512 bib_c_bse_ab = 0.0_dp
513
514 END IF
515
516 CALL timeset(routinen//"_loop", handle2)
517
518 IF (eri_method == do_eri_mme .AND. &
519 (ri_metric%potential_type == do_potential_coulomb .OR. ri_metric%potential_type == do_potential_long) .OR. &
520 eri_method == do_eri_os .AND. ri_metric%potential_type == do_potential_coulomb) THEN
521
522 NULLIFY (mat_munu_local_l)
523 ALLOCATE (mat_munu_local_l(my_group_l_size))
524 DO lll = 1, my_group_l_size
525 NULLIFY (mat_munu_local_l(lll)%matrix)
526 ALLOCATE (mat_munu_local_l(lll)%matrix)
527 CALL dbcsr_copy(mat_munu_local_l(lll)%matrix, mat_munu%matrix)
528 CALL dbcsr_set(mat_munu_local_l(lll)%matrix, 0.0_dp)
529 END DO
530 CALL mp2_eri_3c_integrate(eri_param, ri_metric, para_env_sub, qs_env, &
531 first_c=my_group_l_start, last_c=my_group_l_end, &
532 mat_ab=mat_munu_local_l, &
533 basis_type_a="ORB", basis_type_b="ORB", &
534 basis_type_c="RI_AUX", &
535 sab_nl=sab_orb_sub, eri_method=eri_method)
536
537 DO ispin = 1, nspins
538 DO lll = 1, my_group_l_size
539 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu_local_l(lll), intermed_mat(ispin), &
540 bib_c(ispin)%array(lll, :, :), &
541 mo_coeff_o(ispin)%matrix, mo_coeff_v(ispin)%matrix, &
542 eps_filter, &
543 my_b_virtual_end(ispin), my_b_virtual_start(ispin))
544 END DO
545 CALL contract_b_l(bib_c(ispin)%array, my_lrows, gd_b_virtual(ispin)%sizes, &
546 gd_array%sizes, qs_env%mp2_env%eri_blksize, &
547 ngroup, color_sub, para_env, para_env_sub)
548 END DO
549
550 IF (my_do_gw) THEN
551
552 DO ispin = 1, nspins
553 DO lll = 1, my_group_l_size
554 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu_local_l(lll), intermed_mat_gw(ispin), &
555 bib_c_gw(ispin)%array(lll, :, :), &
556 mo_coeff_gw(ispin)%matrix, mo_coeff_all(ispin)%matrix, eps_filter, &
557 my_b_all_end, my_b_all_start)
558 END DO
559 CALL contract_b_l(bib_c_gw(ispin)%array, my_lrows, gd_b_all%sizes, gd_array%sizes, qs_env%mp2_env%eri_blksize, &
560 ngroup, color_sub, para_env, para_env_sub)
561 END DO
562 END IF
563
564 IF (do_bse) THEN
565
566 ! B^ab_P matrix elements for BSE
567 DO lll = 1, my_group_l_size
568 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu_local_l(lll), intermed_mat_bse_ab, &
569 bib_c_bse_ab(lll, :, :), &
570 mo_coeff_v(1)%matrix, mo_coeff_v(1)%matrix, eps_filter, &
571 my_b_all_end, my_b_all_start)
572 END DO
573 CALL contract_b_l(bib_c_bse_ab, my_lrows, gd_b_virt_bse%sizes, gd_array%sizes, qs_env%mp2_env%eri_blksize, &
574 ngroup, color_sub, para_env, para_env_sub)
575
576 ! B^ij_P matrix elements for BSE
577 DO lll = 1, my_group_l_size
578 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu_local_l(lll), intermed_mat_bse_ij, &
579 bib_c_bse_ij(lll, :, :), &
580 mo_coeff_o(1)%matrix, mo_coeff_o(1)%matrix, eps_filter, &
581 my_b_occ_bse_end, my_b_occ_bse_start)
582 END DO
583 CALL contract_b_l(bib_c_bse_ij, my_lrows, gd_b_occ_bse%sizes, gd_array%sizes, qs_env%mp2_env%eri_blksize, &
584 ngroup, color_sub, para_env, para_env_sub)
585
586 END IF
587
588 DO lll = 1, my_group_l_size
589 CALL dbcsr_release_p(mat_munu_local_l(lll)%matrix)
590 END DO
591 DEALLOCATE (mat_munu_local_l)
592
593 ELSE IF (do_gpw) THEN
594
595 CALL prepare_gpw(qs_env, dft_control, e_cutoff_old, cutoff_old, relative_cutoff_old, para_env_sub, pw_env_sub, &
596 auxbas_pw_pool, poisson_env, task_list_sub, rho_r, rho_g, pot_g, psi_l, sab_orb_sub)
597
598 DO i_counter = 1, my_group_l_size
599
600 CALL mp2_eri_3c_integrate_gpw(mo_coeff(1), psi_l, rho_g, atomic_kind_set, qs_kind_set, cell, dft_control, &
601 particle_set, pw_env_sub, my_lrows(:, i_counter), poisson_env, rho_r, pot_g, &
602 ri_metric, mat_munu, qs_env, task_list_sub)
603
604 DO ispin = 1, nspins
605 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu, intermed_mat(ispin), &
606 bib_c(ispin)%array(i_counter, :, :), &
607 mo_coeff_o(ispin)%matrix, mo_coeff_v(ispin)%matrix, eps_filter, &
608 my_b_virtual_end(ispin), my_b_virtual_start(ispin))
609
610 END DO
611
612 IF (my_do_gw) THEN
613 ! transform (K|mu nu) to (K|nm), n corresponds to corrected GW levels, m is in nmo
614 DO ispin = 1, nspins
615 CALL ao_to_mo_and_store_b(para_env_sub, mat_munu, intermed_mat_gw(ispin), &
616 bib_c_gw(ispin)%array(i_counter, :, :), &
617 mo_coeff_gw(ispin)%matrix, mo_coeff_all(ispin)%matrix, eps_filter, &
618 my_b_all_end, my_b_all_start)
619
620 END DO
621 END IF
622
623 END DO
624
625 CALL cleanup_gpw(qs_env, e_cutoff_old, cutoff_old, relative_cutoff_old, para_env_sub, pw_env_sub, &
626 task_list_sub, auxbas_pw_pool, rho_r, rho_g, pot_g, psi_l)
627 ELSE
628 cpabort("Integration method not implemented!")
629 END IF
630
631 CALL timestop(handle2)
632
633 DEALLOCATE (my_lrows)
634
635 DO ispin = 1, nspins
636 CALL release_intermediate_matrices(intermed_mat(ispin))
637 END DO
638 DEALLOCATE (intermed_mat)
639
640 IF (my_do_gw) THEN
641 DO ispin = 1, nspins
642 CALL release_intermediate_matrices(intermed_mat_gw(ispin))
643 END DO
644 DEALLOCATE (intermed_mat_gw)
645 END IF
646
647 IF (do_bse) THEN
648 CALL release_intermediate_matrices(intermed_mat_bse_ab)
649 CALL release_intermediate_matrices(intermed_mat_bse_ij)
650 END IF
651
652 ! imag. time = low-scaling SOS-MP2, RPA, GW
653 ELSE
654
655 memory_info = qs_env%mp2_env%ri_rpa_im_time%memory_info
656
657 ! we need 3 tensors:
658 ! 1) t_3c_overl_int: 3c overlap integrals, optimized for easy access to integral blocks
659 ! (atomic blocks)
660 ! 2) t_3c_O: 3c overlap integrals, optimized for contraction (split blocks)
661 ! 3) t_3c_M: tensor M, optimized for contraction
662
663 CALL get_qs_env(qs_env, natom=natom, nkind=nkind, dft_control=dft_control)
664
665 pdims_t3c = 0
666 CALL dbt_pgrid_create(para_env, pdims_t3c, pgrid_t3c_overl)
667
668 ! set up basis
669 ALLOCATE (sizes_ri(natom), sizes_ao(natom))
670 ALLOCATE (basis_set_ri_aux(nkind), basis_set_ao(nkind))
671 CALL basis_set_list_setup(basis_set_ri_aux, "RI_AUX", qs_kind_set)
672 CALL get_particle_set(particle_set, qs_kind_set, nsgf=sizes_ri, basis=basis_set_ri_aux)
673 CALL basis_set_list_setup(basis_set_ao, "ORB", qs_kind_set)
674 CALL get_particle_set(particle_set, qs_kind_set, nsgf=sizes_ao, basis=basis_set_ao)
675
676 ! make sure we use the QS%EPS_PGF_ORB
677 qs_section => section_vals_get_subs_vals(qs_env%input, "DFT%QS")
678 CALL section_vals_val_get(qs_section, "EPS_PGF_ORB", n_rep_val=n_rep)
679 IF (n_rep /= 0) THEN
680 CALL section_vals_val_get(qs_section, "EPS_PGF_ORB", r_val=eps_pgf_orb)
681 ELSE
682 CALL section_vals_val_get(qs_section, "EPS_DEFAULT", r_val=eps_pgf_orb)
683 eps_pgf_orb = sqrt(eps_pgf_orb)
684 END IF
685 eps_pgf_orb_old = dft_control%qs_control%eps_pgf_orb
686
687 DO ibasis = 1, SIZE(basis_set_ao)
688 orb_basis => basis_set_ao(ibasis)%gto_basis_set
689 CALL init_interaction_radii_orb_basis(orb_basis, eps_pgf_orb)
690 ri_basis => basis_set_ri_aux(ibasis)%gto_basis_set
691 CALL init_interaction_radii_orb_basis(ri_basis, eps_pgf_orb)
692 END DO
693
694 cut_memory_int = qs_env%mp2_env%ri_rpa_im_time%cut_memory
695 CALL create_tensor_batches(sizes_ri, cut_memory_int, starts_array_mc_int, ends_array_mc_int, &
696 starts_array_mc_block_int, ends_array_mc_block_int)
697
698 DEALLOCATE (starts_array_mc_int, ends_array_mc_int)
699
700 CALL create_3c_tensor(t_3c_overl_int_template, dist_ri, dist_ao_1, dist_ao_2, pgrid_t3c_overl, &
701 sizes_ri, sizes_ao, sizes_ao, map1=[1, 2], map2=[3], &
702 name="O (RI AO | AO)")
703
704 CALL get_qs_env(qs_env, nkind=nkind, particle_set=particle_set)
705 CALL dbt_mp_environ_pgrid(pgrid_t3c_overl, pdims, pcoord)
706 CALL mp_comm_t3c_2%create(pgrid_t3c_overl%mp_comm_2d, 3, pdims)
707 CALL distribution_3d_create(dist_3d, dist_ri, dist_ao_1, dist_ao_2, &
708 nkind, particle_set, mp_comm_t3c_2, own_comm=.true.)
709 DEALLOCATE (dist_ri, dist_ao_1, dist_ao_2)
710
711 CALL build_3c_neighbor_lists(nl_3c, basis_set_ri_aux, basis_set_ao, basis_set_ao, &
712 dist_3d, ri_metric, "RPA_3c_nl", qs_env, &
713 sym_jk=.NOT. do_kpoints_cubic_rpa, own_dist=.true.)
714
715 ! init k points
716 IF (do_kpoints_cubic_rpa) THEN
717 ! set up new kpoint type with periodic images according to eps_grid from MP2 section
718 ! instead of eps_pgf_orb from QS section
719 CALL kpoint_init_cell_index(kpoints, nl_3c%jk_list, para_env, dft_control)
720 IF (unit_nr > 0) WRITE (unit=unit_nr, fmt="(T3,A,T75,i6)") &
721 "3C_OVERLAP_INTEGRALS_INFO| Number of periodic images considered:", dft_control%nimages
722
723 nimg = dft_control%nimages
724 ELSE
725 nimg = 1
726 END IF
727
728 ALLOCATE (t_3c_overl_int(nimg, nimg))
729
730 DO i = 1, SIZE(t_3c_overl_int, 1)
731 DO j = 1, SIZE(t_3c_overl_int, 2)
732 CALL dbt_create(t_3c_overl_int_template, t_3c_overl_int(i, j))
733 END DO
734 END DO
735
736 CALL dbt_destroy(t_3c_overl_int_template)
737
738 ! split blocks to improve load balancing for tensor contraction
739 min_bsize = qs_env%mp2_env%ri_rpa_im_time%min_bsize
740
741 CALL pgf_block_sizes(atomic_kind_set, basis_set_ao, min_bsize, sizes_ao_split)
742 CALL pgf_block_sizes(atomic_kind_set, basis_set_ri_aux, min_bsize, sizes_ri_split)
743
744 pdims_t3c = 0
745 CALL dbt_pgrid_create(para_env, pdims_t3c, pgrid_t3c_m)
746
747 associate(cut_memory => qs_env%mp2_env%ri_rpa_im_time%cut_memory)
748 CALL create_tensor_batches(sizes_ao_split, cut_memory, starts_array_mc, ends_array_mc, &
749 starts_array_mc_block, ends_array_mc_block)
750 CALL create_tensor_batches(sizes_ri_split, cut_memory, &
751 qs_env%mp2_env%ri_rpa_im_time%starts_array_mc_RI, &
752 qs_env%mp2_env%ri_rpa_im_time%ends_array_mc_RI, &
753 qs_env%mp2_env%ri_rpa_im_time%starts_array_mc_block_RI, &
754 qs_env%mp2_env%ri_rpa_im_time%ends_array_mc_block_RI)
755
756 END associate
757 cut_memory = qs_env%mp2_env%ri_rpa_im_time%cut_memory
758
759 CALL create_3c_tensor(t_3c_m, dist_ri, dist_ao_1, dist_ao_2, pgrid_t3c_m, &
760 sizes_ri_split, sizes_ao_split, sizes_ao_split, &
761 map1=[1], map2=[2, 3], &
762 name="M (RI | AO AO)")
763 DEALLOCATE (dist_ri, dist_ao_1, dist_ao_2)
764 CALL dbt_pgrid_destroy(pgrid_t3c_m)
765
766 ALLOCATE (t_3c_o(SIZE(t_3c_overl_int, 1), SIZE(t_3c_overl_int, 2)))
767 ALLOCATE (t_3c_o_compressed(SIZE(t_3c_overl_int, 1), SIZE(t_3c_overl_int, 2), cut_memory))
768 ALLOCATE (t_3c_o_ind(SIZE(t_3c_overl_int, 1), SIZE(t_3c_overl_int, 2), cut_memory))
769 CALL create_3c_tensor(t_3c_o(1, 1), dist_ri, dist_ao_1, dist_ao_2, pgrid_t3c_overl, &
770 sizes_ri_split, sizes_ao_split, sizes_ao_split, &
771 map1=[1, 2], map2=[3], &
772 name="O (RI AO | AO)")
773 DEALLOCATE (dist_ri, dist_ao_1, dist_ao_2)
774 CALL dbt_pgrid_destroy(pgrid_t3c_overl)
775
776 DO i = 1, SIZE(t_3c_o, 1)
777 DO j = 1, SIZE(t_3c_o, 2)
778 IF (i > 1 .OR. j > 1) CALL dbt_create(t_3c_o(1, 1), t_3c_o(i, j))
779 END DO
780 END DO
781
782 ! build integrals in batches and copy to optimized format
783 ! note: integrals are stored in terms of atomic blocks. To avoid a memory bottleneck,
784 ! integrals are calculated in batches and copied to optimized format with subatomic blocks
785
786 DO cm = 1, cut_memory_int
787 CALL build_3c_integrals(t_3c_overl_int, &
788 qs_env%mp2_env%ri_rpa_im_time%eps_filter/2, &
789 qs_env, &
790 nl_3c, &
791 int_eps=qs_env%mp2_env%ri_rpa_im_time%eps_filter/2, &
792 basis_i=basis_set_ri_aux, &
793 basis_j=basis_set_ao, basis_k=basis_set_ao, &
794 potential_parameter=ri_metric, &
795 do_kpoints=do_kpoints_cubic_rpa, &
796 bounds_i=[starts_array_mc_block_int(cm), ends_array_mc_block_int(cm)], desymmetrize=.false.)
797 CALL timeset(routinen//"_copy_3c", handle4)
798 ! copy integral tensor t_3c_overl_int to t_3c_O tensor optimized for contraction
799 DO i = 1, SIZE(t_3c_overl_int, 1)
800 DO j = 1, SIZE(t_3c_overl_int, 2)
801
802 CALL dbt_copy(t_3c_overl_int(i, j), t_3c_o(i, j), order=[1, 3, 2], &
803 summation=.true., move_data=.true.)
804 CALL dbt_clear(t_3c_overl_int(i, j))
805 CALL dbt_filter(t_3c_o(i, j), qs_env%mp2_env%ri_rpa_im_time%eps_filter/2)
806 END DO
807 END DO
808 CALL timestop(handle4)
809 END DO
810
811 DO i = 1, SIZE(t_3c_overl_int, 1)
812 DO j = 1, SIZE(t_3c_overl_int, 2)
813 CALL dbt_destroy(t_3c_overl_int(i, j))
814 END DO
815 END DO
816 DEALLOCATE (t_3c_overl_int)
817
818 CALL timeset(routinen//"_copy_3c", handle4)
819 ! desymmetrize
820 CALL dbt_create(t_3c_o(1, 1), t_3c_tmp)
821 DO jcell = 1, nimg
822 DO kcell = 1, jcell
823 CALL dbt_copy(t_3c_o(jcell, kcell), t_3c_tmp)
824 CALL dbt_copy(t_3c_tmp, t_3c_o(kcell, jcell), order=[1, 3, 2], summation=.true., move_data=.true.)
825 CALL dbt_filter(t_3c_o(kcell, jcell), qs_env%mp2_env%ri_rpa_im_time%eps_filter)
826 END DO
827 END DO
828 DO jcell = 1, nimg
829 DO kcell = jcell + 1, nimg
830 CALL dbt_copy(t_3c_o(jcell, kcell), t_3c_tmp)
831 CALL dbt_copy(t_3c_tmp, t_3c_o(kcell, jcell), order=[1, 3, 2], summation=.false., move_data=.true.)
832 CALL dbt_filter(t_3c_o(kcell, jcell), qs_env%mp2_env%ri_rpa_im_time%eps_filter)
833 END DO
834 END DO
835
836 CALL dbt_get_info(t_3c_o(1, 1), nfull_total=bounds_3c)
837 CALL get_tensor_occupancy(t_3c_o(1, 1), nze, occ)
838 memory_3c = 0.0_dp
839
840 bounds(:, 1) = [1, bounds_3c(1)]
841 bounds(:, 3) = [1, bounds_3c(3)]
842 DO i = 1, SIZE(t_3c_o, 1)
843 DO j = 1, SIZE(t_3c_o, 2)
844 DO i_mem = 1, cut_memory
845 bounds(:, 2) = [starts_array_mc(i_mem), ends_array_mc(i_mem)]
846 CALL dbt_copy(t_3c_o(i, j), t_3c_tmp, bounds=bounds)
847
848 CALL alloc_containers(t_3c_o_compressed(i, j, i_mem), 1)
849 CALL compress_tensor(t_3c_tmp, t_3c_o_ind(i, j, i_mem)%ind, &
850 t_3c_o_compressed(i, j, i_mem), &
851 qs_env%mp2_env%ri_rpa_im_time%eps_compress, memory_3c)
852 END DO
853 CALL dbt_clear(t_3c_o(i, j))
854 END DO
855 END DO
856
857 CALL para_env%sum(memory_3c)
858
859 compression_factor = real(nze, dp)*1.0e-06*8.0_dp/memory_3c
860
861 IF (unit_nr > 0) THEN
862 WRITE (unit=unit_nr, fmt="((T3,A,T66,F11.2,A4))") &
863 "MEMORY_INFO| Memory for 3-center integrals (compressed):", memory_3c, ' MiB'
864
865 WRITE (unit=unit_nr, fmt="((T3,A,T60,F21.2))") &
866 "MEMORY_INFO| Compression factor: ", compression_factor
867 END IF
868
869 CALL dbt_destroy(t_3c_tmp)
870
871 CALL timestop(handle4)
872
873 DO ibasis = 1, SIZE(basis_set_ao)
874 orb_basis => basis_set_ao(ibasis)%gto_basis_set
875 CALL init_interaction_radii_orb_basis(orb_basis, eps_pgf_orb_old)
876 ri_basis => basis_set_ri_aux(ibasis)%gto_basis_set
877 CALL init_interaction_radii_orb_basis(ri_basis, eps_pgf_orb_old)
878 END DO
879
880 DEALLOCATE (basis_set_ri_aux, basis_set_ao)
881
882 CALL neighbor_list_3c_destroy(nl_3c)
883
884 END IF
885
886 CALL timestop(handle)
887
888 END SUBROUTINE mp2_ri_gpw_compute_in
889
890! **************************************************************************************************
891!> \brief Contract (P|ai) = (R|P) x (R|ai)
892!> \param BIb_C (R|ai)
893!> \param my_Lrows (R|P)
894!> \param sizes_B number of a (virtual) indices per subgroup process
895!> \param sizes_L number of P / R (auxiliary) indices per subgroup
896!> \param blk_size ...
897!> \param ngroup how many subgroups (NG)
898!> \param igroup subgroup color
899!> \param mp_comm communicator
900!> \param para_env_sub ...
901! **************************************************************************************************
902 SUBROUTINE contract_b_l(BIb_C, my_Lrows, sizes_B, sizes_L, blk_size, ngroup, igroup, mp_comm, para_env_sub)
903 REAL(kind=dp), DIMENSION(:, :, :), INTENT(INOUT) :: bib_c
904 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: my_lrows
905 INTEGER, DIMENSION(:), INTENT(IN) :: sizes_b, sizes_l
906 INTEGER, DIMENSION(2), INTENT(IN) :: blk_size
907 INTEGER, INTENT(IN) :: ngroup, igroup
908
909 CLASS(mp_comm_type), INTENT(IN) :: mp_comm
910 TYPE(mp_para_env_type), INTENT(IN) :: para_env_sub
911
912 CHARACTER(LEN=*), PARAMETER :: routinen = 'contract_B_L'
913 LOGICAL, PARAMETER :: debug = .false.
914
915 INTEGER :: check_proc, handle, i, iend, ii, ioff, &
916 istart, loc_a, loc_p, nblk_per_thread
917 INTEGER, ALLOCATABLE, DIMENSION(:) :: block_ind_l_p, block_ind_l_r
918 INTEGER, DIMENSION(1) :: dist_b_i, map_b_1, map_l_1, map_l_2, &
919 sizes_i
920 INTEGER, DIMENSION(2) :: map_b_2, pdims_l
921 INTEGER, DIMENSION(3) :: pdims_b
922 LOGICAL :: found
923 INTEGER, DIMENSION(ngroup) :: dist_l_p, dist_l_r
924 INTEGER, DIMENSION(para_env_sub%num_pe) :: dist_b_a
925 TYPE(dbt_distribution_type) :: dist_b, dist_l
926 TYPE(dbt_pgrid_type) :: mp_comm_b, mp_comm_l
927 TYPE(dbt_type) :: tb_in, tb_in_split, tb_out, &
928 tb_out_split, tl, tl_split
929
930 CALL timeset(routinen, handle)
931
932 sizes_i(1) = SIZE(bib_c, 3)
933
934 associate(nproc => para_env_sub%num_pe, iproc => para_env_sub%mepos, iproc_glob => mp_comm%mepos)
935
936 ! local block index for R/P and a
937 loc_p = igroup + 1; loc_a = iproc + 1
938
939 cpassert(SIZE(sizes_l) .EQ. ngroup)
940 cpassert(SIZE(sizes_b) .EQ. nproc)
941 cpassert(sizes_l(loc_p) .EQ. SIZE(bib_c, 1))
942 cpassert(sizes_l(loc_p) .EQ. SIZE(my_lrows, 2))
943 cpassert(sizes_b(loc_a) .EQ. SIZE(bib_c, 2))
944
945 ! Tensor distributions as follows:
946 ! Process grid NG x Nw
947 ! Each process has coordinates (np, nw)
948 ! tB_in: (R|ai): R distributed (np), a distributed (nw)
949 ! tB_out: (P|ai): P distributed (np), a distributed (nw)
950 ! tL: (R|P): R distributed (nw), P distributed (np)
951
952 ! define mappings between tensor index and matrix index:
953 ! (R|ai) and (P|ai):
954 map_b_1 = [1] ! index 1 (R or P) maps to 1st matrix index (np distributed)
955 map_b_2 = [2, 3] ! indices 2, 3 (a, i) map to 2nd matrix index (nw distributed)
956 ! (R|P):
957 map_l_1 = [2] ! index 2 (P) maps to 1st matrix index (np distributed)
958 map_l_2 = [1] ! index 1 (R) maps to 2nd matrix index (nw distributed)
959
960 ! derive nd process grid that is compatible with distributions and 2d process grid
961 ! (R|ai) / (P|ai) on process grid NG x Nw x 1
962 ! (R|P) on process grid NG x Nw
963 pdims_b = [ngroup, nproc, 1]
964 pdims_l = [nproc, ngroup]
965
966 CALL dbt_pgrid_create(mp_comm, pdims_b, mp_comm_b)
967 CALL dbt_pgrid_create(mp_comm, pdims_l, mp_comm_l)
968
969 ! setup distribution vectors such that distribution matches parallel data layout of BIb_C and my_Lrows
970 dist_b_i = [0]
971 dist_b_a = (/(i, i=0, nproc - 1)/)
972 dist_l_r = (/(modulo(i, nproc), i=0, ngroup - 1)/) ! R index is replicated in my_Lrows, we impose a cyclic distribution
973 dist_l_p = (/(i, i=0, ngroup - 1)/)
974
975 ! create distributions and tensors
976 CALL dbt_distribution_new(dist_b, mp_comm_b, dist_l_p, dist_b_a, dist_b_i)
977 CALL dbt_distribution_new(dist_l, mp_comm_l, dist_l_r, dist_l_p)
978
979 CALL dbt_create(tb_in, "(R|ai)", dist_b, map_b_1, map_b_2, sizes_l, sizes_b, sizes_i)
980 CALL dbt_create(tb_out, "(P|ai)", dist_b, map_b_1, map_b_2, sizes_l, sizes_b, sizes_i)
981 CALL dbt_create(tl, "(R|P)", dist_l, map_l_1, map_l_2, sizes_l, sizes_l)
982
983 IF (debug) THEN
984 ! check that tensor distribution is correct
985 CALL dbt_get_stored_coordinates(tb_in, [loc_p, loc_a, 1], check_proc)
986 cpassert(check_proc .EQ. iproc_glob)
987 END IF
988
989 ! reserve (R|ai) block
990!$OMP PARALLEL DEFAULT(NONE) SHARED(tB_in,loc_P,loc_a)
991 CALL dbt_reserve_blocks(tb_in, [loc_p], [loc_a], [1])
992!$OMP END PARALLEL
993
994 ! reserve (R|P) blocks
995 ! in my_Lrows, R index is replicated. For (R|P), we distribute quadratic blocks cyclically over
996 ! the processes in a subgroup.
997 ! There are NG blocks, so each process holds at most NG/Nw+1 blocks.
998 ALLOCATE (block_ind_l_r(ngroup/nproc + 1))
999 ALLOCATE (block_ind_l_p(ngroup/nproc + 1))
1000 block_ind_l_r(:) = 0; block_ind_l_p(:) = 0
1001 ii = 0
1002 DO i = 1, ngroup
1003 CALL dbt_get_stored_coordinates(tl, [i, loc_p], check_proc)
1004 IF (check_proc == iproc_glob) THEN
1005 ii = ii + 1
1006 block_ind_l_r(ii) = i
1007 block_ind_l_p(ii) = loc_p
1008 END IF
1009 END DO
1010
1011!TODO: Parallelize creation of block list.
1012!$OMP PARALLEL DEFAULT(NONE) SHARED(tL,block_ind_L_R,block_ind_L_P,ii) &
1013!$OMP PRIVATE(nblk_per_thread,istart,iend)
1014 nblk_per_thread = ii/omp_get_num_threads() + 1
1015 istart = omp_get_thread_num()*nblk_per_thread + 1
1016 iend = min(istart + nblk_per_thread, ii)
1017 CALL dbt_reserve_blocks(tl, block_ind_l_r(istart:iend), block_ind_l_p(istart:iend))
1018!$OMP END PARALLEL
1019
1020 ! insert (R|ai) block
1021 CALL dbt_put_block(tb_in, [loc_p, loc_a, 1], shape(bib_c), bib_c)
1022
1023 ! insert (R|P) blocks
1024 ioff = 0
1025 DO i = 1, ngroup
1026 istart = ioff + 1; iend = ioff + sizes_l(i)
1027 ioff = ioff + sizes_l(i)
1028 CALL dbt_get_stored_coordinates(tl, [i, loc_p], check_proc)
1029 IF (check_proc == iproc_glob) THEN
1030 CALL dbt_put_block(tl, [i, loc_p], [sizes_l(i), sizes_l(loc_p)], my_lrows(istart:iend, :))
1031 END IF
1032 END DO
1033 END associate
1034
1035 CALL dbt_split_blocks(tb_in, tb_in_split, [blk_size(2), blk_size(1), blk_size(1)])
1036 CALL dbt_split_blocks(tl, tl_split, [blk_size(2), blk_size(2)])
1037 CALL dbt_split_blocks(tb_out, tb_out_split, [blk_size(2), blk_size(1), blk_size(1)])
1038
1039 ! contract
1040 CALL dbt_contract(alpha=1.0_dp, tensor_1=tb_in_split, tensor_2=tl_split, &
1041 beta=0.0_dp, tensor_3=tb_out_split, &
1042 contract_1=[1], notcontract_1=[2, 3], &
1043 contract_2=[1], notcontract_2=[2], &
1044 map_1=[2, 3], map_2=[1], optimize_dist=.true.)
1045
1046 ! retrieve local block of contraction result (P|ai)
1047 CALL dbt_copy(tb_out_split, tb_out)
1048
1049 CALL dbt_get_block(tb_out, [loc_p, loc_a, 1], shape(bib_c), bib_c, found)
1050 cpassert(found)
1051
1052 ! cleanup
1053 CALL dbt_destroy(tb_in)
1054 CALL dbt_destroy(tb_in_split)
1055 CALL dbt_destroy(tb_out)
1056 CALL dbt_destroy(tb_out_split)
1057 CALL dbt_destroy(tl)
1058 CALL dbt_destroy(tl_split)
1059
1060 CALL dbt_distribution_destroy(dist_b)
1061 CALL dbt_distribution_destroy(dist_l)
1062
1063 CALL dbt_pgrid_destroy(mp_comm_b)
1064 CALL dbt_pgrid_destroy(mp_comm_l)
1065
1066 CALL timestop(handle)
1067
1068 END SUBROUTINE contract_b_l
1069
1070! **************************************************************************************************
1071!> \brief Encapsulate building of intermediate matrices matrix_ia_jnu(_beta
1072!> matrix_ia_jb(_beta),fm_BIb_jb(_beta),matrix_in_jnu(for G0W0) and
1073!> fm_BIb_all(for G0W0)
1074!> \param intermed_mat ...
1075!> \param mo_coeff_templ ...
1076!> \param size_1 ...
1077!> \param size_2 ...
1078!> \param matrix_name_2 ...
1079!> \param blacs_env_sub ...
1080!> \param para_env_sub ...
1081!> \author Jan Wilhelm
1082! **************************************************************************************************
1083 SUBROUTINE create_intermediate_matrices(intermed_mat, mo_coeff_templ, size_1, size_2, &
1084 matrix_name_2, blacs_env_sub, para_env_sub)
1085
1086 TYPE(intermediate_matrix_type), INTENT(OUT) :: intermed_mat
1087 TYPE(dbcsr_type), INTENT(INOUT) :: mo_coeff_templ
1088 INTEGER, INTENT(IN) :: size_1, size_2
1089 CHARACTER(LEN=*), INTENT(IN) :: matrix_name_2
1090 TYPE(cp_blacs_env_type), POINTER :: blacs_env_sub
1091 TYPE(mp_para_env_type), POINTER :: para_env_sub
1092
1093 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_intermediate_matrices'
1094
1095 INTEGER :: handle, ncol_local, nfullcols_total, &
1096 nfullrows_total, nrow_local
1097 INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
1098 TYPE(cp_fm_struct_type), POINTER :: fm_struct
1099
1100 CALL timeset(routinen, handle)
1101
1102 ! initialize and create the matrix (K|jnu)
1103 CALL dbcsr_create(intermed_mat%matrix_ia_jnu, template=mo_coeff_templ)
1104
1105 ! Allocate Sparse matrices: (K|jb)
1106 CALL cp_dbcsr_m_by_n_from_template(intermed_mat%matrix_ia_jb, template=mo_coeff_templ, m=size_2, n=size_1, &
1107 sym=dbcsr_type_no_symmetry)
1108
1109 ! set all to zero in such a way that the memory is actually allocated
1110 CALL dbcsr_set(intermed_mat%matrix_ia_jnu, 0.0_dp)
1111 CALL dbcsr_set(intermed_mat%matrix_ia_jb, 0.0_dp)
1112
1113 ! create the analogous of matrix_ia_jb in fm type
1114 NULLIFY (fm_struct)
1115 CALL dbcsr_get_info(intermed_mat%matrix_ia_jb, nfullrows_total=nfullrows_total, nfullcols_total=nfullcols_total)
1116 CALL cp_fm_struct_create(fm_struct, context=blacs_env_sub, nrow_global=nfullrows_total, &
1117 ncol_global=nfullcols_total, para_env=para_env_sub)
1118 CALL cp_fm_create(intermed_mat%fm_BIb_jb, fm_struct, name="fm_BIb_jb_"//matrix_name_2)
1119
1120 CALL copy_dbcsr_to_fm(intermed_mat%matrix_ia_jb, intermed_mat%fm_BIb_jb)
1121 CALL cp_fm_struct_release(fm_struct)
1122
1123 CALL cp_fm_get_info(matrix=intermed_mat%fm_BIb_jb, &
1124 nrow_local=nrow_local, &
1125 ncol_local=ncol_local, &
1126 row_indices=row_indices, &
1127 col_indices=col_indices)
1128
1129 intermed_mat%max_row_col_local = max(nrow_local, ncol_local)
1130 CALL para_env_sub%max(intermed_mat%max_row_col_local)
1131
1132 ALLOCATE (intermed_mat%local_col_row_info(0:intermed_mat%max_row_col_local, 2))
1133 intermed_mat%local_col_row_info = 0
1134 ! 0,1 nrows
1135 intermed_mat%local_col_row_info(0, 1) = nrow_local
1136 intermed_mat%local_col_row_info(1:nrow_local, 1) = row_indices(1:nrow_local)
1137 ! 0,2 ncols
1138 intermed_mat%local_col_row_info(0, 2) = ncol_local
1139 intermed_mat%local_col_row_info(1:ncol_local, 2) = col_indices(1:ncol_local)
1140
1141 intermed_mat%descr = matrix_name_2
1142
1143 CALL timestop(handle)
1144
1145 END SUBROUTINE create_intermediate_matrices
1146
1147! **************************************************************************************************
1148!> \brief Encapsulate ERI postprocessing: AO to MO transformation and store in B matrix.
1149!> \param para_env ...
1150!> \param mat_munu ...
1151!> \param intermed_mat ...
1152!> \param BIb_jb ...
1153!> \param mo_coeff_o ...
1154!> \param mo_coeff_v ...
1155!> \param eps_filter ...
1156!> \param my_B_end ...
1157!> \param my_B_start ...
1158! **************************************************************************************************
1159 SUBROUTINE ao_to_mo_and_store_b(para_env, mat_munu, intermed_mat, BIb_jb, &
1160 mo_coeff_o, mo_coeff_v, eps_filter, &
1161 my_B_end, my_B_start)
1162 TYPE(mp_para_env_type), INTENT(IN) :: para_env
1163 TYPE(dbcsr_p_type), INTENT(IN) :: mat_munu
1164 TYPE(intermediate_matrix_type), INTENT(INOUT) :: intermed_mat
1165 REAL(kind=dp), DIMENSION(:, :), INTENT(OUT) :: bib_jb
1166 TYPE(dbcsr_type), POINTER :: mo_coeff_o, mo_coeff_v
1167 REAL(kind=dp), INTENT(IN) :: eps_filter
1168 INTEGER, INTENT(IN) :: my_b_end, my_b_start
1169
1170 CHARACTER(LEN=*), PARAMETER :: routinen = 'ao_to_mo_and_store_B'
1171
1172 INTEGER :: handle
1173
1174 CALL timeset(routinen//"_mult_"//trim(intermed_mat%descr), handle)
1175
1176 CALL dbcsr_multiply("N", "N", 1.0_dp, mat_munu%matrix, mo_coeff_o, &
1177 0.0_dp, intermed_mat%matrix_ia_jnu, filter_eps=eps_filter)
1178 CALL dbcsr_multiply("T", "N", 1.0_dp, intermed_mat%matrix_ia_jnu, mo_coeff_v, &
1179 0.0_dp, intermed_mat%matrix_ia_jb, filter_eps=eps_filter)
1180 CALL timestop(handle)
1181
1182 CALL timeset(routinen//"_E_Ex_"//trim(intermed_mat%descr), handle)
1183 CALL copy_dbcsr_to_fm(intermed_mat%matrix_ia_jb, intermed_mat%fm_BIb_jb)
1184
1185 CALL grep_my_integrals(para_env, intermed_mat%fm_BIb_jb, bib_jb, intermed_mat%max_row_col_local, &
1186 intermed_mat%local_col_row_info, &
1187 my_b_end, my_b_start)
1188
1189 CALL timestop(handle)
1190 END SUBROUTINE ao_to_mo_and_store_b
1191
1192! **************************************************************************************************
1193!> \brief ...
1194!> \param intermed_mat ...
1195! **************************************************************************************************
1196 SUBROUTINE release_intermediate_matrices(intermed_mat)
1197 TYPE(intermediate_matrix_type), INTENT(INOUT) :: intermed_mat
1198
1199 CALL dbcsr_release(intermed_mat%matrix_ia_jnu)
1200 CALL dbcsr_release(intermed_mat%matrix_ia_jb)
1201 CALL cp_fm_release(intermed_mat%fm_BIb_jb)
1202 DEALLOCATE (intermed_mat%local_col_row_info)
1203
1204 END SUBROUTINE
1205
1206! **************************************************************************************************
1207!> \brief ...
1208!> \param qs_env ...
1209!> \param kpoints ...
1210!> \param unit_nr ...
1211! **************************************************************************************************
1212 SUBROUTINE compute_kpoints(qs_env, kpoints, unit_nr)
1213
1214 TYPE(qs_environment_type), POINTER :: qs_env
1215 TYPE(kpoint_type), POINTER :: kpoints
1216 INTEGER :: unit_nr
1217
1218 CHARACTER(LEN=*), PARAMETER :: routinen = 'compute_kpoints'
1219
1220 INTEGER :: handle, i, i_dim, ix, iy, iz, nkp, &
1221 nkp_extra, nkp_orig
1222 INTEGER, DIMENSION(3) :: nkp_grid, nkp_grid_extra, periodic
1223 LOGICAL :: do_extrapolate_kpoints
1224 TYPE(cell_type), POINTER :: cell
1225 TYPE(dft_control_type), POINTER :: dft_control
1226 TYPE(mp_para_env_type), POINTER :: para_env
1227 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1228 POINTER :: sab_orb
1229
1230 CALL timeset(routinen, handle)
1231
1232 NULLIFY (cell, dft_control, para_env)
1233 CALL get_qs_env(qs_env=qs_env, cell=cell, para_env=para_env, dft_control=dft_control, sab_orb=sab_orb)
1234 CALL get_cell(cell=cell, periodic=periodic)
1235
1236 ! general because we augment a Monkhorst-Pack mesh by additional points in the BZ
1237 kpoints%kp_scheme = "GENERAL"
1238 kpoints%symmetry = .false.
1239 kpoints%verbose = .false.
1240 kpoints%full_grid = .true.
1241 kpoints%use_real_wfn = .false.
1242 kpoints%eps_geo = 1.e-6_dp
1243 nkp_grid(1:3) = qs_env%mp2_env%ri_rpa_im_time%kp_grid(1:3)
1244 do_extrapolate_kpoints = qs_env%mp2_env%ri_rpa_im_time%do_extrapolate_kpoints
1245
1246 DO i_dim = 1, 3
1247 IF (periodic(i_dim) == 1) THEN
1248 cpassert(modulo(nkp_grid(i_dim), 2) == 0)
1249 END IF
1250 IF (periodic(i_dim) == 0) THEN
1251 cpassert(nkp_grid(i_dim) == 1)
1252 END IF
1253 END DO
1254
1255 nkp_orig = nkp_grid(1)*nkp_grid(2)*nkp_grid(3)/2
1256
1257 IF (do_extrapolate_kpoints) THEN
1258
1259 cpassert(qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method == kp_weights_w_uniform)
1260
1261 DO i_dim = 1, 3
1262 IF (periodic(i_dim) == 1) nkp_grid_extra(i_dim) = nkp_grid(i_dim) + 2
1263 IF (periodic(i_dim) == 0) nkp_grid_extra(i_dim) = 1
1264 END DO
1265
1266 qs_env%mp2_env%ri_rpa_im_time%kp_grid_extra(1:3) = nkp_grid_extra(1:3)
1267
1268 nkp_extra = nkp_grid_extra(1)*nkp_grid_extra(2)*nkp_grid_extra(3)/2
1269
1270 ELSE
1271
1272 nkp_grid_extra(1:3) = 0
1273 nkp_extra = 0
1274
1275 END IF
1276
1277 nkp = nkp_orig + nkp_extra
1278
1279 qs_env%mp2_env%ri_rpa_im_time%nkp_orig = nkp_orig
1280 qs_env%mp2_env%ri_rpa_im_time%nkp_extra = nkp_extra
1281
1282 ALLOCATE (kpoints%xkp(3, nkp), kpoints%wkp(nkp))
1283
1284 kpoints%nkp_grid(1:3) = nkp_grid(1:3)
1285 kpoints%nkp = nkp
1286
1287 ALLOCATE (qs_env%mp2_env%ri_rpa_im_time%wkp_V(nkp))
1288 IF (do_extrapolate_kpoints) THEN
1289 kpoints%wkp(1:nkp_orig) = 1.0_dp/real(nkp_orig, kind=dp) &
1290 /(1.0_dp - sqrt(real(nkp_extra, kind=dp)/real(nkp_orig, kind=dp)))
1291 kpoints%wkp(nkp_orig + 1:nkp) = 1.0_dp/real(nkp_extra, kind=dp) &
1292 /(1.0_dp - sqrt(real(nkp_orig, kind=dp)/real(nkp_extra, kind=dp)))
1293 qs_env%mp2_env%ri_rpa_im_time%wkp_V(1:nkp_orig) = 0.0_dp
1294 qs_env%mp2_env%ri_rpa_im_time%wkp_V(nkp_orig + 1:nkp) = 1.0_dp/real(nkp_extra, kind=dp)
1295 ELSE
1296 kpoints%wkp(:) = 1.0_dp/real(nkp, kind=dp)
1297 qs_env%mp2_env%ri_rpa_im_time%wkp_V(:) = kpoints%wkp(:)
1298 END IF
1299
1300 i = 0
1301 DO ix = 1, nkp_grid(1)
1302 DO iy = 1, nkp_grid(2)
1303 DO iz = 1, nkp_grid(3)
1304
1305 IF (i == nkp_orig) cycle
1306 i = i + 1
1307
1308 kpoints%xkp(1, i) = real(2*ix - nkp_grid(1) - 1, kind=dp)/(2._dp*real(nkp_grid(1), kind=dp))
1309 kpoints%xkp(2, i) = real(2*iy - nkp_grid(2) - 1, kind=dp)/(2._dp*real(nkp_grid(2), kind=dp))
1310 kpoints%xkp(3, i) = real(2*iz - nkp_grid(3) - 1, kind=dp)/(2._dp*real(nkp_grid(3), kind=dp))
1311
1312 END DO
1313 END DO
1314 END DO
1315
1316 DO ix = 1, nkp_grid_extra(1)
1317 DO iy = 1, nkp_grid_extra(2)
1318 DO iz = 1, nkp_grid_extra(3)
1319
1320 i = i + 1
1321 IF (i > nkp) cycle
1322
1323 kpoints%xkp(1, i) = real(2*ix - nkp_grid_extra(1) - 1, kind=dp)/(2._dp*real(nkp_grid_extra(1), kind=dp))
1324 kpoints%xkp(2, i) = real(2*iy - nkp_grid_extra(2) - 1, kind=dp)/(2._dp*real(nkp_grid_extra(2), kind=dp))
1325 kpoints%xkp(3, i) = real(2*iz - nkp_grid_extra(3) - 1, kind=dp)/(2._dp*real(nkp_grid_extra(3), kind=dp))
1326
1327 END DO
1328 END DO
1329 END DO
1330
1331 CALL kpoint_init_cell_index(kpoints, sab_orb, para_env, dft_control)
1332
1333 CALL set_qs_env(qs_env, kpoints=kpoints)
1334
1335 IF (unit_nr > 0) THEN
1336
1337 IF (do_extrapolate_kpoints) THEN
1338 WRITE (unit=unit_nr, fmt="(T3,A,T69,3I4)") "KPOINT_INFO| K-point mesh for V (leading to Sigma^x):", nkp_grid(1:3)
1339 WRITE (unit=unit_nr, fmt="(T3,A,T69)") "KPOINT_INFO| K-point extrapolation for W^c is used (W^c leads to Sigma^c):"
1340 WRITE (unit=unit_nr, fmt="(T3,A,T69,3I4)") "KPOINT_INFO| K-point mesh 1 for W^c:", nkp_grid(1:3)
1341 WRITE (unit=unit_nr, fmt="(T3,A,T69,3I4)") "KPOINT_INFO| K-point mesh 2 for W^c:", nkp_grid_extra(1:3)
1342 ELSE
1343 WRITE (unit=unit_nr, fmt="(T3,A,T69,3I4)") "KPOINT_INFO| K-point mesh for V and W:", nkp_grid(1:3)
1344 WRITE (unit=unit_nr, fmt="(T3,A,T75,I6)") "KPOINT_INFO| Number of kpoints for V and W:", nkp
1345 END IF
1346
1347 SELECT CASE (qs_env%mp2_env%ri_rpa_im_time%kpoint_weights_W_method)
1348 CASE (kp_weights_w_tailored)
1349 WRITE (unit=unit_nr, fmt="(T3,A,T81)") &
1350 "KPOINT_INFO| K-point weights for W: TAILORED"
1351 CASE (kp_weights_w_auto)
1352 WRITE (unit=unit_nr, fmt="(T3,A,T81)") &
1353 "KPOINT_INFO| K-point weights for W: AUTO"
1354 CASE (kp_weights_w_uniform)
1355 WRITE (unit=unit_nr, fmt="(T3,A,T81)") &
1356 "KPOINT_INFO| K-point weights for W: UNIFORM"
1357 END SELECT
1358
1359 END IF
1360
1361 CALL timestop(handle)
1362
1363 END SUBROUTINE compute_kpoints
1364
1365! **************************************************************************************************
1366!> \brief ...
1367!> \param para_env_sub ...
1368!> \param fm_BIb_jb ...
1369!> \param BIb_jb ...
1370!> \param max_row_col_local ...
1371!> \param local_col_row_info ...
1372!> \param my_B_virtual_end ...
1373!> \param my_B_virtual_start ...
1374! **************************************************************************************************
1375 SUBROUTINE grep_my_integrals(para_env_sub, fm_BIb_jb, BIb_jb, max_row_col_local, &
1376 local_col_row_info, &
1377 my_B_virtual_end, my_B_virtual_start)
1378 TYPE(mp_para_env_type), INTENT(IN) :: para_env_sub
1379 TYPE(cp_fm_type), INTENT(IN) :: fm_bib_jb
1380 REAL(kind=dp), DIMENSION(:, :), INTENT(OUT) :: bib_jb
1381 INTEGER, INTENT(IN) :: max_row_col_local
1382 INTEGER, ALLOCATABLE, DIMENSION(:, :), INTENT(IN) :: local_col_row_info
1383 INTEGER, INTENT(IN) :: my_b_virtual_end, my_b_virtual_start
1384
1385 INTEGER :: i_global, iib, j_global, jjb, ncol_rec, &
1386 nrow_rec, proc_receive, proc_send, &
1387 proc_shift
1388 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: rec_col_row_info
1389 INTEGER, DIMENSION(:), POINTER :: col_indices_rec, row_indices_rec
1390 REAL(kind=dp), DIMENSION(:, :), POINTER :: local_bi, rec_bi
1391
1392 ALLOCATE (rec_col_row_info(0:max_row_col_local, 2))
1393
1394 rec_col_row_info(:, :) = local_col_row_info
1395
1396 nrow_rec = rec_col_row_info(0, 1)
1397 ncol_rec = rec_col_row_info(0, 2)
1398
1399 ALLOCATE (row_indices_rec(nrow_rec))
1400 row_indices_rec = rec_col_row_info(1:nrow_rec, 1)
1401
1402 ALLOCATE (col_indices_rec(ncol_rec))
1403 col_indices_rec = rec_col_row_info(1:ncol_rec, 2)
1404
1405 ! accumulate data on BIb_jb buffer starting from myself
1406 DO jjb = 1, ncol_rec
1407 j_global = col_indices_rec(jjb)
1408 IF (j_global >= my_b_virtual_start .AND. j_global <= my_b_virtual_end) THEN
1409 DO iib = 1, nrow_rec
1410 i_global = row_indices_rec(iib)
1411 bib_jb(j_global - my_b_virtual_start + 1, i_global) = fm_bib_jb%local_data(iib, jjb)
1412 END DO
1413 END IF
1414 END DO
1415
1416 DEALLOCATE (row_indices_rec)
1417 DEALLOCATE (col_indices_rec)
1418
1419 IF (para_env_sub%num_pe > 1) THEN
1420 ALLOCATE (local_bi(nrow_rec, ncol_rec))
1421 local_bi(1:nrow_rec, 1:ncol_rec) = fm_bib_jb%local_data(1:nrow_rec, 1:ncol_rec)
1422
1423 DO proc_shift = 1, para_env_sub%num_pe - 1
1424 proc_send = modulo(para_env_sub%mepos + proc_shift, para_env_sub%num_pe)
1425 proc_receive = modulo(para_env_sub%mepos - proc_shift, para_env_sub%num_pe)
1426
1427 ! first exchange information on the local data
1428 rec_col_row_info = 0
1429 CALL para_env_sub%sendrecv(local_col_row_info, proc_send, rec_col_row_info, proc_receive)
1430 nrow_rec = rec_col_row_info(0, 1)
1431 ncol_rec = rec_col_row_info(0, 2)
1432
1433 ALLOCATE (row_indices_rec(nrow_rec))
1434 row_indices_rec = rec_col_row_info(1:nrow_rec, 1)
1435
1436 ALLOCATE (col_indices_rec(ncol_rec))
1437 col_indices_rec = rec_col_row_info(1:ncol_rec, 2)
1438
1439 ALLOCATE (rec_bi(nrow_rec, ncol_rec))
1440 rec_bi = 0.0_dp
1441
1442 ! then send and receive the real data
1443 CALL para_env_sub%sendrecv(local_bi, proc_send, rec_bi, proc_receive)
1444
1445 ! accumulate the received data on BIb_jb buffer
1446 DO jjb = 1, ncol_rec
1447 j_global = col_indices_rec(jjb)
1448 IF (j_global >= my_b_virtual_start .AND. j_global <= my_b_virtual_end) THEN
1449 DO iib = 1, nrow_rec
1450 i_global = row_indices_rec(iib)
1451 bib_jb(j_global - my_b_virtual_start + 1, i_global) = rec_bi(iib, jjb)
1452 END DO
1453 END IF
1454 END DO
1455
1456 DEALLOCATE (col_indices_rec)
1457 DEALLOCATE (row_indices_rec)
1458 DEALLOCATE (rec_bi)
1459 END DO
1460
1461 DEALLOCATE (local_bi)
1462 END IF
1463
1464 DEALLOCATE (rec_col_row_info)
1465
1466 END SUBROUTINE grep_my_integrals
1467
1468END MODULE mp2_integrals
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Define the atomic kind types and their sub types.
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public delben2013
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_m_by_n_from_template(matrix, template, m, n, sym, data_type)
Utility function to create an arbitrary shaped dbcsr matrix with the same processor grid as the templ...
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
Interface to Minimax-Ewald method for periodic ERI's to be used in CP2K.
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_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp)
creates a new full matrix with the given structure
various routines to log and control the output. The idea is that decisions about where to log should ...
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_from_cp2k(value, unit_str, defaults, power)
converts from the internal cp2k units to the given unit
Definition cp_units.F:1179
This is the start of a dbt_api, all publically needed functions are exported here....
Definition dbt_api.F:17
Types to describe group distributions.
Types and set/get functions for HFX.
Definition hfx_types.F:15
subroutine, public alloc_containers(data, bin_size)
...
Definition hfx_types.F:2906
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_eri_os
integer, parameter, public kp_weights_w_auto
integer, parameter, public kp_weights_w_uniform
integer, parameter, public do_eri_mme
integer, parameter, public do_potential_truncated
integer, parameter, public do_potential_id
integer, parameter, public do_potential_coulomb
integer, parameter, public do_potential_short
integer, parameter, public kp_weights_w_tailored
integer, parameter, public do_potential_long
integer, parameter, public do_eri_gpw
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
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
integer, parameter, public default_string_length
Definition kinds.F:57
Routines needed for kpoint calculation.
subroutine, public kpoint_init_cell_index(kpoint, sab_nl, para_env, dft_control)
Generates the mapping of cell indices and linear RS index CELL (0,0,0) is always mapped to index 1.
Types and basic routines needed for a kpoint calculation.
2- and 3-center electron repulsion integral routines based on libint2 Currently available operators: ...
pure logical function, public compare_potential_types(potential1, potential2)
Helper function to compare libint_potential_types.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:106
Interface to the message passing library MPI.
Routines to calculate 2c- and 3c-integrals for RI with GPW.
Definition mp2_eri_gpw.F:13
subroutine, public prepare_gpw(qs_env, dft_control, e_cutoff_old, cutoff_old, relative_cutoff_old, para_env_sub, pw_env_sub, auxbas_pw_pool, poisson_env, task_list_sub, rho_r, rho_g, pot_g, psi_l, sab_orb_sub)
Prepares GPW calculation for RI-MP2/RI-RPA.
subroutine, public mp2_eri_3c_integrate_gpw(mo_coeff, psi_l, rho_g, atomic_kind_set, qs_kind_set, cell, dft_control, particle_set, pw_env_sub, external_vector, poisson_env, rho_r, pot_g, potential_parameter, mat_munu, qs_env, task_list_sub)
...
subroutine, public cleanup_gpw(qs_env, e_cutoff_old, cutoff_old, relative_cutoff_old, para_env_sub, pw_env_sub, task_list_sub, auxbas_pw_pool, rho_r, rho_g, pot_g, psi_l)
Cleanup GPW integration for RI-MP2/RI-RPA.
Interface to direct methods for electron repulsion integrals for MP2.
Definition mp2_eri.F:12
subroutine, public mp2_eri_3c_integrate(param, potential_parameter, para_env, qs_env, first_c, last_c, mat_ab, basis_type_a, basis_type_b, basis_type_c, sab_nl, eri_method, pabc, force_a, force_b, force_c, mat_dabc, mat_adbc, mat_abdc)
high-level integration routine for 3c integrals (ab|c) over CP2K basis sets. For each local function ...
Definition mp2_eri.F:680
Routines to calculate and distribute 2c- and 3c- integrals for RI.
subroutine, public compute_kpoints(qs_env, kpoints, unit_nr)
...
subroutine, public mp2_ri_gpw_compute_in(bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, gd_array, gd_b_virtual, dimen_ri, dimen_ri_red, qs_env, para_env, para_env_sub, color_sub, cell, particle_set, atomic_kind_set, qs_kind_set, mo_coeff, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, nmo, homo, mat_munu, sab_orb_sub, mo_coeff_o, mo_coeff_v, mo_coeff_all, mo_coeff_gw, eps_filter, unit_nr, mp2_memory, calc_pq_cond_num, calc_forces, blacs_env_sub, my_do_gw, do_bse, gd_b_all, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, gw_corr_lev_occ, gw_corr_lev_virt, do_im_time, do_kpoints_cubic_rpa, kpoints, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, ri_metric, gd_b_occ_bse, gd_b_virt_bse)
with ri mp2 gpw
Framework for 2c-integrals for RI.
Definition mp2_ri_2c.F:14
subroutine, public get_2c_integrals(qs_env, eri_method, eri_param, para_env, para_env_sub, mp2_memory, my_lrows, my_vrows, fm_matrix_pq, ngroup, color_sub, dimen_ri, dimen_ri_red, kpoints, my_group_l_size, my_group_l_start, my_group_l_end, gd_array, calc_pq_cond_num, do_svd, potential, ri_metric, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, do_im_time, do_kpoints, mp2_eps_pgf_orb_s, qs_kind_set, sab_orb_sub, calc_forces, unit_nr)
...
Definition mp2_ri_2c.F:162
Types needed for MP2 calculations.
Definition mp2_types.F:14
Define methods related to particle_type.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis)
Get the components of a particle set.
Define the data structure for the particle information.
container for various plainwaves related things
functions related to the poisson solver on regular grids
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
subroutine, public set_qs_env(qs_env, super_cell, mos, qmmm, qmmm_periodic, ewald_env, ewald_pw, mpools, rho_external, external_vxc, mask, scf_control, rel_control, qs_charges, ks_env, ks_qmmm_env, wf_history, scf_env, active_space, input, oce, rho_atom_set, rho0_atom_set, rho0_mpole, run_rtp, rtp, rhoz_set, rhoz_tot, ecoul_1c, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, efield, linres_control, xas_env, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, ls_scf_env, do_transport, transport_env, lri_env, lri_density, exstate_env, ec_env, dispersion_env, gcp_env, mp2_env, bs_env, kg_env, force, kpoints, wanniercentres, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, rhs)
Set the QUICKSTEP environment.
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_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, 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, rhs)
Get the QUICKSTEP environment.
Some utility functions for the calculation of integrals.
subroutine, public basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
Set up an easy accessible list of the basis sets for all kinds.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public init_interaction_radii_orb_basis(orb_basis_set, eps_pgf_orb, eps_pgf_short)
...
Define the quickstep kind type and their sub types.
Define the neighbor list data types and the corresponding functionality.
Utility methods to build 3-center integral tensors of various types.
subroutine, public distribution_3d_create(dist_3d, dist1, dist2, dist3, nkind, particle_set, mp_comm_3d, own_comm)
Create a 3d distribution.
subroutine, public pgf_block_sizes(atomic_kind_set, basis, min_blk_size, pgf_blk_sizes)
...
subroutine, public create_tensor_batches(sizes, nbatches, starts_array, ends_array, starts_array_block, ends_array_block)
...
subroutine, public create_3c_tensor(t3c, dist_1, dist_2, dist_3, pgrid, sizes_1, sizes_2, sizes_3, map1, map2, name)
...
Utility methods to build 3-center integral tensors of various types.
Definition qs_tensors.F:11
subroutine, public compress_tensor(tensor, blk_indices, compressed, eps, memory)
...
subroutine, public neighbor_list_3c_destroy(ijk_list)
Destroy 3c neighborlist.
Definition qs_tensors.F:383
subroutine, public get_tensor_occupancy(tensor, nze, occ)
...
subroutine, public build_3c_neighbor_lists(ijk_list, basis_i, basis_j, basis_k, dist_3d, potential_parameter, name, qs_env, sym_ij, sym_jk, sym_ik, molecular, op_pos, own_dist)
Build a 3-center neighbor list.
Definition qs_tensors.F:282
subroutine, public build_3c_integrals(t3c, filter_eps, qs_env, nl_3c, basis_i, basis_j, basis_k, potential_parameter, int_eps, op_pos, do_kpoints, do_hfx_kpoints, desymmetrize, bounds_i, bounds_j, bounds_k, ri_range, img_to_ri_cell)
Build 3-center integral tensor.
types for task lists
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
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
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
contained for different pw related things
environment for the poisson solver
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.