Line data Source code
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
10 : !> \par History
11 : !> 01.2026 Maximilian Graml: add more bounds to exploit sparsity in 3c integrals, fixes
12 : !> \author Jan Wilhelm
13 : !> \date 07.2023
14 : ! **************************************************************************************************
15 : MODULE gw_utils
16 : USE atomic_kind_types, ONLY: atomic_kind_type,&
17 : get_atomic_kind_set
18 : USE basis_set_types, ONLY: get_gto_basis_set,&
19 : gto_basis_set_type
20 : USE bibliography, ONLY: Graml2024,&
21 : cite_reference
22 : USE cell_types, ONLY: cell_type,&
23 : pbc,&
24 : scaled_to_real
25 : USE cp_blacs_env, ONLY: cp_blacs_env_create,&
26 : cp_blacs_env_release,&
27 : cp_blacs_env_type
28 : USE cp_cfm_types, ONLY: cp_cfm_create,&
29 : cp_cfm_release,&
30 : cp_cfm_to_cfm,&
31 : cp_cfm_to_fm,&
32 : cp_cfm_type
33 : USE cp_control_types, ONLY: dft_control_type
34 : USE cp_dbcsr_api, ONLY: &
35 : dbcsr_create, dbcsr_distribution_release, dbcsr_distribution_type, dbcsr_p_type, &
36 : dbcsr_release, dbcsr_set, dbcsr_type, dbcsr_type_no_symmetry, dbcsr_type_symmetric
37 : USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
38 : copy_fm_to_dbcsr,&
39 : cp_dbcsr_dist2d_to_dist,&
40 : dbcsr_allocate_matrix_set,&
41 : dbcsr_deallocate_matrix_set
42 : USE cp_files, ONLY: close_file,&
43 : open_file
44 : USE cp_fm_basic_linalg, ONLY: cp_fm_scale_and_add
45 : USE cp_fm_struct, ONLY: cp_fm_struct_create,&
46 : cp_fm_struct_release,&
47 : cp_fm_struct_type
48 : USE cp_fm_types, ONLY: cp_fm_create,&
49 : cp_fm_get_diag,&
50 : cp_fm_release,&
51 : cp_fm_set_all,&
52 : cp_fm_type
53 : USE cp_log_handling, ONLY: cp_get_default_logger,&
54 : cp_logger_type
55 : USE cp_output_handling, ONLY: cp_print_key_generate_filename
56 : USE dbt_api, ONLY: &
57 : dbt_clear, dbt_create, dbt_destroy, dbt_filter, dbt_iterator_blocks_left, &
58 : dbt_iterator_next_block, dbt_iterator_start, dbt_iterator_stop, dbt_iterator_type, &
59 : dbt_mp_environ_pgrid, dbt_pgrid_create, dbt_pgrid_destroy, dbt_pgrid_type, dbt_type
60 : USE distribution_2d_types, ONLY: distribution_2d_type
61 : USE gw_communication, ONLY: fm_to_local_array
62 : USE gw_integrals, ONLY: build_3c_integral_block
63 : USE input_constants, ONLY: do_potential_truncated,&
64 : large_cell_Gamma,&
65 : large_cell_Gamma_ri_rs,&
66 : ri_rpa_g0w0_crossing_newton,&
67 : rtp_method_bse,&
68 : small_cell_full_kp,&
69 : xc_none
70 : USE input_section_types, ONLY: section_vals_get,&
71 : section_vals_get_subs_vals,&
72 : section_vals_type,&
73 : section_vals_val_get,&
74 : section_vals_val_set
75 : USE kinds, ONLY: default_path_length,&
76 : dp,&
77 : int_8
78 : USE kpoint_k_r_trafo_simple, ONLY: rs_to_kp
79 : USE kpoint_types, ONLY: get_kpoint_info,&
80 : kpoint_create,&
81 : kpoint_type
82 : USE libint_2c_3c, ONLY: libint_potential_type
83 : USE libint_wrapper, ONLY: cp_libint_static_cleanup,&
84 : cp_libint_static_init
85 : USE machine, ONLY: m_memory,&
86 : m_walltime
87 : USE mathconstants, ONLY: gaussi,&
88 : z_one,&
89 : z_zero
90 : USE mathlib, ONLY: diag_complex,&
91 : gcd
92 : USE message_passing, ONLY: mp_cart_type,&
93 : mp_para_env_type
94 : USE minimax_exp, ONLY: get_exp_minimax_coeff
95 : USE minimax_exp_gw, ONLY: get_exp_minimax_coeff_gw
96 : USE minimax_rpa, ONLY: get_rpa_minimax_coeff,&
97 : get_rpa_minimax_coeff_larger_grid
98 : USE mp2_gpw, ONLY: create_mat_munu
99 : USE mp2_grids, ONLY: get_l_sq_wghts_cos_tf_t_to_w,&
100 : get_l_sq_wghts_cos_tf_w_to_t,&
101 : get_l_sq_wghts_sin_tf_t_to_w
102 : USE mp2_ri_2c, ONLY: trunc_coulomb_for_exchange
103 : USE parallel_gemm_api, ONLY: parallel_gemm
104 : USE particle_methods, ONLY: get_particle_set
105 : USE particle_types, ONLY: particle_type
106 : USE physcon, ONLY: angstrom,&
107 : evolt
108 : USE post_scf_bandstructure_types, ONLY: post_scf_bandstructure_type
109 : USE post_scf_bandstructure_utils, ONLY: rsmat_to_kp
110 : USE qs_energy_types, ONLY: qs_energy_type
111 : USE qs_environment_types, ONLY: get_qs_env,&
112 : qs_env_part_release,&
113 : qs_environment_type
114 : USE qs_integral_utils, ONLY: basis_set_list_setup
115 : USE qs_interactions, ONLY: init_interaction_radii_orb_basis
116 : USE qs_kind_types, ONLY: get_qs_kind,&
117 : qs_kind_type
118 : USE qs_ks_methods, ONLY: qs_ks_build_kohn_sham_matrix
119 : USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type,&
120 : release_neighbor_list_sets
121 : USE qs_tensors, ONLY: build_2c_integrals,&
122 : build_2c_neighbor_lists,&
123 : build_3c_integrals,&
124 : build_3c_neighbor_lists,&
125 : get_tensor_occupancy,&
126 : neighbor_list_3c_destroy
127 : USE qs_tensors_types, ONLY: create_2c_tensor,&
128 : create_3c_tensor,&
129 : distribution_3d_create,&
130 : distribution_3d_type,&
131 : neighbor_list_3c_type
132 : USE rpa_gw, ONLY: continuation_pade
133 : #include "base/base_uses.f90"
134 :
135 : IMPLICIT NONE
136 :
137 : PRIVATE
138 :
139 : PUBLIC :: create_and_init_bs_env_for_gw, de_init_bs_env, get_i_j_atoms, &
140 : compute_xkp, time_to_freq, analyt_conti_and_print, &
141 : add_R, is_cell_in_index_to_cell, get_V_tr_R, power
142 :
143 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'gw_utils'
144 :
145 : CONTAINS
146 :
147 : ! **************************************************************************************************
148 : !> \brief ...
149 : !> \param qs_env ...
150 : !> \param bs_env ...
151 : !> \param bs_sec ...
152 : ! **************************************************************************************************
153 42 : SUBROUTINE create_and_init_bs_env_for_gw(qs_env, bs_env, bs_sec)
154 : TYPE(qs_environment_type), POINTER :: qs_env
155 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
156 : TYPE(section_vals_type), POINTER :: bs_sec
157 :
158 : CHARACTER(LEN=*), PARAMETER :: routineN = 'create_and_init_bs_env_for_gw'
159 :
160 : INTEGER :: handle
161 :
162 42 : CALL timeset(routineN, handle)
163 :
164 42 : CALL cite_reference(Graml2024)
165 :
166 42 : CALL read_gw_input_parameters(bs_env, bs_sec)
167 :
168 42 : CALL print_header_and_input_parameters(bs_env)
169 :
170 42 : CALL setup_AO_and_RI_basis_set(qs_env, bs_env)
171 :
172 42 : CALL get_RI_basis_and_basis_function_indices(qs_env, bs_env)
173 :
174 42 : CALL set_heuristic_parameters(bs_env, qs_env)
175 :
176 42 : CALL cp_libint_static_init()
177 :
178 42 : CALL setup_kpoints_chi_eps_W(bs_env, bs_env%kpoints_chi_eps_W)
179 :
180 42 : IF (bs_env%small_cell_full_kp_or_large_cell_Gamma == small_cell_full_kp) THEN
181 8 : CALL setup_cells_3c(qs_env, bs_env)
182 : END IF
183 :
184 42 : CALL set_parallelization_parameters(qs_env, bs_env)
185 :
186 42 : CALL allocate_matrices(qs_env, bs_env)
187 :
188 42 : CALL compute_V_xc(qs_env, bs_env)
189 :
190 42 : CALL create_tensors(qs_env, bs_env)
191 :
192 76 : SELECT CASE (bs_env%small_cell_full_kp_or_large_cell_Gamma)
193 : CASE (large_cell_Gamma, large_cell_Gamma_ri_rs)
194 :
195 34 : CALL allocate_GW_eigenvalues(bs_env)
196 :
197 34 : CALL check_sparsity_3c(qs_env, bs_env)
198 :
199 34 : CALL set_sparsity_parallelization_parameters(bs_env)
200 :
201 34 : CALL check_for_restart_files(qs_env, bs_env)
202 :
203 : CASE (small_cell_full_kp)
204 :
205 8 : CALL compute_3c_integrals(qs_env, bs_env)
206 :
207 8 : CALL setup_cells_Delta_R(bs_env)
208 :
209 8 : CALL setup_parallelization_Delta_R(bs_env)
210 :
211 8 : CALL allocate_matrices_small_cell_full_kp(qs_env, bs_env)
212 :
213 8 : CALL trafo_V_xc_R_to_kp(qs_env, bs_env)
214 :
215 50 : CALL heuristic_RI_regularization(qs_env, bs_env)
216 :
217 : END SELECT
218 :
219 42 : CALL setup_time_and_frequency_minimax_grid(bs_env)
220 :
221 : ! free memory in qs_env; only if one is not calculating the LDOS because
222 : ! we need real-space grid operations in pw_env, task_list for the LDOS
223 : ! Recommendation in case of memory issues: first perform GW calculation without calculating
224 : ! LDOS (to safe memor). Then, use GW restart files
225 : ! in a subsequent calculation to calculate the LDOS
226 : ! Marek : TODO - boolean that does not interfere with RTP init but sets this to correct value
227 : IF (.NOT. bs_env%do_ldos .AND. .FALSE.) THEN
228 : CALL qs_env_part_release(qs_env)
229 : END IF
230 :
231 42 : CALL timestop(handle)
232 :
233 42 : END SUBROUTINE create_and_init_bs_env_for_gw
234 :
235 : ! **************************************************************************************************
236 : !> \brief ...
237 : !> \param bs_env ...
238 : ! **************************************************************************************************
239 42 : SUBROUTINE de_init_bs_env(bs_env)
240 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
241 :
242 : CHARACTER(LEN=*), PARAMETER :: routineN = 'de_init_bs_env'
243 :
244 : INTEGER :: handle
245 :
246 42 : CALL timeset(routineN, handle)
247 : ! deallocate quantities here which:
248 : ! 1. cannot be deallocated in bs_env_release due to circular dependencies
249 : ! 2. consume a lot of memory and should not be kept until the quantity is
250 : ! deallocated in bs_env_release
251 :
252 42 : IF (ASSOCIATED(bs_env%nl_3c%ij_list) .AND. (bs_env%rtp_method == rtp_method_bse)) THEN
253 14 : IF (bs_env%unit_nr > 0) WRITE (bs_env%unit_nr, *) "Retaining nl_3c for RTBSE"
254 : ELSE
255 28 : CALL neighbor_list_3c_destroy(bs_env%nl_3c)
256 : END IF
257 :
258 42 : CALL cp_libint_static_cleanup()
259 :
260 42 : CALL timestop(handle)
261 :
262 42 : END SUBROUTINE de_init_bs_env
263 :
264 : ! **************************************************************************************************
265 : !> \brief ...
266 : !> \param bs_env ...
267 : !> \param bs_sec ...
268 : ! **************************************************************************************************
269 42 : SUBROUTINE read_gw_input_parameters(bs_env, bs_sec)
270 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
271 : TYPE(section_vals_type), POINTER :: bs_sec
272 :
273 : CHARACTER(LEN=*), PARAMETER :: routineN = 'read_gw_input_parameters'
274 :
275 : INTEGER :: handle
276 : TYPE(section_vals_type), POINTER :: gw_sec
277 :
278 42 : CALL timeset(routineN, handle)
279 :
280 42 : NULLIFY (gw_sec)
281 42 : gw_sec => section_vals_get_subs_vals(bs_sec, "GW")
282 :
283 42 : CALL section_vals_val_get(gw_sec, "NUM_TIME_FREQ_POINTS", i_val=bs_env%num_time_freq_points)
284 42 : CALL section_vals_val_get(gw_sec, "EPS_FILTER", r_val=bs_env%eps_filter)
285 42 : CALL section_vals_val_get(gw_sec, "REGULARIZATION_RI", r_val=bs_env%input_regularization_RI)
286 42 : CALL section_vals_val_get(gw_sec, "REGULARIZATION_MINIMAX", r_val=bs_env%input_regularization_minimax)
287 42 : CALL section_vals_val_get(gw_sec, "CUTOFF_RADIUS_RI", r_val=bs_env%ri_metric%cutoff_radius)
288 42 : CALL section_vals_val_get(gw_sec, "MEMORY_PER_PROC", r_val=bs_env%input_memory_per_proc_GB)
289 42 : CALL section_vals_val_get(gw_sec, "APPROX_KP_EXTRAPOL", l_val=bs_env%approx_kp_extrapol)
290 42 : CALL section_vals_val_get(gw_sec, "SIZE_LATTICE_SUM", i_val=bs_env%size_lattice_sum_V)
291 42 : CALL section_vals_val_get(gw_sec, "KPOINTS_W", i_vals=bs_env%nkp_grid_chi_eps_W_input)
292 42 : CALL section_vals_val_get(gw_sec, "HEDIN_SHIFT", l_val=bs_env%do_hedin_shift)
293 42 : CALL section_vals_val_get(gw_sec, "FREQ_MAX_FIT", r_val=bs_env%freq_max_fit)
294 42 : CALL section_vals_val_get(gw_sec, "PRINT%PRINT_DBT_CONTRACT", l_val=bs_env%print_contract)
295 42 : CALL section_vals_val_get(gw_sec, "PRINT%PRINT_DBT_CONTRACT_VERBOSE", l_val=bs_env%print_contract_verbose)
296 :
297 42 : IF (bs_env%print_contract) THEN
298 0 : bs_env%unit_nr_contract = bs_env%unit_nr
299 : ELSE
300 42 : bs_env%unit_nr_contract = 0
301 : END IF
302 42 : CALL timestop(handle)
303 :
304 42 : END SUBROUTINE read_gw_input_parameters
305 :
306 : ! **************************************************************************************************
307 : !> \brief ...
308 : !> \param qs_env ...
309 : !> \param bs_env ...
310 : ! **************************************************************************************************
311 42 : SUBROUTINE setup_AO_and_RI_basis_set(qs_env, bs_env)
312 : TYPE(qs_environment_type), POINTER :: qs_env
313 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
314 :
315 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_AO_and_RI_basis_set'
316 :
317 : INTEGER :: handle, natom, nkind
318 42 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
319 42 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
320 :
321 42 : CALL timeset(routineN, handle)
322 :
323 : CALL get_qs_env(qs_env, &
324 : qs_kind_set=qs_kind_set, &
325 : particle_set=particle_set, &
326 42 : natom=natom, nkind=nkind)
327 :
328 : ! set up basis
329 168 : ALLOCATE (bs_env%sizes_RI(natom), bs_env%sizes_AO(natom))
330 304 : ALLOCATE (bs_env%basis_set_RI(nkind), bs_env%basis_set_AO(nkind))
331 :
332 42 : CALL basis_set_list_setup(bs_env%basis_set_RI, "RI_AUX", qs_kind_set)
333 42 : CALL basis_set_list_setup(bs_env%basis_set_AO, "ORB", qs_kind_set)
334 :
335 : CALL get_particle_set(particle_set, qs_kind_set, nsgf=bs_env%sizes_RI, &
336 42 : basis=bs_env%basis_set_RI)
337 : CALL get_particle_set(particle_set, qs_kind_set, nsgf=bs_env%sizes_AO, &
338 42 : basis=bs_env%basis_set_AO)
339 :
340 42 : CALL timestop(handle)
341 :
342 42 : END SUBROUTINE setup_AO_and_RI_basis_set
343 :
344 : ! **************************************************************************************************
345 : !> \brief ...
346 : !> \param qs_env ...
347 : !> \param bs_env ...
348 : ! **************************************************************************************************
349 42 : SUBROUTINE get_RI_basis_and_basis_function_indices(qs_env, bs_env)
350 : TYPE(qs_environment_type), POINTER :: qs_env
351 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
352 :
353 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_RI_basis_and_basis_function_indices'
354 :
355 : INTEGER :: handle, i_RI, iatom, ikind, iset, &
356 : max_AO_bf_per_atom, n_ao_test, n_atom, &
357 : n_kind, n_RI, nset, nsgf, u
358 42 : INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of
359 42 : INTEGER, DIMENSION(:), POINTER :: l_max, l_min, nsgf_set
360 42 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
361 : TYPE(gto_basis_set_type), POINTER :: basis_set_a
362 42 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
363 :
364 42 : CALL timeset(routineN, handle)
365 :
366 : ! determine RI basis set size
367 42 : CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, qs_kind_set=qs_kind_set)
368 :
369 42 : n_kind = SIZE(qs_kind_set)
370 42 : n_atom = bs_env%n_atom
371 :
372 42 : CALL get_atomic_kind_set(atomic_kind_set, kind_of=kind_of)
373 :
374 110 : DO ikind = 1, n_kind
375 : CALL get_qs_kind(qs_kind=qs_kind_set(ikind), basis_set=basis_set_a, &
376 68 : basis_type="RI_AUX")
377 110 : IF (.NOT. ASSOCIATED(basis_set_a)) THEN
378 : CALL cp_abort(__LOCATION__, &
379 0 : "At least one RI_AUX basis set was not explicitly invoked in &KIND-section.")
380 : END IF
381 : END DO
382 :
383 126 : ALLOCATE (bs_env%i_RI_start_from_atom(n_atom))
384 84 : ALLOCATE (bs_env%i_RI_end_from_atom(n_atom))
385 84 : ALLOCATE (bs_env%i_ao_start_from_atom(n_atom))
386 84 : ALLOCATE (bs_env%i_ao_end_from_atom(n_atom))
387 :
388 42 : n_RI = 0
389 144 : DO iatom = 1, n_atom
390 102 : bs_env%i_RI_start_from_atom(iatom) = n_RI + 1
391 102 : ikind = kind_of(iatom)
392 102 : CALL get_qs_kind(qs_kind=qs_kind_set(ikind), nsgf=nsgf, basis_type="RI_AUX")
393 102 : n_RI = n_RI + nsgf
394 144 : bs_env%i_RI_end_from_atom(iatom) = n_RI
395 : END DO
396 42 : bs_env%n_RI = n_RI
397 :
398 42 : max_AO_bf_per_atom = 0
399 42 : n_ao_test = 0
400 144 : DO iatom = 1, n_atom
401 102 : bs_env%i_ao_start_from_atom(iatom) = n_ao_test + 1
402 102 : ikind = kind_of(iatom)
403 102 : CALL get_qs_kind(qs_kind=qs_kind_set(ikind), nsgf=nsgf, basis_type="ORB")
404 102 : n_ao_test = n_ao_test + nsgf
405 102 : bs_env%i_ao_end_from_atom(iatom) = n_ao_test
406 144 : max_AO_bf_per_atom = MAX(max_AO_bf_per_atom, nsgf)
407 : END DO
408 42 : CPASSERT(n_ao_test == bs_env%n_ao)
409 42 : bs_env%max_AO_bf_per_atom = max_AO_bf_per_atom
410 :
411 126 : ALLOCATE (bs_env%l_RI(n_RI))
412 42 : i_RI = 0
413 144 : DO iatom = 1, n_atom
414 102 : ikind = kind_of(iatom)
415 :
416 102 : nset = bs_env%basis_set_RI(ikind)%gto_basis_set%nset
417 102 : l_max => bs_env%basis_set_RI(ikind)%gto_basis_set%lmax
418 102 : l_min => bs_env%basis_set_RI(ikind)%gto_basis_set%lmin
419 102 : nsgf_set => bs_env%basis_set_RI(ikind)%gto_basis_set%nsgf_set
420 :
421 456 : DO iset = 1, nset
422 312 : CPASSERT(l_max(iset) == l_min(iset))
423 940 : bs_env%l_RI(i_RI + 1:i_RI + nsgf_set(iset)) = l_max(iset)
424 414 : i_RI = i_RI + nsgf_set(iset)
425 : END DO
426 :
427 : END DO
428 42 : CPASSERT(i_RI == n_RI)
429 :
430 42 : u = bs_env%unit_nr
431 :
432 42 : IF (u > 0) THEN
433 21 : WRITE (u, FMT="(T2,A)") " "
434 21 : WRITE (u, FMT="(T2,2A,T75,I8)") "Number of auxiliary Gaussian basis functions ", &
435 42 : "for χ, ε, W", n_RI
436 : END IF
437 :
438 42 : CALL timestop(handle)
439 :
440 84 : END SUBROUTINE get_RI_basis_and_basis_function_indices
441 :
442 : ! **************************************************************************************************
443 : !> \brief ...
444 : !> \param bs_env ...
445 : !> \param kpoints ...
446 : ! **************************************************************************************************
447 42 : SUBROUTINE setup_kpoints_chi_eps_W(bs_env, kpoints)
448 :
449 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
450 : TYPE(kpoint_type), POINTER :: kpoints
451 :
452 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_kpoints_chi_eps_W'
453 :
454 : INTEGER :: handle, i_dim, n_dim, nkp, nkp_extra, &
455 : nkp_orig, u
456 : INTEGER, DIMENSION(3) :: nkp_grid, nkp_grid_extra, periodic
457 : REAL(KIND=dp) :: exp_s_p, n_dim_inv
458 :
459 42 : CALL timeset(routineN, handle)
460 :
461 : ! routine adapted from mp2_integrals.F
462 42 : NULLIFY (kpoints)
463 42 : CALL kpoint_create(kpoints)
464 :
465 42 : kpoints%kp_scheme = "GENERAL"
466 :
467 168 : periodic(1:3) = bs_env%periodic(1:3)
468 :
469 42 : CPASSERT(SIZE(bs_env%nkp_grid_chi_eps_W_input) == 3)
470 :
471 : IF (bs_env%nkp_grid_chi_eps_W_input(1) > 0 .AND. &
472 42 : bs_env%nkp_grid_chi_eps_W_input(2) > 0 .AND. &
473 : bs_env%nkp_grid_chi_eps_W_input(3) > 0) THEN
474 : ! 1. k-point mesh for χ, ε, W from input
475 0 : DO i_dim = 1, 3
476 0 : SELECT CASE (periodic(i_dim))
477 : CASE (0)
478 0 : nkp_grid(i_dim) = 1
479 0 : nkp_grid_extra(i_dim) = 1
480 : CASE (1)
481 0 : nkp_grid(i_dim) = bs_env%nkp_grid_chi_eps_W_input(i_dim)
482 0 : nkp_grid_extra(i_dim) = nkp_grid(i_dim)*2
483 : CASE DEFAULT
484 0 : CPABORT("Error in periodicity.")
485 : END SELECT
486 : END DO
487 :
488 : ELSE IF (bs_env%nkp_grid_chi_eps_W_input(1) == -1 .AND. &
489 42 : bs_env%nkp_grid_chi_eps_W_input(2) == -1 .AND. &
490 : bs_env%nkp_grid_chi_eps_W_input(3) == -1) THEN
491 : ! 2. automatic k-point mesh for χ, ε, W
492 :
493 168 : DO i_dim = 1, 3
494 :
495 126 : CPASSERT(periodic(i_dim) == 0 .OR. periodic(i_dim) == 1)
496 :
497 42 : SELECT CASE (periodic(i_dim))
498 : CASE (0)
499 90 : nkp_grid(i_dim) = 1
500 90 : nkp_grid_extra(i_dim) = 1
501 : CASE (1)
502 56 : SELECT CASE (bs_env%small_cell_full_kp_or_large_cell_Gamma)
503 : CASE (large_cell_Gamma, large_cell_Gamma_ri_rs)
504 20 : nkp_grid(i_dim) = 4
505 20 : nkp_grid_extra(i_dim) = 6
506 : CASE (small_cell_full_kp)
507 16 : nkp_grid(i_dim) = bs_env%kpoints_scf_desymm%nkp_grid(i_dim)*4
508 36 : nkp_grid_extra(i_dim) = bs_env%kpoints_scf_desymm%nkp_grid(i_dim)*8
509 : END SELECT
510 : CASE DEFAULT
511 126 : CPABORT("Error in periodicity.")
512 : END SELECT
513 :
514 : END DO
515 :
516 : ELSE
517 :
518 0 : CPABORT("An error occured when setting up the k-mesh for W.")
519 :
520 : END IF
521 :
522 42 : nkp_orig = MAX(nkp_grid(1)*nkp_grid(2)*nkp_grid(3)/2, 1)
523 :
524 42 : nkp_extra = nkp_grid_extra(1)*nkp_grid_extra(2)*nkp_grid_extra(3)/2
525 :
526 42 : nkp = nkp_orig + nkp_extra
527 :
528 168 : kpoints%nkp_grid(1:3) = nkp_grid(1:3)
529 42 : kpoints%nkp = nkp
530 :
531 168 : bs_env%nkp_grid_chi_eps_W_orig(1:3) = nkp_grid(1:3)
532 168 : bs_env%nkp_grid_chi_eps_W_extra(1:3) = nkp_grid_extra(1:3)
533 42 : bs_env%nkp_chi_eps_W_orig = nkp_orig
534 42 : bs_env%nkp_chi_eps_W_extra = nkp_extra
535 42 : bs_env%nkp_chi_eps_W_orig_plus_extra = nkp
536 :
537 210 : ALLOCATE (kpoints%xkp(3, nkp), kpoints%wkp(nkp))
538 126 : ALLOCATE (bs_env%wkp_no_extra(nkp), bs_env%wkp_s_p(nkp))
539 :
540 42 : CALL compute_xkp(kpoints%xkp, 1, nkp_orig, nkp_grid)
541 42 : CALL compute_xkp(kpoints%xkp, nkp_orig + 1, nkp, nkp_grid_extra)
542 :
543 168 : n_dim = SUM(periodic)
544 42 : IF (n_dim == 0) THEN
545 : ! molecules
546 24 : kpoints%wkp(1) = 1.0_dp
547 24 : bs_env%wkp_s_p(1) = 1.0_dp
548 24 : bs_env%wkp_no_extra(1) = 1.0_dp
549 : ELSE
550 :
551 18 : n_dim_inv = 1.0_dp/REAL(n_dim, KIND=dp)
552 :
553 : ! k-point weights are chosen to automatically extrapolate the k-point mesh
554 18 : CALL compute_wkp(kpoints%wkp(1:nkp_orig), nkp_orig, nkp_extra, n_dim_inv)
555 18 : CALL compute_wkp(kpoints%wkp(nkp_orig + 1:nkp), nkp_extra, nkp_orig, n_dim_inv)
556 :
557 1122 : bs_env%wkp_no_extra(1:nkp_orig) = 0.0_dp
558 4294 : bs_env%wkp_no_extra(nkp_orig + 1:nkp) = 1.0_dp/REAL(nkp_extra, KIND=dp)
559 :
560 18 : IF (n_dim == 3) THEN
561 : ! W_PQ(k) for an s-function P and a p-function Q diverges as 1/k at k=0
562 : ! (instead of 1/k^2 for P and Q both being s-functions).
563 0 : exp_s_p = 2.0_dp*n_dim_inv
564 0 : CALL compute_wkp(bs_env%wkp_s_p(1:nkp_orig), nkp_orig, nkp_extra, exp_s_p)
565 0 : CALL compute_wkp(bs_env%wkp_s_p(nkp_orig + 1:nkp), nkp_extra, nkp_orig, exp_s_p)
566 : ELSE
567 5398 : bs_env%wkp_s_p(1:nkp) = bs_env%wkp_no_extra(1:nkp)
568 : END IF
569 :
570 : END IF
571 :
572 42 : IF (bs_env%approx_kp_extrapol) THEN
573 2 : bs_env%wkp_orig = 1.0_dp/REAL(nkp_orig, KIND=dp)
574 : END IF
575 :
576 : ! heuristic parameter: how many k-points for χ, ε, and W are used simultaneously
577 : ! (less simultaneous k-points: less memory, but more computational effort because of
578 : ! recomputation of V(k))
579 42 : bs_env%nkp_chi_eps_W_batch = 4
580 :
581 : bs_env%num_chi_eps_W_batches = (bs_env%nkp_chi_eps_W_orig_plus_extra - 1)/ &
582 42 : bs_env%nkp_chi_eps_W_batch + 1
583 :
584 42 : u = bs_env%unit_nr
585 :
586 42 : IF (u > 0) THEN
587 21 : WRITE (u, FMT="(T2,A)") " "
588 21 : WRITE (u, FMT="(T2,1A,T71,3I4)") "K-point mesh 1 for χ, ε, W", nkp_grid(1:3)
589 21 : WRITE (u, FMT="(T2,2A,T71,3I4)") "K-point mesh 2 for χ, ε, W ", &
590 42 : "(for k-point extrapolation of W)", nkp_grid_extra(1:3)
591 21 : WRITE (u, FMT="(T2,A,T80,L)") "Approximate the k-point extrapolation", &
592 42 : bs_env%approx_kp_extrapol
593 : END IF
594 :
595 42 : CALL timestop(handle)
596 :
597 42 : END SUBROUTINE setup_kpoints_chi_eps_W
598 :
599 : ! **************************************************************************************************
600 : !> \brief ...
601 : !> \param xkp ...
602 : !> \param ikp_start ...
603 : !> \param ikp_end ...
604 : !> \param grid ...
605 : ! **************************************************************************************************
606 84 : SUBROUTINE compute_xkp(xkp, ikp_start, ikp_end, grid)
607 :
608 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: xkp
609 : INTEGER :: ikp_start, ikp_end
610 : INTEGER, DIMENSION(3) :: grid
611 :
612 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_xkp'
613 :
614 : INTEGER :: handle, i, ix, iy, iz
615 :
616 84 : CALL timeset(routineN, handle)
617 :
618 84 : i = ikp_start
619 292 : DO ix = 1, grid(1)
620 3456 : DO iy = 1, grid(2)
621 14180 : DO iz = 1, grid(3)
622 :
623 10808 : IF (i > ikp_end) CYCLE
624 :
625 5404 : xkp(1, i) = REAL(2*ix - grid(1) - 1, KIND=dp)/(2._dp*REAL(grid(1), KIND=dp))
626 5404 : xkp(2, i) = REAL(2*iy - grid(2) - 1, KIND=dp)/(2._dp*REAL(grid(2), KIND=dp))
627 5404 : xkp(3, i) = REAL(2*iz - grid(3) - 1, KIND=dp)/(2._dp*REAL(grid(3), KIND=dp))
628 13972 : i = i + 1
629 :
630 : END DO
631 : END DO
632 : END DO
633 :
634 84 : CALL timestop(handle)
635 :
636 84 : END SUBROUTINE compute_xkp
637 :
638 : ! **************************************************************************************************
639 : !> \brief ...
640 : !> \param wkp ...
641 : !> \param nkp_1 ...
642 : !> \param nkp_2 ...
643 : !> \param exponent ...
644 : ! **************************************************************************************************
645 36 : SUBROUTINE compute_wkp(wkp, nkp_1, nkp_2, exponent)
646 : REAL(KIND=dp), DIMENSION(:) :: wkp
647 : INTEGER :: nkp_1, nkp_2
648 : REAL(KIND=dp) :: exponent
649 :
650 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_wkp'
651 :
652 : INTEGER :: handle
653 : REAL(KIND=dp) :: nkp_ratio
654 :
655 36 : CALL timeset(routineN, handle)
656 :
657 36 : nkp_ratio = REAL(nkp_2, KIND=dp)/REAL(nkp_1, KIND=dp)
658 :
659 5416 : wkp(:) = 1.0_dp/REAL(nkp_1, KIND=dp)/(1.0_dp - nkp_ratio**exponent)
660 :
661 36 : CALL timestop(handle)
662 :
663 36 : END SUBROUTINE compute_wkp
664 :
665 : ! **************************************************************************************************
666 : !> \brief ...
667 : !> \param qs_env ...
668 : !> \param bs_env ...
669 : ! **************************************************************************************************
670 42 : SUBROUTINE allocate_matrices(qs_env, bs_env)
671 : TYPE(qs_environment_type), POINTER :: qs_env
672 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
673 :
674 : CHARACTER(LEN=*), PARAMETER :: routineN = 'allocate_matrices'
675 :
676 : INTEGER :: handle, i_t
677 : TYPE(cp_blacs_env_type), POINTER :: blacs_env, blacs_env_tensor
678 : TYPE(cp_fm_struct_type), POINTER :: fm_struct, fm_struct_RI_global
679 : TYPE(mp_para_env_type), POINTER :: para_env
680 :
681 42 : CALL timeset(routineN, handle)
682 :
683 42 : CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
684 :
685 42 : fm_struct => bs_env%fm_ks_Gamma(1)%matrix_struct
686 :
687 42 : CALL cp_fm_create(bs_env%fm_Gocc, fm_struct)
688 42 : CALL cp_fm_create(bs_env%fm_Gvir, fm_struct)
689 :
690 42 : NULLIFY (fm_struct_RI_global)
691 : CALL cp_fm_struct_create(fm_struct_RI_global, context=blacs_env, nrow_global=bs_env%n_RI, &
692 42 : ncol_global=bs_env%n_RI, para_env=para_env)
693 42 : CALL cp_fm_create(bs_env%fm_RI_RI, fm_struct_RI_global)
694 42 : CALL cp_fm_create(bs_env%fm_chi_Gamma_freq, fm_struct_RI_global)
695 42 : CALL cp_fm_create(bs_env%fm_W_MIC_freq, fm_struct_RI_global)
696 42 : IF (bs_env%approx_kp_extrapol) THEN
697 2 : CALL cp_fm_create(bs_env%fm_W_MIC_freq_1_extra, fm_struct_RI_global)
698 2 : CALL cp_fm_create(bs_env%fm_W_MIC_freq_1_no_extra, fm_struct_RI_global)
699 2 : CALL cp_fm_set_all(bs_env%fm_W_MIC_freq_1_extra, 0.0_dp)
700 2 : CALL cp_fm_set_all(bs_env%fm_W_MIC_freq_1_no_extra, 0.0_dp)
701 : END IF
702 42 : CALL cp_fm_struct_release(fm_struct_RI_global)
703 :
704 : ! create blacs_env for subgroups of tensor operations
705 42 : NULLIFY (blacs_env_tensor)
706 42 : CALL cp_blacs_env_create(blacs_env=blacs_env_tensor, para_env=bs_env%para_env_tensor)
707 :
708 : ! allocate dbcsr matrices in the tensor subgroup; actually, one only needs a small
709 : ! subset of blocks in the tensor subgroup, however, all atomic blocks are allocated.
710 : ! One might think of creating a dbcsr matrix with only the blocks that are needed
711 : ! in the tensor subgroup
712 : CALL create_mat_munu(bs_env%mat_ao_ao_tensor, qs_env, bs_env%eps_atom_grid_2d_mat, &
713 42 : blacs_env_tensor, do_ri_aux_basis=.FALSE.)
714 :
715 : CALL create_mat_munu(bs_env%mat_RI_RI_tensor, qs_env, bs_env%eps_atom_grid_2d_mat, &
716 42 : blacs_env_tensor, do_ri_aux_basis=.TRUE.)
717 :
718 : CALL create_mat_munu(bs_env%mat_RI_RI, qs_env, bs_env%eps_atom_grid_2d_mat, &
719 42 : blacs_env, do_ri_aux_basis=.TRUE.)
720 :
721 42 : CALL cp_blacs_env_release(blacs_env_tensor)
722 :
723 42 : NULLIFY (bs_env%mat_chi_Gamma_tau)
724 42 : CALL dbcsr_allocate_matrix_set(bs_env%mat_chi_Gamma_tau, bs_env%num_time_freq_points)
725 :
726 572 : DO i_t = 1, bs_env%num_time_freq_points
727 530 : ALLOCATE (bs_env%mat_chi_Gamma_tau(i_t)%matrix)
728 572 : CALL dbcsr_create(bs_env%mat_chi_Gamma_tau(i_t)%matrix, template=bs_env%mat_RI_RI%matrix)
729 : END DO
730 :
731 42 : CALL timestop(handle)
732 :
733 42 : END SUBROUTINE allocate_matrices
734 :
735 : ! **************************************************************************************************
736 : !> \brief ...
737 : !> \param bs_env ...
738 : ! **************************************************************************************************
739 34 : SUBROUTINE allocate_GW_eigenvalues(bs_env)
740 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
741 :
742 : CHARACTER(LEN=*), PARAMETER :: routineN = 'allocate_GW_eigenvalues'
743 :
744 : INTEGER :: handle
745 :
746 34 : CALL timeset(routineN, handle)
747 :
748 170 : ALLOCATE (bs_env%eigenval_G0W0(bs_env%n_ao, bs_env%nkp_bs_and_DOS, bs_env%n_spin))
749 170 : ALLOCATE (bs_env%eigenval_HF(bs_env%n_ao, bs_env%nkp_bs_and_DOS, bs_env%n_spin))
750 :
751 34 : CALL timestop(handle)
752 :
753 34 : END SUBROUTINE allocate_GW_eigenvalues
754 :
755 : ! **************************************************************************************************
756 : !> \brief ...
757 : !> \param qs_env ...
758 : !> \param bs_env ...
759 : ! **************************************************************************************************
760 42 : SUBROUTINE create_tensors(qs_env, bs_env)
761 : TYPE(qs_environment_type), POINTER :: qs_env
762 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
763 :
764 : CHARACTER(LEN=*), PARAMETER :: routineN = 'create_tensors'
765 :
766 : INTEGER :: handle
767 :
768 42 : CALL timeset(routineN, handle)
769 :
770 42 : CALL init_interaction_radii(bs_env)
771 :
772 : ! split blocks does not improve load balancing/efficienfy for tensor contraction, so we go
773 : ! with the standard atomic blocks
774 : CALL create_3c_t(bs_env%t_RI_AO__AO, bs_env%para_env_tensor, "(RI AO | AO)", [1, 2], [3], &
775 : bs_env%sizes_RI, bs_env%sizes_AO, &
776 42 : create_nl_3c=.TRUE., nl_3c=bs_env%nl_3c, qs_env=qs_env)
777 : CALL create_3c_t(bs_env%t_RI__AO_AO, bs_env%para_env_tensor, "(RI | AO AO)", [1], [2, 3], &
778 42 : bs_env%sizes_RI, bs_env%sizes_AO)
779 :
780 42 : CALL create_2c_t(bs_env, bs_env%sizes_RI, bs_env%sizes_AO)
781 :
782 42 : CALL timestop(handle)
783 :
784 42 : END SUBROUTINE create_tensors
785 :
786 : ! **************************************************************************************************
787 : !> \brief ...
788 : !> \param qs_env ...
789 : !> \param bs_env ...
790 : ! **************************************************************************************************
791 34 : SUBROUTINE check_sparsity_3c(qs_env, bs_env)
792 : TYPE(qs_environment_type), POINTER :: qs_env
793 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
794 :
795 : CHARACTER(LEN=*), PARAMETER :: routineN = 'check_sparsity_3c'
796 :
797 : INTEGER :: handle, n_atom_step, RI_atom
798 : INTEGER(int_8) :: non_zero_elements_sum, nze
799 : REAL(dp) :: max_dist_AO_atoms, occ, occupation_sum
800 : REAL(KIND=dp) :: t1, t2
801 34 : TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :) :: t_3c_global_array
802 :
803 : !TYPE(dbt_type) :: t_3c_global
804 :
805 : !TYPE(neighbor_list_3c_type) :: nl_3c_global
806 :
807 34 : CALL timeset(routineN, handle)
808 :
809 : ! check the sparsity of 3c integral tensor (µν|P); calculate maximum distance between
810 : ! AO atoms µ, ν where at least a single integral (µν|P) is larger than the filter threshold
811 :
812 306 : ALLOCATE (t_3c_global_array(1, 1))
813 34 : CALL dbt_create(bs_env%t_RI_AO__AO, t_3c_global_array(1, 1))
814 :
815 : ! Allocate arrays to store min/max indices for overlap with other AO/RI functions on each atom
816 : ! (Filled during loop via get_i_j_atom_ranges)
817 136 : ALLOCATE (bs_env%min_RI_idx_from_AO_AO_atom(bs_env%n_atom, bs_env%n_atom))
818 136 : ALLOCATE (bs_env%max_RI_idx_from_AO_AO_atom(bs_env%n_atom, bs_env%n_atom))
819 136 : ALLOCATE (bs_env%min_AO_idx_from_RI_AO_atom(bs_env%n_atom, bs_env%n_atom))
820 136 : ALLOCATE (bs_env%max_AO_idx_from_RI_AO_atom(bs_env%n_atom, bs_env%n_atom))
821 310 : bs_env%min_RI_idx_from_AO_AO_atom(:, :) = bs_env%n_RI
822 310 : bs_env%max_RI_idx_from_AO_AO_atom(:, :) = 1
823 310 : bs_env%min_AO_idx_from_RI_AO_atom(:, :) = bs_env%n_AO
824 310 : bs_env%max_AO_idx_from_RI_AO_atom(:, :) = 1
825 :
826 34 : CALL bs_env%para_env%sync()
827 34 : t1 = m_walltime()
828 :
829 34 : occupation_sum = 0.0_dp
830 34 : non_zero_elements_sum = 0
831 34 : max_dist_AO_atoms = 0.0_dp
832 34 : n_atom_step = INT(SQRT(REAL(bs_env%n_atom, KIND=dp)))
833 : ! do not compute full 3c integrals at once because it may cause out of memory
834 114 : DO RI_atom = 1, bs_env%n_atom, n_atom_step
835 :
836 : CALL build_3c_integrals(t_3c_global_array, &
837 : bs_env%eps_filter, &
838 : qs_env, &
839 : bs_env%nl_3c, &
840 : int_eps=bs_env%eps_filter, &
841 : basis_i=bs_env%basis_set_RI, &
842 : basis_j=bs_env%basis_set_AO, &
843 : basis_k=bs_env%basis_set_AO, &
844 : bounds_i=[RI_atom, MIN(RI_atom + n_atom_step - 1, bs_env%n_atom)], &
845 : potential_parameter=bs_env%ri_metric, &
846 240 : desymmetrize=.FALSE.)
847 :
848 80 : CALL dbt_filter(t_3c_global_array(1, 1), bs_env%eps_filter)
849 :
850 80 : CALL bs_env%para_env%sync()
851 :
852 80 : CALL get_tensor_occupancy(t_3c_global_array(1, 1), nze, occ)
853 80 : non_zero_elements_sum = non_zero_elements_sum + nze
854 80 : occupation_sum = occupation_sum + occ
855 :
856 80 : CALL get_max_dist_AO_atoms(t_3c_global_array(1, 1), max_dist_AO_atoms, qs_env)
857 :
858 : ! Extract indices per block
859 80 : CALL get_i_j_atom_ranges(t_3c_global_array(1, 1), bs_env)
860 :
861 194 : CALL dbt_clear(t_3c_global_array(1, 1))
862 :
863 : END DO
864 :
865 34 : t2 = m_walltime()
866 :
867 : ! Sync/max for max_dist_AO_atoms is done inside each get_max_dist_AO_atoms
868 34 : bs_env%max_dist_AO_atoms = max_dist_AO_atoms
869 : ! occupation_sum is a global quantity, also needs no sync here
870 34 : bs_env%occupation_3c_int = occupation_sum
871 :
872 34 : CALL bs_env%para_env%min(bs_env%min_RI_idx_from_AO_AO_atom)
873 34 : CALL bs_env%para_env%max(bs_env%max_RI_idx_from_AO_AO_atom)
874 34 : CALL bs_env%para_env%min(bs_env%min_AO_idx_from_RI_AO_atom)
875 34 : CALL bs_env%para_env%max(bs_env%max_AO_idx_from_RI_AO_atom)
876 :
877 34 : CALL dbt_destroy(t_3c_global_array(1, 1))
878 68 : DEALLOCATE (t_3c_global_array)
879 :
880 34 : IF (bs_env%unit_nr > 0) THEN
881 17 : WRITE (bs_env%unit_nr, '(T2,A)') ''
882 : WRITE (bs_env%unit_nr, '(T2,A,F27.1,A)') &
883 17 : 'Computed 3-center integrals (µν|P), execution time', t2 - t1, ' s'
884 17 : WRITE (bs_env%unit_nr, '(T2,A,F48.3,A)') 'Percentage of non-zero (µν|P)', &
885 34 : bs_env%occupation_3c_int*100, ' %'
886 17 : WRITE (bs_env%unit_nr, '(T2,A,F33.1,A)') 'Max. distance between µ,ν in non-zero (µν|P)', &
887 34 : bs_env%max_dist_AO_atoms*angstrom, ' A'
888 17 : WRITE (bs_env%unit_nr, '(T2,2A,I20,A)') 'Required memory if storing all 3-center ', &
889 34 : 'integrals (µν|P)', INT(REAL(non_zero_elements_sum, KIND=dp)*8.0E-9_dp), ' GB'
890 : END IF
891 :
892 34 : CALL timestop(handle)
893 :
894 68 : END SUBROUTINE check_sparsity_3c
895 :
896 : ! **************************************************************************************************
897 : !> \brief ...
898 : !> \param t_3c ...
899 : !> \param bs_env ...
900 : ! **************************************************************************************************
901 80 : SUBROUTINE get_i_j_atom_ranges(t_3c, bs_env)
902 : TYPE(dbt_type) :: t_3c
903 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
904 :
905 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_i_j_atom_ranges'
906 :
907 : INTEGER :: handle, idx_AO_end, idx_AO_start, &
908 : idx_RI_end, idx_RI_start
909 : INTEGER, DIMENSION(3) :: atom_ind
910 : TYPE(dbt_iterator_type) :: iter
911 :
912 80 : CALL timeset(routineN, handle)
913 :
914 : ! Loop over blocks in 3c, for given min_atom: RI_min/max index from min_atom
915 : !$OMP PARALLEL DEFAULT(NONE) &
916 : !$OMP SHARED(t_3c, bs_env) &
917 : !$OMP PRIVATE(iter, atom_ind, &
918 80 : !$OMP idx_RI_start, idx_RI_end, idx_AO_start, idx_AO_end)
919 :
920 : CALL dbt_iterator_start(iter, t_3c)
921 : DO WHILE (dbt_iterator_blocks_left(iter))
922 : CALL dbt_iterator_next_block(iter, atom_ind)
923 :
924 : ! Pre-fetch indices to avoid referencing 'bs_env' twice inside the ATOMIC blocks
925 : idx_RI_start = bs_env%i_RI_start_from_atom(atom_ind(1))
926 : idx_RI_end = bs_env%i_RI_end_from_atom(atom_ind(1))
927 :
928 : idx_AO_start = bs_env%i_ao_start_from_atom(atom_ind(2))
929 : idx_AO_end = bs_env%i_ao_end_from_atom(atom_ind(2))
930 :
931 : ! Update values safely inside ATOMIC blocks, otherwise race conditions occur
932 : !$OMP ATOMIC UPDATE
933 : bs_env%min_RI_idx_from_AO_AO_atom(atom_ind(2), atom_ind(3)) = &
934 : MIN(bs_env%min_RI_idx_from_AO_AO_atom(atom_ind(2), atom_ind(3)), idx_RI_start)
935 : !$OMP ATOMIC UPDATE
936 : bs_env%max_RI_idx_from_AO_AO_atom(atom_ind(2), atom_ind(3)) = &
937 : MAX(bs_env%max_RI_idx_from_AO_AO_atom(atom_ind(2), atom_ind(3)), idx_RI_end)
938 :
939 : !$OMP ATOMIC UPDATE
940 : bs_env%min_AO_idx_from_RI_AO_atom(atom_ind(1), atom_ind(3)) = &
941 : MIN(bs_env%min_AO_idx_from_RI_AO_atom(atom_ind(1), atom_ind(3)), idx_AO_start)
942 : !$OMP ATOMIC UPDATE
943 : bs_env%max_AO_idx_from_RI_AO_atom(atom_ind(1), atom_ind(3)) = &
944 : MAX(bs_env%max_AO_idx_from_RI_AO_atom(atom_ind(1), atom_ind(3)), idx_AO_end)
945 :
946 : END DO
947 : CALL dbt_iterator_stop(iter)
948 : !$OMP END PARALLEL
949 :
950 80 : CALL timestop(handle)
951 :
952 80 : END SUBROUTINE get_i_j_atom_ranges
953 :
954 : ! **************************************************************************************************
955 : !> \brief ...
956 : !> \param bs_env ...
957 : !> \param sizes_RI ...
958 : !> \param sizes_AO ...
959 : ! **************************************************************************************************
960 42 : SUBROUTINE create_2c_t(bs_env, sizes_RI, sizes_AO)
961 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
962 : INTEGER, ALLOCATABLE, DIMENSION(:) :: sizes_RI, sizes_AO
963 :
964 : CHARACTER(LEN=*), PARAMETER :: routineN = 'create_2c_t'
965 :
966 : INTEGER :: handle
967 42 : INTEGER, ALLOCATABLE, DIMENSION(:) :: dist_1, dist_2
968 : INTEGER, DIMENSION(2) :: pdims_2d
969 126 : TYPE(dbt_pgrid_type) :: pgrid_2d
970 :
971 42 : CALL timeset(routineN, handle)
972 :
973 : ! inspired from rpa_im_time.F / hfx_types.F
974 :
975 42 : pdims_2d = 0
976 42 : CALL dbt_pgrid_create(bs_env%para_env_tensor, pdims_2d, pgrid_2d)
977 :
978 : CALL create_2c_tensor(bs_env%t_G, dist_1, dist_2, pgrid_2d, sizes_AO, sizes_AO, &
979 42 : name="(AO | AO)")
980 42 : DEALLOCATE (dist_1, dist_2)
981 : CALL create_2c_tensor(bs_env%t_chi, dist_1, dist_2, pgrid_2d, sizes_RI, sizes_RI, &
982 42 : name="(RI | RI)")
983 42 : DEALLOCATE (dist_1, dist_2)
984 : CALL create_2c_tensor(bs_env%t_W, dist_1, dist_2, pgrid_2d, sizes_RI, sizes_RI, &
985 42 : name="(RI | RI)")
986 42 : DEALLOCATE (dist_1, dist_2)
987 42 : CALL dbt_pgrid_destroy(pgrid_2d)
988 :
989 42 : CALL timestop(handle)
990 :
991 42 : END SUBROUTINE create_2c_t
992 :
993 : ! **************************************************************************************************
994 : !> \brief ...
995 : !> \param tensor ...
996 : !> \param para_env ...
997 : !> \param tensor_name ...
998 : !> \param map1 ...
999 : !> \param map2 ...
1000 : !> \param sizes_RI ...
1001 : !> \param sizes_AO ...
1002 : !> \param create_nl_3c ...
1003 : !> \param nl_3c ...
1004 : !> \param qs_env ...
1005 : ! **************************************************************************************************
1006 84 : SUBROUTINE create_3c_t(tensor, para_env, tensor_name, map1, map2, sizes_RI, sizes_AO, &
1007 : create_nl_3c, nl_3c, qs_env)
1008 : TYPE(dbt_type) :: tensor
1009 : TYPE(mp_para_env_type), POINTER :: para_env
1010 : CHARACTER(LEN=12) :: tensor_name
1011 : INTEGER, DIMENSION(:) :: map1, map2
1012 : INTEGER, ALLOCATABLE, DIMENSION(:) :: sizes_RI, sizes_AO
1013 : LOGICAL, OPTIONAL :: create_nl_3c
1014 : TYPE(neighbor_list_3c_type), OPTIONAL :: nl_3c
1015 : TYPE(qs_environment_type), OPTIONAL, POINTER :: qs_env
1016 :
1017 : CHARACTER(LEN=*), PARAMETER :: routineN = 'create_3c_t'
1018 :
1019 : INTEGER :: handle, nkind
1020 84 : INTEGER, ALLOCATABLE, DIMENSION(:) :: dist_AO_1, dist_AO_2, dist_RI
1021 : INTEGER, DIMENSION(3) :: pcoord, pdims, pdims_3d
1022 : LOGICAL :: my_create_nl_3c
1023 252 : TYPE(dbt_pgrid_type) :: pgrid_3d
1024 : TYPE(distribution_3d_type) :: dist_3d
1025 84 : TYPE(mp_cart_type) :: mp_comm_t3c_2
1026 84 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1027 :
1028 84 : CALL timeset(routineN, handle)
1029 :
1030 84 : pdims_3d = 0
1031 84 : CALL dbt_pgrid_create(para_env, pdims_3d, pgrid_3d)
1032 : CALL create_3c_tensor(tensor, dist_RI, dist_AO_1, dist_AO_2, &
1033 : pgrid_3d, sizes_RI, sizes_AO, sizes_AO, &
1034 84 : map1=map1, map2=map2, name=tensor_name)
1035 :
1036 84 : IF (PRESENT(create_nl_3c)) THEN
1037 42 : my_create_nl_3c = create_nl_3c
1038 : ELSE
1039 : my_create_nl_3c = .FALSE.
1040 : END IF
1041 :
1042 42 : IF (my_create_nl_3c) THEN
1043 42 : CALL get_qs_env(qs_env, nkind=nkind, particle_set=particle_set)
1044 42 : CALL dbt_mp_environ_pgrid(pgrid_3d, pdims, pcoord)
1045 42 : CALL mp_comm_t3c_2%create(pgrid_3d%mp_comm_2d, 3, pdims)
1046 : CALL distribution_3d_create(dist_3d, dist_RI, dist_AO_1, dist_AO_2, &
1047 42 : nkind, particle_set, mp_comm_t3c_2, own_comm=.TRUE.)
1048 :
1049 : CALL build_3c_neighbor_lists(nl_3c, &
1050 : qs_env%bs_env%basis_set_RI, &
1051 : qs_env%bs_env%basis_set_AO, &
1052 : qs_env%bs_env%basis_set_AO, &
1053 : dist_3d, qs_env%bs_env%ri_metric, &
1054 42 : "GW_3c_nl", qs_env, own_dist=.TRUE.)
1055 : END IF
1056 :
1057 84 : DEALLOCATE (dist_RI, dist_AO_1, dist_AO_2)
1058 84 : CALL dbt_pgrid_destroy(pgrid_3d)
1059 :
1060 84 : CALL timestop(handle)
1061 :
1062 168 : END SUBROUTINE create_3c_t
1063 :
1064 : ! **************************************************************************************************
1065 : !> \brief ...
1066 : !> \param bs_env ...
1067 : ! **************************************************************************************************
1068 42 : SUBROUTINE init_interaction_radii(bs_env)
1069 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1070 :
1071 : CHARACTER(LEN=*), PARAMETER :: routineN = 'init_interaction_radii'
1072 :
1073 : INTEGER :: handle, ibasis
1074 : TYPE(gto_basis_set_type), POINTER :: orb_basis, ri_basis
1075 :
1076 42 : CALL timeset(routineN, handle)
1077 :
1078 110 : DO ibasis = 1, SIZE(bs_env%basis_set_AO)
1079 :
1080 68 : orb_basis => bs_env%basis_set_AO(ibasis)%gto_basis_set
1081 68 : CALL init_interaction_radii_orb_basis(orb_basis, bs_env%eps_filter)
1082 :
1083 68 : ri_basis => bs_env%basis_set_RI(ibasis)%gto_basis_set
1084 110 : CALL init_interaction_radii_orb_basis(ri_basis, bs_env%eps_filter)
1085 :
1086 : END DO
1087 :
1088 42 : CALL timestop(handle)
1089 :
1090 42 : END SUBROUTINE init_interaction_radii
1091 :
1092 : ! **************************************************************************************************
1093 : !> \brief ...
1094 : !> \param t_3c_int ...
1095 : !> \param max_dist_AO_atoms ...
1096 : !> \param qs_env ...
1097 : ! **************************************************************************************************
1098 80 : SUBROUTINE get_max_dist_AO_atoms(t_3c_int, max_dist_AO_atoms, qs_env)
1099 : TYPE(dbt_type) :: t_3c_int
1100 : REAL(KIND=dp), INTENT(INOUT) :: max_dist_AO_atoms
1101 : TYPE(qs_environment_type), POINTER :: qs_env
1102 :
1103 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_max_dist_AO_atoms'
1104 :
1105 : INTEGER :: atom_1, atom_2, handle, num_cells
1106 : INTEGER, DIMENSION(3) :: atom_ind
1107 80 : INTEGER, DIMENSION(:, :), POINTER :: index_to_cell
1108 : REAL(KIND=dp) :: abs_rab
1109 : REAL(KIND=dp), DIMENSION(3) :: rab
1110 : TYPE(cell_type), POINTER :: cell
1111 : TYPE(dbt_iterator_type) :: iter
1112 : TYPE(mp_para_env_type), POINTER :: para_env
1113 80 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1114 :
1115 80 : CALL timeset(routineN, handle)
1116 :
1117 80 : NULLIFY (cell, particle_set, para_env)
1118 80 : CALL get_qs_env(qs_env, cell=cell, particle_set=particle_set, para_env=para_env)
1119 :
1120 : ! max_dist_AO_atoms is compared to earlier steps in the loop with step n_atom_step
1121 : ! do not initialize/overwrite here
1122 :
1123 : ! IMPORTANT: Use thread-local copy for max_dist_AO_atoms via REDUCTION to avoid race conditions
1124 : !$OMP PARALLEL DEFAULT(NONE) &
1125 : !$OMP SHARED(t_3c_int, num_cells, index_to_cell, particle_set, cell) &
1126 : !$OMP PRIVATE(iter, atom_ind, rab, abs_rab, atom_1, atom_2) &
1127 80 : !$OMP REDUCTION(MAX:max_dist_AO_atoms)
1128 :
1129 : CALL dbt_iterator_start(iter, t_3c_int)
1130 : DO WHILE (dbt_iterator_blocks_left(iter))
1131 : CALL dbt_iterator_next_block(iter, atom_ind)
1132 :
1133 : atom_1 = atom_ind(2)
1134 : atom_2 = atom_ind(3)
1135 : rab = pbc(particle_set(atom_1)%r(1:3), particle_set(atom_2)%r(1:3), cell)
1136 : abs_rab = SQRT(rab(1)**2 + rab(2)**2 + rab(3)**2)
1137 :
1138 : ! Reduction takes care of using a thread-local copy
1139 : max_dist_AO_atoms = MAX(max_dist_AO_atoms, abs_rab)
1140 :
1141 : END DO
1142 : CALL dbt_iterator_stop(iter)
1143 : !$OMP END PARALLEL
1144 :
1145 80 : CALL para_env%max(max_dist_AO_atoms)
1146 :
1147 80 : CALL timestop(handle)
1148 :
1149 80 : END SUBROUTINE get_max_dist_AO_atoms
1150 :
1151 : ! **************************************************************************************************
1152 : !> \brief ...
1153 : !> \param bs_env ...
1154 : ! **************************************************************************************************
1155 34 : SUBROUTINE set_sparsity_parallelization_parameters(bs_env)
1156 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1157 :
1158 : CHARACTER(LEN=*), PARAMETER :: routineN = 'set_sparsity_parallelization_parameters'
1159 :
1160 : INTEGER :: handle, i_ivl, IL_ivl, j_ivl, n_atom_per_IL_ivl, n_atom_per_ivl, n_intervals_i, &
1161 : n_intervals_inner_loop_atoms, n_intervals_j, u
1162 : INTEGER(KIND=int_8) :: input_memory_per_proc
1163 :
1164 34 : CALL timeset(routineN, handle)
1165 :
1166 : ! heuristic parameter to prevent out of memory
1167 34 : bs_env%safety_factor_memory = 0.10_dp
1168 :
1169 34 : input_memory_per_proc = INT(bs_env%input_memory_per_proc_GB*1.0E9_dp, KIND=int_8)
1170 :
1171 : ! choose atomic range for λ ("i_atom"), ν ("j_atom") in
1172 : ! M_λνP(iτ) = sum_µ (µν|P) G^occ_µλ(i|τ|,k=0)
1173 : ! N_νλQ(iτ) = sum_σ (σλ|Q) G^vir_σν(i|τ|,k=0)
1174 : ! such that M and N fit into the memory
1175 : n_atom_per_ivl = INT(SQRT(bs_env%safety_factor_memory*input_memory_per_proc &
1176 : *bs_env%group_size_tensor/24/bs_env%n_RI &
1177 34 : /SQRT(bs_env%occupation_3c_int)))/bs_env%max_AO_bf_per_atom
1178 :
1179 34 : n_intervals_i = (bs_env%n_atom_i - 1)/n_atom_per_ivl + 1
1180 34 : n_intervals_j = (bs_env%n_atom_j - 1)/n_atom_per_ivl + 1
1181 :
1182 34 : bs_env%n_atom_per_interval_ij = n_atom_per_ivl
1183 34 : bs_env%n_intervals_i = n_intervals_i
1184 34 : bs_env%n_intervals_j = n_intervals_j
1185 :
1186 102 : ALLOCATE (bs_env%i_atom_intervals(2, n_intervals_i))
1187 102 : ALLOCATE (bs_env%j_atom_intervals(2, n_intervals_j))
1188 :
1189 68 : DO i_ivl = 1, n_intervals_i
1190 34 : bs_env%i_atom_intervals(1, i_ivl) = (i_ivl - 1)*n_atom_per_ivl + bs_env%atoms_i(1)
1191 : bs_env%i_atom_intervals(2, i_ivl) = MIN(i_ivl*n_atom_per_ivl + bs_env%atoms_i(1) - 1, &
1192 68 : bs_env%atoms_i(2))
1193 : END DO
1194 :
1195 68 : DO j_ivl = 1, n_intervals_j
1196 34 : bs_env%j_atom_intervals(1, j_ivl) = (j_ivl - 1)*n_atom_per_ivl + bs_env%atoms_j(1)
1197 : bs_env%j_atom_intervals(2, j_ivl) = MIN(j_ivl*n_atom_per_ivl + bs_env%atoms_j(1) - 1, &
1198 68 : bs_env%atoms_j(2))
1199 : END DO
1200 :
1201 136 : ALLOCATE (bs_env%skip_Sigma_occ(n_intervals_i, n_intervals_j))
1202 102 : ALLOCATE (bs_env%skip_Sigma_vir(n_intervals_i, n_intervals_j))
1203 102 : bs_env%skip_Sigma_occ(:, :) = .FALSE.
1204 102 : bs_env%skip_Sigma_vir(:, :) = .FALSE.
1205 34 : bs_env%n_skip_chi = 0
1206 :
1207 102 : ALLOCATE (bs_env%skip_chi(n_intervals_i, n_intervals_j))
1208 102 : bs_env%skip_chi(:, :) = .FALSE.
1209 34 : bs_env%n_skip_sigma = 0
1210 :
1211 : ! choose atomic range for µ and σ ("inner loop (IL) atom") in
1212 : ! M_λνP(iτ) = sum_µ (µν|P) G^occ_µλ(i|τ|,k=0)
1213 : ! N_νλQ(iτ) = sum_σ (σλ|Q) G^vir_σν(i|τ|,k=0)
1214 : n_atom_per_IL_ivl = MIN(INT(bs_env%safety_factor_memory*input_memory_per_proc &
1215 : *bs_env%group_size_tensor/n_atom_per_ivl &
1216 : /bs_env%max_AO_bf_per_atom &
1217 : /bs_env%n_RI/8/SQRT(bs_env%occupation_3c_int) &
1218 34 : /bs_env%max_AO_bf_per_atom), bs_env%n_atom)
1219 :
1220 34 : n_intervals_inner_loop_atoms = (bs_env%n_atom - 1)/n_atom_per_IL_ivl + 1
1221 :
1222 34 : bs_env%n_atom_per_IL_interval = n_atom_per_IL_ivl
1223 34 : bs_env%n_intervals_inner_loop_atoms = n_intervals_inner_loop_atoms
1224 :
1225 102 : ALLOCATE (bs_env%inner_loop_atom_intervals(2, n_intervals_inner_loop_atoms))
1226 68 : DO IL_ivl = 1, n_intervals_inner_loop_atoms
1227 34 : bs_env%inner_loop_atom_intervals(1, IL_ivl) = (IL_ivl - 1)*n_atom_per_IL_ivl + 1
1228 68 : bs_env%inner_loop_atom_intervals(2, IL_ivl) = MIN(IL_ivl*n_atom_per_IL_ivl, bs_env%n_atom)
1229 : END DO
1230 :
1231 34 : u = bs_env%unit_nr
1232 34 : IF (u > 0) THEN
1233 17 : WRITE (u, '(T2,A)') ''
1234 17 : WRITE (u, '(T2,A,I33)') 'Number of i and j atoms in M_λνP(τ), N_νλQ(τ):', n_atom_per_ivl
1235 17 : WRITE (u, '(T2,A,I18)') 'Number of inner loop atoms for µ in M_λνP = sum_µ (µν|P) G_µλ', &
1236 34 : n_atom_per_IL_ivl
1237 : END IF
1238 :
1239 34 : CALL timestop(handle)
1240 :
1241 34 : END SUBROUTINE set_sparsity_parallelization_parameters
1242 :
1243 : ! **************************************************************************************************
1244 : !> \brief ...
1245 : !> \param qs_env ...
1246 : !> \param bs_env ...
1247 : ! **************************************************************************************************
1248 34 : SUBROUTINE check_for_restart_files(qs_env, bs_env)
1249 : TYPE(qs_environment_type), POINTER :: qs_env
1250 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1251 :
1252 : CHARACTER(LEN=*), PARAMETER :: routineN = 'check_for_restart_files'
1253 :
1254 : CHARACTER(LEN=9) :: frmt
1255 : CHARACTER(len=default_path_length) :: f_chi, f_s_n, f_s_p, f_s_x, f_w_t, &
1256 : prefix, project_name, Z_lP_name
1257 : INTEGER :: handle, i_spin, i_t_or_w, ind, n_spin, &
1258 : num_time_freq_points
1259 : LOGICAL :: chi_exists, Sigma_neg_time_exists, &
1260 : Sigma_pos_time_exists, &
1261 : Sigma_x_spin_exists, W_time_exists, &
1262 : Z_lP_exists
1263 : TYPE(cp_logger_type), POINTER :: logger
1264 : TYPE(section_vals_type), POINTER :: input, print_key
1265 :
1266 34 : CALL timeset(routineN, handle)
1267 :
1268 34 : num_time_freq_points = bs_env%num_time_freq_points
1269 34 : n_spin = bs_env%n_spin
1270 :
1271 102 : ALLOCATE (bs_env%read_chi(num_time_freq_points))
1272 68 : ALLOCATE (bs_env%calc_chi(num_time_freq_points))
1273 136 : ALLOCATE (bs_env%Sigma_c_exists(num_time_freq_points, n_spin))
1274 :
1275 34 : CALL get_qs_env(qs_env, input=input)
1276 :
1277 34 : logger => cp_get_default_logger()
1278 34 : print_key => section_vals_get_subs_vals(input, 'PROPERTIES%BANDSTRUCTURE%GW%PRINT%RESTART')
1279 : project_name = cp_print_key_generate_filename(logger, print_key, extension="", &
1280 34 : my_local=.FALSE.)
1281 34 : WRITE (prefix, '(2A)') TRIM(project_name), "-RESTART_"
1282 34 : bs_env%prefix = prefix
1283 :
1284 34 : bs_env%all_W_exist = .TRUE.
1285 :
1286 508 : DO i_t_or_w = 1, num_time_freq_points
1287 :
1288 474 : IF (i_t_or_w < 10) THEN
1289 294 : WRITE (frmt, '(A)') '(3A,I1,A)'
1290 294 : WRITE (f_chi, frmt) TRIM(prefix), bs_env%chi_name, "_0", i_t_or_w, ".matrix"
1291 294 : WRITE (f_W_t, frmt) TRIM(prefix), bs_env%W_time_name, "_0", i_t_or_w, ".matrix"
1292 180 : ELSE IF (i_t_or_w < 100) THEN
1293 180 : WRITE (frmt, '(A)') '(3A,I2,A)'
1294 180 : WRITE (f_chi, frmt) TRIM(prefix), bs_env%chi_name, "_", i_t_or_w, ".matrix"
1295 180 : WRITE (f_W_t, frmt) TRIM(prefix), bs_env%W_time_name, "_", i_t_or_w, ".matrix"
1296 : ELSE
1297 0 : CPABORT('Please implement more than 99 time/frequency points.')
1298 : END IF
1299 :
1300 474 : INQUIRE (file=TRIM(f_chi), exist=chi_exists)
1301 474 : INQUIRE (file=TRIM(f_W_t), exist=W_time_exists)
1302 :
1303 474 : bs_env%read_chi(i_t_or_w) = chi_exists
1304 474 : bs_env%calc_chi(i_t_or_w) = .NOT. chi_exists
1305 :
1306 474 : bs_env%all_W_exist = bs_env%all_W_exist .AND. W_time_exists
1307 :
1308 : ! the self-energy is spin-dependent
1309 1042 : DO i_spin = 1, n_spin
1310 :
1311 534 : ind = i_t_or_w + (i_spin - 1)*num_time_freq_points
1312 :
1313 534 : IF (ind < 10) THEN
1314 294 : WRITE (frmt, '(A)') '(3A,I1,A)'
1315 294 : WRITE (f_S_p, frmt) TRIM(prefix), bs_env%Sigma_p_name, "_0", ind, ".matrix"
1316 294 : WRITE (f_S_n, frmt) TRIM(prefix), bs_env%Sigma_n_name, "_0", ind, ".matrix"
1317 240 : ELSE IF (i_t_or_w < 100) THEN
1318 240 : WRITE (frmt, '(A)') '(3A,I2,A)'
1319 240 : WRITE (f_S_p, frmt) TRIM(prefix), bs_env%Sigma_p_name, "_", ind, ".matrix"
1320 240 : WRITE (f_S_n, frmt) TRIM(prefix), bs_env%Sigma_n_name, "_", ind, ".matrix"
1321 : END IF
1322 :
1323 534 : INQUIRE (file=TRIM(f_S_p), exist=Sigma_pos_time_exists)
1324 534 : INQUIRE (file=TRIM(f_S_n), exist=Sigma_neg_time_exists)
1325 :
1326 : bs_env%Sigma_c_exists(i_t_or_w, i_spin) = Sigma_pos_time_exists .AND. &
1327 1422 : Sigma_neg_time_exists
1328 :
1329 : END DO
1330 :
1331 : END DO
1332 :
1333 : ! Marek : In the RTBSE run, check also for zero frequency W
1334 34 : IF (bs_env%rtp_method == rtp_method_bse) THEN
1335 14 : WRITE (f_W_t, '(3A,I1,A)') TRIM(prefix), "W_freq_rtp", "_0", 0, ".matrix"
1336 14 : INQUIRE (file=TRIM(f_W_t), exist=W_time_exists)
1337 24 : bs_env%all_W_exist = bs_env%all_W_exist .AND. W_time_exists
1338 : END IF
1339 :
1340 : ! Check for Restart Z_lP file
1341 34 : IF (bs_env%do_gw_ri_rs) THEN
1342 10 : WRITE (Z_lP_name, '(3A)') TRIM(prefix), "Z_lP", ".matrix"
1343 10 : INQUIRE (file=TRIM(Z_lP_name), exist=Z_lP_exists)
1344 10 : bs_env%ri_rs%Z_lP_exists = Z_lP_exists
1345 : END IF
1346 :
1347 34 : IF (bs_env%all_W_exist) THEN
1348 106 : bs_env%read_chi(:) = .FALSE.
1349 106 : bs_env%calc_chi(:) = .FALSE.
1350 : END IF
1351 :
1352 34 : bs_env%Sigma_x_exists = .TRUE.
1353 74 : DO i_spin = 1, n_spin
1354 40 : WRITE (f_S_x, '(3A,I1,A)') TRIM(prefix), bs_env%Sigma_x_name, "_0", i_spin, ".matrix"
1355 40 : INQUIRE (file=TRIM(f_S_x), exist=Sigma_x_spin_exists)
1356 106 : bs_env%Sigma_x_exists = bs_env%Sigma_x_exists .AND. Sigma_x_spin_exists
1357 : END DO
1358 :
1359 : ! If any restart files are read, check if the SCF converged in 1 step.
1360 : ! This is important because a re-iterated SCF can lead to spurious GW results
1361 : IF (ANY(bs_env%read_chi(:)) &
1362 : .OR. ANY(bs_env%Sigma_c_exists) &
1363 : .OR. bs_env%all_W_exist &
1364 954 : .OR. bs_env%Sigma_x_exists &
1365 : ) THEN
1366 :
1367 6 : IF (qs_env%scf_env%iter_count /= 1) THEN
1368 : CALL cp_warn(__LOCATION__, "SCF needed more than 1 step, "// &
1369 6 : "which might lead to spurious GW results when using GW restart files. ")
1370 : END IF
1371 : END IF
1372 :
1373 34 : CALL timestop(handle)
1374 :
1375 34 : END SUBROUTINE check_for_restart_files
1376 :
1377 : ! **************************************************************************************************
1378 : !> \brief ...
1379 : !> \param qs_env ...
1380 : !> \param bs_env ...
1381 : ! **************************************************************************************************
1382 42 : SUBROUTINE set_parallelization_parameters(qs_env, bs_env)
1383 : TYPE(qs_environment_type), POINTER :: qs_env
1384 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1385 :
1386 : CHARACTER(LEN=*), PARAMETER :: routineN = 'set_parallelization_parameters'
1387 :
1388 : INTEGER :: color_sub, dummy_1, dummy_2, handle, &
1389 : num_pe, num_t_groups, u
1390 : TYPE(mp_para_env_type), POINTER :: para_env
1391 :
1392 42 : CALL timeset(routineN, handle)
1393 :
1394 42 : CALL get_qs_env(qs_env, para_env=para_env)
1395 :
1396 42 : num_pe = para_env%num_pe
1397 : ! if not already set, use all processors for the group (for large-cell GW, performance
1398 : ! seems to be best for a single group with all MPI processes per group)
1399 42 : IF (bs_env%group_size_tensor < 0 .OR. bs_env%group_size_tensor > num_pe) &
1400 34 : bs_env%group_size_tensor = num_pe
1401 :
1402 : ! group_size_tensor must divide num_pe without rest; otherwise everything will be complicated
1403 42 : IF (MODULO(num_pe, bs_env%group_size_tensor) /= 0) THEN
1404 0 : CALL find_good_group_size(num_pe, bs_env%group_size_tensor)
1405 : END IF
1406 :
1407 : ! para_env_tensor for tensor subgroups
1408 42 : color_sub = para_env%mepos/bs_env%group_size_tensor
1409 42 : bs_env%tensor_group_color = color_sub
1410 :
1411 42 : ALLOCATE (bs_env%para_env_tensor)
1412 42 : CALL bs_env%para_env_tensor%from_split(para_env, color_sub)
1413 :
1414 42 : num_t_groups = para_env%num_pe/bs_env%group_size_tensor
1415 42 : bs_env%num_tensor_groups = num_t_groups
1416 :
1417 : CALL get_i_j_atoms(bs_env%atoms_i, bs_env%atoms_j, bs_env%n_atom_i, bs_env%n_atom_j, &
1418 42 : color_sub, bs_env)
1419 :
1420 126 : ALLOCATE (bs_env%atoms_i_t_group(2, num_t_groups))
1421 84 : ALLOCATE (bs_env%atoms_j_t_group(2, num_t_groups))
1422 92 : DO color_sub = 0, num_t_groups - 1
1423 : CALL get_i_j_atoms(bs_env%atoms_i_t_group(1:2, color_sub + 1), &
1424 : bs_env%atoms_j_t_group(1:2, color_sub + 1), &
1425 92 : dummy_1, dummy_2, color_sub, bs_env)
1426 : END DO
1427 :
1428 42 : u = bs_env%unit_nr
1429 42 : IF (u > 0) THEN
1430 21 : WRITE (u, '(T2,A,I47)') 'Group size for tensor operations', bs_env%group_size_tensor
1431 21 : IF (bs_env%group_size_tensor > 1 .AND. bs_env%n_atom < 5) THEN
1432 17 : WRITE (u, '(T2,A)') 'The requested group size is > 1 which can lead to bad performance.'
1433 17 : WRITE (u, '(T2,A)') 'Using more memory per MPI process might improve performance.'
1434 17 : WRITE (u, '(T2,A)') '(Also increase MEMORY_PER_PROC when using more memory per process.)'
1435 : END IF
1436 : END IF
1437 :
1438 42 : CALL timestop(handle)
1439 :
1440 42 : END SUBROUTINE set_parallelization_parameters
1441 :
1442 : ! **************************************************************************************************
1443 : !> \brief ...
1444 : !> \param num_pe ...
1445 : !> \param group_size ...
1446 : ! **************************************************************************************************
1447 0 : SUBROUTINE find_good_group_size(num_pe, group_size)
1448 :
1449 : INTEGER :: num_pe, group_size
1450 :
1451 : CHARACTER(LEN=*), PARAMETER :: routineN = 'find_good_group_size'
1452 :
1453 : INTEGER :: group_size_minus, group_size_orig, &
1454 : group_size_plus, handle, i_diff
1455 :
1456 0 : CALL timeset(routineN, handle)
1457 :
1458 0 : group_size_orig = group_size
1459 :
1460 0 : DO i_diff = 1, num_pe
1461 :
1462 0 : group_size_minus = group_size - i_diff
1463 :
1464 0 : IF (MODULO(num_pe, group_size_minus) == 0 .AND. group_size_minus > 0) THEN
1465 0 : group_size = group_size_minus
1466 0 : EXIT
1467 : END IF
1468 :
1469 0 : group_size_plus = group_size + i_diff
1470 :
1471 0 : IF (MODULO(num_pe, group_size_plus) == 0 .AND. group_size_plus <= num_pe) THEN
1472 0 : group_size = group_size_plus
1473 0 : EXIT
1474 : END IF
1475 :
1476 : END DO
1477 :
1478 0 : IF (group_size_orig == group_size) CPABORT("Group size error")
1479 :
1480 0 : CALL timestop(handle)
1481 :
1482 0 : END SUBROUTINE find_good_group_size
1483 :
1484 : ! **************************************************************************************************
1485 : !> \brief ...
1486 : !> \param atoms_i ...
1487 : !> \param atoms_j ...
1488 : !> \param n_atom_i ...
1489 : !> \param n_atom_j ...
1490 : !> \param color_sub ...
1491 : !> \param bs_env ...
1492 : ! **************************************************************************************************
1493 92 : SUBROUTINE get_i_j_atoms(atoms_i, atoms_j, n_atom_i, n_atom_j, color_sub, bs_env)
1494 :
1495 : INTEGER, DIMENSION(2) :: atoms_i, atoms_j
1496 : INTEGER :: n_atom_i, n_atom_j, color_sub
1497 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1498 :
1499 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_i_j_atoms'
1500 :
1501 : INTEGER :: handle, i_atoms_per_group, i_group, &
1502 : ipcol, ipcol_loop, iprow, iprow_loop, &
1503 : j_atoms_per_group, npcol, nprow
1504 :
1505 92 : CALL timeset(routineN, handle)
1506 :
1507 : ! create a square mesh of tensor groups for iatom and jatom; code from blacs_env_create
1508 92 : CALL square_mesh(nprow, npcol, bs_env%num_tensor_groups)
1509 :
1510 92 : i_group = 0
1511 184 : DO ipcol_loop = 0, npcol - 1
1512 300 : DO iprow_loop = 0, nprow - 1
1513 116 : IF (i_group == color_sub) THEN
1514 92 : iprow = iprow_loop
1515 92 : ipcol = ipcol_loop
1516 : END IF
1517 208 : i_group = i_group + 1
1518 : END DO
1519 : END DO
1520 :
1521 92 : IF (MODULO(bs_env%n_atom, nprow) == 0) THEN
1522 74 : i_atoms_per_group = bs_env%n_atom/nprow
1523 : ELSE
1524 18 : i_atoms_per_group = bs_env%n_atom/nprow + 1
1525 : END IF
1526 :
1527 92 : IF (MODULO(bs_env%n_atom, npcol) == 0) THEN
1528 92 : j_atoms_per_group = bs_env%n_atom/npcol
1529 : ELSE
1530 0 : j_atoms_per_group = bs_env%n_atom/npcol + 1
1531 : END IF
1532 :
1533 92 : atoms_i(1) = iprow*i_atoms_per_group + 1
1534 92 : atoms_i(2) = MIN((iprow + 1)*i_atoms_per_group, bs_env%n_atom)
1535 92 : n_atom_i = atoms_i(2) - atoms_i(1) + 1
1536 :
1537 92 : atoms_j(1) = ipcol*j_atoms_per_group + 1
1538 92 : atoms_j(2) = MIN((ipcol + 1)*j_atoms_per_group, bs_env%n_atom)
1539 92 : n_atom_j = atoms_j(2) - atoms_j(1) + 1
1540 :
1541 92 : CALL timestop(handle)
1542 :
1543 92 : END SUBROUTINE get_i_j_atoms
1544 :
1545 : ! **************************************************************************************************
1546 : !> \brief ...
1547 : !> \param nprow ...
1548 : !> \param npcol ...
1549 : !> \param nproc ...
1550 : ! **************************************************************************************************
1551 92 : SUBROUTINE square_mesh(nprow, npcol, nproc)
1552 : INTEGER :: nprow, npcol, nproc
1553 :
1554 : CHARACTER(LEN=*), PARAMETER :: routineN = 'square_mesh'
1555 :
1556 : INTEGER :: gcd_max, handle, ipe, jpe
1557 :
1558 92 : CALL timeset(routineN, handle)
1559 :
1560 92 : gcd_max = -1
1561 208 : DO ipe = 1, CEILING(SQRT(REAL(nproc, dp)))
1562 116 : jpe = nproc/ipe
1563 116 : IF (ipe*jpe /= nproc) CYCLE
1564 208 : IF (gcd(ipe, jpe) >= gcd_max) THEN
1565 116 : nprow = ipe
1566 116 : npcol = jpe
1567 116 : gcd_max = gcd(ipe, jpe)
1568 : END IF
1569 : END DO
1570 :
1571 92 : CALL timestop(handle)
1572 :
1573 92 : END SUBROUTINE square_mesh
1574 :
1575 : ! **************************************************************************************************
1576 : !> \brief ...
1577 : !> \param bs_env ...
1578 : !> \param qs_env ...
1579 : ! **************************************************************************************************
1580 42 : SUBROUTINE set_heuristic_parameters(bs_env, qs_env)
1581 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1582 : TYPE(qs_environment_type), OPTIONAL, POINTER :: qs_env
1583 :
1584 : CHARACTER(LEN=*), PARAMETER :: routineN = 'set_heuristic_parameters'
1585 :
1586 : INTEGER :: handle, u
1587 : LOGICAL :: do_BvK_cell
1588 :
1589 42 : CALL timeset(routineN, handle)
1590 :
1591 : ! for generating numerically stable minimax Fourier integration weights
1592 42 : bs_env%num_points_per_magnitude = 200
1593 :
1594 42 : IF (bs_env%input_regularization_minimax > -1.0E-12_dp) THEN
1595 0 : bs_env%regularization_minimax = bs_env%input_regularization_minimax
1596 : ELSE
1597 : ! for periodic systems and for 20 minimax points, we use a regularized minimax mesh
1598 : ! (from experience: regularized minimax meshes converges faster for periodic systems
1599 : ! and for 20 pts)
1600 168 : IF (SUM(bs_env%periodic) /= 0 .OR. bs_env%num_time_freq_points >= 20) THEN
1601 32 : bs_env%regularization_minimax = 1.0E-6_dp
1602 : ELSE
1603 10 : bs_env%regularization_minimax = 0.0_dp
1604 : END IF
1605 : END IF
1606 :
1607 42 : bs_env%stabilize_exp = 70.0_dp
1608 42 : bs_env%eps_atom_grid_2d_mat = 1.0E-50_dp
1609 :
1610 : ! use a 16-parameter Padé fit
1611 42 : bs_env%nparam_pade = 16
1612 :
1613 : ! resolution of the identity with the truncated Coulomb metric, cutoff radius 3 Angström
1614 42 : bs_env%ri_metric%potential_type = do_potential_truncated
1615 42 : bs_env%ri_metric%omega = 0.0_dp
1616 : ! cutoff radius is specified in the input
1617 42 : bs_env%ri_metric%filename = "t_c_g.dat"
1618 :
1619 42 : bs_env%eps_eigval_mat_RI = 0.0_dp
1620 :
1621 42 : IF (bs_env%input_regularization_RI > -1.0E-12_dp) THEN
1622 0 : bs_env%regularization_RI = bs_env%input_regularization_RI
1623 : ELSE
1624 : ! default case:
1625 :
1626 : ! 1. for periodic systems, we use the regularized resolution of the identity per default
1627 42 : bs_env%regularization_RI = 1.0E-2_dp
1628 :
1629 : ! 2. for molecules, no regularization is necessary
1630 168 : IF (SUM(bs_env%periodic) == 0) bs_env%regularization_RI = 0.0_dp
1631 :
1632 : END IF
1633 :
1634 : ! truncated Coulomb operator for exchange self-energy
1635 : ! (see details in Guidon, VandeVondele, Hutter, JCTC 5, 3010 (2009) and references therein)
1636 42 : do_BvK_cell = bs_env%small_cell_full_kp_or_large_cell_Gamma == small_cell_full_kp
1637 : CALL trunc_coulomb_for_exchange(qs_env, bs_env%trunc_coulomb, &
1638 : rel_cutoff_trunc_coulomb_ri_x=0.5_dp, &
1639 : cell_grid=bs_env%cell_grid_scf_desymm, &
1640 42 : do_BvK_cell=do_BvK_cell)
1641 :
1642 : ! for small-cell GW, we need more cells than normally used by the filter bs_env%eps_filter
1643 : ! (in particular for computing the self-energy because of higher number of cells needed)
1644 42 : bs_env%heuristic_filter_factor = 1.0E-4
1645 :
1646 42 : u = bs_env%unit_nr
1647 42 : IF (u > 0) THEN
1648 21 : WRITE (u, FMT="(T2,2A,F21.1,A)") "Cutoff radius for the truncated Coulomb ", &
1649 42 : "operator in Σ^x:", bs_env%trunc_coulomb%cutoff_radius*angstrom, " Å"
1650 21 : WRITE (u, FMT="(T2,2A,F15.1,A)") "Cutoff radius for the truncated Coulomb ", &
1651 42 : "operator in RI metric:", bs_env%ri_metric%cutoff_radius*angstrom, " Å"
1652 21 : WRITE (u, FMT="(T2,A,ES48.1)") "Regularization parameter of RI ", bs_env%regularization_RI
1653 21 : WRITE (u, FMT="(T2,A,ES38.1)") "Regularization parameter of minimax grids", &
1654 42 : bs_env%regularization_minimax
1655 21 : WRITE (u, FMT="(T2,A,I53)") "Lattice sum size for V(k):", bs_env%size_lattice_sum_V
1656 : END IF
1657 :
1658 42 : CALL timestop(handle)
1659 :
1660 42 : END SUBROUTINE set_heuristic_parameters
1661 :
1662 : ! **************************************************************************************************
1663 : !> \brief ...
1664 : !> \param bs_env ...
1665 : ! **************************************************************************************************
1666 42 : SUBROUTINE print_header_and_input_parameters(bs_env)
1667 :
1668 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1669 :
1670 : CHARACTER(LEN=*), PARAMETER :: routineN = 'print_header_and_input_parameters'
1671 :
1672 : INTEGER :: handle, u
1673 :
1674 42 : CALL timeset(routineN, handle)
1675 :
1676 42 : u = bs_env%unit_nr
1677 :
1678 42 : IF (u > 0) THEN
1679 21 : WRITE (u, '(T2,A)') ' '
1680 21 : WRITE (u, '(T2,A)') REPEAT('-', 79)
1681 21 : WRITE (u, '(T2,A,A78)') '-', '-'
1682 21 : WRITE (u, '(T2,A,A46,A32)') '-', 'GW CALCULATION', '-'
1683 21 : WRITE (u, '(T2,A,A78)') '-', '-'
1684 21 : WRITE (u, '(T2,A)') REPEAT('-', 79)
1685 21 : WRITE (u, '(T2,A)') ' '
1686 21 : WRITE (u, '(T2,A,I45)') 'Input: Number of time/freq. points', bs_env%num_time_freq_points
1687 21 : WRITE (u, "(T2,A,F44.1,A)") 'Input: ω_max for fitting Σ(iω) (eV)', bs_env%freq_max_fit*evolt
1688 21 : WRITE (u, '(T2,A,ES27.1)') 'Input: Filter threshold for sparse tensor operations', &
1689 42 : bs_env%eps_filter
1690 21 : WRITE (u, "(T2,A,L55)") 'Input: Apply Hedin shift', bs_env%do_hedin_shift
1691 21 : WRITE (u, '(T2,A,F37.1,A)') 'Input: Available memory per MPI process', &
1692 42 : bs_env%input_memory_per_proc_GB, ' GB'
1693 21 : IF (bs_env%do_gw_ri_rs) THEN
1694 5 : WRITE (u, FMT="(T2,A,ES44.1)") " "
1695 5 : WRITE (u, FMT="(T2,A,ES37.1)") "INPUT: Regularization parameter for RI-RS ", bs_env%ri_rs%tikhonov
1696 5 : WRITE (u, FMT="(T2,A,F31.1,A)") "INPUT: Cutoff radius for grid points in RI-RS ", &
1697 10 : bs_env%ri_rs%cutoff_radius_ri_rs*angstrom, " Å"
1698 5 : WRITE (u, FMT="(T2,A,I42)") "INPUT: Chunk size for DBCSR matrices ", bs_env%ri_rs%chunk_size_dbcsr
1699 5 : WRITE (u, FMT="(T2,A,ES44.1)") " "
1700 : END IF
1701 : END IF
1702 :
1703 42 : CALL timestop(handle)
1704 :
1705 42 : END SUBROUTINE print_header_and_input_parameters
1706 :
1707 : ! **************************************************************************************************
1708 : !> \brief ...
1709 : !> \param qs_env ...
1710 : !> \param bs_env ...
1711 : ! **************************************************************************************************
1712 84 : SUBROUTINE compute_V_xc(qs_env, bs_env)
1713 : TYPE(qs_environment_type), POINTER :: qs_env
1714 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1715 :
1716 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_V_xc'
1717 :
1718 : INTEGER :: handle, img, ispin, myfun, nimages
1719 : LOGICAL :: hf_present
1720 : REAL(KIND=dp) :: energy_ex, energy_exc, energy_total, &
1721 : myfraction
1722 42 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mat_ks_without_v_xc
1723 42 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks_kp
1724 : TYPE(dft_control_type), POINTER :: dft_control
1725 : TYPE(qs_energy_type), POINTER :: energy
1726 : TYPE(section_vals_type), POINTER :: hf_section, input, xc_section
1727 :
1728 42 : CALL timeset(routineN, handle)
1729 :
1730 42 : CALL get_qs_env(qs_env, input=input, energy=energy, dft_control=dft_control)
1731 :
1732 : ! previously, dft_control%nimages set to # neighbor cells, revert for Γ-only KS matrix
1733 42 : nimages = dft_control%nimages
1734 42 : dft_control%nimages = bs_env%nimages_scf
1735 :
1736 : ! we need to reset XC functional, therefore, get XC input
1737 42 : xc_section => section_vals_get_subs_vals(input, "DFT%XC")
1738 42 : CALL section_vals_val_get(xc_section, "XC_FUNCTIONAL%_SECTION_PARAMETERS_", i_val=myfun)
1739 42 : CALL section_vals_val_set(xc_section, "XC_FUNCTIONAL%_SECTION_PARAMETERS_", i_val=xc_none)
1740 : ! IF (ASSOCIATED(section_vals_get_subs_vals(xc_section, "HF", can_return_null=.TRUE.))) THEN
1741 42 : hf_section => section_vals_get_subs_vals(input, "DFT%XC%HF", can_return_null=.TRUE.)
1742 42 : hf_present = .FALSE.
1743 42 : IF (ASSOCIATED(hf_section)) THEN
1744 42 : CALL section_vals_get(hf_section, explicit=hf_present)
1745 : END IF
1746 42 : IF (hf_present) THEN
1747 : ! Special case for handling hfx
1748 4 : CALL section_vals_val_get(xc_section, "HF%FRACTION", r_val=myfraction)
1749 4 : CALL section_vals_val_set(xc_section, "HF%FRACTION", r_val=0.0_dp)
1750 : END IF
1751 :
1752 : ! save the energy before the energy gets updated
1753 42 : energy_total = energy%total
1754 42 : energy_exc = energy%exc
1755 42 : energy_ex = energy%ex
1756 :
1757 76 : SELECT CASE (bs_env%small_cell_full_kp_or_large_cell_Gamma)
1758 : CASE (large_cell_Gamma, large_cell_Gamma_ri_rs)
1759 :
1760 34 : NULLIFY (mat_ks_without_v_xc)
1761 34 : CALL dbcsr_allocate_matrix_set(mat_ks_without_v_xc, bs_env%n_spin)
1762 :
1763 74 : DO ispin = 1, bs_env%n_spin
1764 40 : ALLOCATE (mat_ks_without_v_xc(ispin)%matrix)
1765 74 : IF (hf_present) THEN
1766 : CALL dbcsr_create(mat_ks_without_v_xc(ispin)%matrix, template=bs_env%mat_ao_ao%matrix, &
1767 6 : matrix_type=dbcsr_type_symmetric)
1768 : ELSE
1769 34 : CALL dbcsr_create(mat_ks_without_v_xc(ispin)%matrix, template=bs_env%mat_ao_ao%matrix)
1770 : END IF
1771 : END DO
1772 :
1773 : ! calculate KS-matrix without XC
1774 : CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., just_energy=.FALSE., &
1775 34 : ext_ks_matrix=mat_ks_without_v_xc)
1776 :
1777 74 : DO ispin = 1, bs_env%n_spin
1778 : ! transfer dbcsr matrix to fm
1779 40 : CALL cp_fm_create(bs_env%fm_V_xc_Gamma(ispin), bs_env%fm_s_Gamma%matrix_struct)
1780 40 : CALL copy_dbcsr_to_fm(mat_ks_without_v_xc(ispin)%matrix, bs_env%fm_V_xc_Gamma(ispin))
1781 :
1782 : ! v_xc = h_ks - h_ks(v_xc = 0)
1783 : CALL cp_fm_scale_and_add(alpha=-1.0_dp, matrix_a=bs_env%fm_V_xc_Gamma(ispin), &
1784 74 : beta=1.0_dp, matrix_b=bs_env%fm_ks_Gamma(ispin))
1785 : END DO
1786 :
1787 34 : CALL dbcsr_deallocate_matrix_set(mat_ks_without_v_xc)
1788 :
1789 : CASE (small_cell_full_kp)
1790 :
1791 : ! calculate KS-matrix without XC
1792 8 : CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., just_energy=.FALSE.)
1793 8 : CALL get_qs_env(qs_env=qs_env, matrix_ks_kp=matrix_ks_kp)
1794 :
1795 288 : ALLOCATE (bs_env%fm_V_xc_R(dft_control%nimages, bs_env%n_spin))
1796 58 : DO ispin = 1, bs_env%n_spin
1797 264 : DO img = 1, dft_control%nimages
1798 : ! safe fm_V_xc_R in fm_matrix because saving in dbcsr matrix caused trouble...
1799 248 : CALL copy_dbcsr_to_fm(matrix_ks_kp(ispin, img)%matrix, bs_env%fm_work_mo(1))
1800 : CALL cp_fm_create(bs_env%fm_V_xc_R(img, ispin), bs_env%fm_work_mo(1)%matrix_struct, &
1801 248 : set_zero=.TRUE.)
1802 : ! store h_ks(v_xc = 0) in fm_V_xc_R
1803 : CALL cp_fm_scale_and_add(alpha=1.0_dp, matrix_a=bs_env%fm_V_xc_R(img, ispin), &
1804 256 : beta=1.0_dp, matrix_b=bs_env%fm_work_mo(1))
1805 : END DO
1806 : END DO
1807 :
1808 : END SELECT
1809 :
1810 : ! set back the energy
1811 42 : energy%total = energy_total
1812 42 : energy%exc = energy_exc
1813 42 : energy%ex = energy_ex
1814 :
1815 : ! set back nimages
1816 42 : dft_control%nimages = nimages
1817 :
1818 : ! set the DFT functional and HF fraction back
1819 : CALL section_vals_val_set(xc_section, "XC_FUNCTIONAL%_SECTION_PARAMETERS_", &
1820 42 : i_val=myfun)
1821 42 : IF (hf_present) THEN
1822 : CALL section_vals_val_set(xc_section, "HF%FRACTION", &
1823 4 : r_val=myfraction)
1824 : END IF
1825 :
1826 42 : IF (bs_env%small_cell_full_kp_or_large_cell_Gamma == small_cell_full_kp) THEN
1827 : ! calculate KS-matrix again with XC
1828 8 : CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., just_energy=.FALSE.)
1829 16 : DO ispin = 1, bs_env%n_spin
1830 264 : DO img = 1, dft_control%nimages
1831 : ! store h_ks in fm_work_mo
1832 248 : CALL copy_dbcsr_to_fm(matrix_ks_kp(ispin, img)%matrix, bs_env%fm_work_mo(1))
1833 : ! v_xc = h_ks - h_ks(v_xc = 0)
1834 : CALL cp_fm_scale_and_add(alpha=-1.0_dp, matrix_a=bs_env%fm_V_xc_R(img, ispin), &
1835 256 : beta=1.0_dp, matrix_b=bs_env%fm_work_mo(1))
1836 : END DO
1837 : END DO
1838 : END IF
1839 :
1840 42 : CALL timestop(handle)
1841 :
1842 42 : END SUBROUTINE compute_V_xc
1843 :
1844 : ! **************************************************************************************************
1845 : !> \brief ...
1846 : !> \param bs_env ...
1847 : ! **************************************************************************************************
1848 42 : SUBROUTINE setup_time_and_frequency_minimax_grid(bs_env)
1849 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1850 :
1851 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_time_and_frequency_minimax_grid'
1852 :
1853 : INTEGER :: handle, homo, i_w, ierr, ispin, j_w, &
1854 : n_mo, num_time_freq_points, u
1855 : REAL(KIND=dp) :: E_max, E_max_ispin, E_min, E_min_ispin, &
1856 : E_range, max_error_min
1857 42 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: points_and_weights
1858 :
1859 42 : CALL timeset(routineN, handle)
1860 :
1861 42 : n_mo = bs_env%n_ao
1862 42 : num_time_freq_points = bs_env%num_time_freq_points
1863 :
1864 126 : ALLOCATE (bs_env%imag_freq_points(num_time_freq_points))
1865 84 : ALLOCATE (bs_env%imag_time_points(num_time_freq_points))
1866 84 : ALLOCATE (bs_env%imag_time_weights_freq_zero(num_time_freq_points))
1867 168 : ALLOCATE (bs_env%weights_cos_t_to_w(num_time_freq_points, num_time_freq_points))
1868 126 : ALLOCATE (bs_env%weights_cos_w_to_t(num_time_freq_points, num_time_freq_points))
1869 126 : ALLOCATE (bs_env%weights_sin_t_to_w(num_time_freq_points, num_time_freq_points))
1870 :
1871 : ! minimum and maximum difference between eigenvalues of unoccupied and an occupied MOs
1872 42 : E_min = 1000.0_dp
1873 42 : E_max = -1000.0_dp
1874 90 : DO ispin = 1, bs_env%n_spin
1875 48 : homo = bs_env%n_occ(ispin)
1876 88 : SELECT CASE (bs_env%small_cell_full_kp_or_large_cell_Gamma)
1877 : CASE (large_cell_Gamma, large_cell_Gamma_ri_rs)
1878 : E_min_ispin = bs_env%eigenval_scf_Gamma(homo + 1, ispin) - &
1879 40 : bs_env%eigenval_scf_Gamma(homo, ispin)
1880 : E_max_ispin = bs_env%eigenval_scf_Gamma(n_mo, ispin) - &
1881 40 : bs_env%eigenval_scf_Gamma(1, ispin)
1882 : CASE (small_cell_full_kp)
1883 : E_min_ispin = MINVAL(bs_env%eigenval_scf(homo + 1, :, ispin)) - &
1884 388 : MAXVAL(bs_env%eigenval_scf(homo, :, ispin))
1885 : E_max_ispin = MAXVAL(bs_env%eigenval_scf(n_mo, :, ispin)) - &
1886 436 : MINVAL(bs_env%eigenval_scf(1, :, ispin))
1887 : END SELECT
1888 48 : E_min = MIN(E_min, E_min_ispin)
1889 90 : E_max = MAX(E_max, E_max_ispin)
1890 : END DO
1891 :
1892 42 : E_range = E_max/E_min
1893 :
1894 126 : ALLOCATE (points_and_weights(2*num_time_freq_points))
1895 :
1896 : ! frequency points
1897 42 : IF (num_time_freq_points <= 20) THEN
1898 42 : CALL get_rpa_minimax_coeff(num_time_freq_points, E_range, points_and_weights, ierr, .FALSE.)
1899 : ELSE
1900 0 : CALL get_rpa_minimax_coeff_larger_grid(num_time_freq_points, E_range, points_and_weights)
1901 : END IF
1902 :
1903 : ! one needs to scale the minimax grids, see Azizi, Wilhelm, Golze, Panades-Barrueta,
1904 : ! Giantomassi, Rinke, Draxl, Gonze et al., 2 publications
1905 572 : bs_env%imag_freq_points(:) = points_and_weights(1:num_time_freq_points)*E_min
1906 :
1907 : ! determine number of fit points in the interval [0,ω_max] for virt, or [-ω_max,0] for occ
1908 42 : bs_env%num_freq_points_fit = 0
1909 572 : DO i_w = 1, num_time_freq_points
1910 572 : IF (bs_env%imag_freq_points(i_w) < bs_env%freq_max_fit) THEN
1911 162 : bs_env%num_freq_points_fit = bs_env%num_freq_points_fit + 1
1912 : END IF
1913 : END DO
1914 :
1915 : ! iω values for the analytic continuation Σ^c_n(iω,k) -> Σ^c_n(ϵ,k)
1916 126 : ALLOCATE (bs_env%imag_freq_points_fit(bs_env%num_freq_points_fit))
1917 42 : j_w = 0
1918 572 : DO i_w = 1, num_time_freq_points
1919 572 : IF (bs_env%imag_freq_points(i_w) < bs_env%freq_max_fit) THEN
1920 162 : j_w = j_w + 1
1921 162 : bs_env%imag_freq_points_fit(j_w) = bs_env%imag_freq_points(i_w)
1922 : END IF
1923 : END DO
1924 :
1925 : ! reset the number of Padé parameters if smaller than the number of
1926 : ! imaginary-frequency points for the fit
1927 42 : IF (bs_env%num_freq_points_fit < bs_env%nparam_pade) THEN
1928 42 : bs_env%nparam_pade = bs_env%num_freq_points_fit
1929 : END IF
1930 :
1931 : ! time points
1932 42 : IF (num_time_freq_points <= 20) THEN
1933 42 : CALL get_exp_minimax_coeff(num_time_freq_points, E_range, points_and_weights)
1934 : ELSE
1935 0 : CALL get_exp_minimax_coeff_gw(num_time_freq_points, E_range, points_and_weights)
1936 : END IF
1937 :
1938 572 : bs_env%imag_time_points(:) = points_and_weights(1:num_time_freq_points)/(2.0_dp*E_min)
1939 572 : bs_env%imag_time_weights_freq_zero(:) = points_and_weights(num_time_freq_points + 1:)/(E_min)
1940 :
1941 42 : DEALLOCATE (points_and_weights)
1942 :
1943 42 : u = bs_env%unit_nr
1944 42 : IF (u > 0) THEN
1945 21 : WRITE (u, '(T2,A)') ''
1946 21 : WRITE (u, '(T2,A,F55.2)') 'SCF direct band gap (eV)', E_min*evolt
1947 21 : WRITE (u, '(T2,A,F53.2)') 'Max. SCF eigval diff. (eV)', E_max*evolt
1948 21 : WRITE (u, '(T2,A,F55.2)') 'E-Range for minimax grid', E_range
1949 21 : WRITE (u, '(T2,A,I27)') 'Number of Padé parameters for analytic continuation:', &
1950 42 : bs_env%nparam_pade
1951 21 : WRITE (u, '(T2,A)') ''
1952 : END IF
1953 :
1954 : ! in minimax grids, Fourier transforms t -> w and w -> t are split using
1955 : ! e^(iwt) = cos(wt) + i sin(wt); we thus calculate weights for trafos with a cos and
1956 : ! sine prefactor; details in Azizi, Wilhelm, Golze, Giantomassi, Panades-Barrueta,
1957 : ! Rinke, Draxl, Gonze et al., 2 publications
1958 :
1959 : ! cosine transform weights imaginary time to imaginary frequency
1960 : CALL get_l_sq_wghts_cos_tf_t_to_w(num_time_freq_points, &
1961 : bs_env%imag_time_points, &
1962 : bs_env%weights_cos_t_to_w, &
1963 : bs_env%imag_freq_points, &
1964 : E_min, E_max, max_error_min, &
1965 : bs_env%num_points_per_magnitude, &
1966 42 : bs_env%regularization_minimax)
1967 :
1968 : ! cosine transform weights imaginary frequency to imaginary time
1969 : CALL get_l_sq_wghts_cos_tf_w_to_t(num_time_freq_points, &
1970 : bs_env%imag_time_points, &
1971 : bs_env%weights_cos_w_to_t, &
1972 : bs_env%imag_freq_points, &
1973 : E_min, E_max, max_error_min, &
1974 : bs_env%num_points_per_magnitude, &
1975 42 : bs_env%regularization_minimax)
1976 :
1977 : ! sine transform weights imaginary time to imaginary frequency
1978 : CALL get_l_sq_wghts_sin_tf_t_to_w(num_time_freq_points, &
1979 : bs_env%imag_time_points, &
1980 : bs_env%weights_sin_t_to_w, &
1981 : bs_env%imag_freq_points, &
1982 : E_min, E_max, max_error_min, &
1983 : bs_env%num_points_per_magnitude, &
1984 42 : bs_env%regularization_minimax)
1985 :
1986 42 : CALL timestop(handle)
1987 :
1988 84 : END SUBROUTINE setup_time_and_frequency_minimax_grid
1989 :
1990 : ! **************************************************************************************************
1991 : !> \brief ...
1992 : !> \param qs_env ...
1993 : !> \param bs_env ...
1994 : ! **************************************************************************************************
1995 8 : SUBROUTINE setup_cells_3c(qs_env, bs_env)
1996 :
1997 : TYPE(qs_environment_type), POINTER :: qs_env
1998 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
1999 :
2000 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_cells_3c'
2001 :
2002 : INTEGER :: atom_i, atom_j, atom_k, block_count, handle, i, i_cell_x, i_cell_x_max, &
2003 : i_cell_x_min, i_size, ikind, img, j, j_cell, j_cell_max, j_cell_y, j_cell_y_max, &
2004 : j_cell_y_min, j_size, k_cell, k_cell_max, k_cell_z, k_cell_z_max, k_cell_z_min, k_size, &
2005 : nimage_pairs_3c, nimages_3c, nimages_3c_max, nkind, u
2006 : INTEGER(KIND=int_8) :: mem_occ_per_proc
2007 8 : INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of, n_other_3c_images_max
2008 8 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: index_to_cell_3c_max, nblocks_3c_max
2009 : INTEGER, DIMENSION(3) :: cell_index, n_max
2010 : REAL(KIND=dp) :: avail_mem_per_proc_GB, cell_dist, cell_radius_3c, dij, dik, djk, eps, &
2011 : exp_min_ao, exp_min_RI, frobenius_norm, mem_3c_GB, mem_occ_per_proc_GB, radius_ao, &
2012 : radius_ao_product, radius_RI
2013 8 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: exp_ao_kind, exp_RI_kind, &
2014 8 : radius_ao_kind, &
2015 8 : radius_ao_product_kind, radius_RI_kind
2016 8 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: int_3c
2017 : REAL(KIND=dp), DIMENSION(3) :: rij, rik, rjk, vec_cell_j, vec_cell_k
2018 8 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: exp_ao, exp_RI
2019 8 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
2020 : TYPE(cell_type), POINTER :: cell
2021 8 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
2022 :
2023 8 : CALL timeset(routineN, handle)
2024 :
2025 8 : CALL get_qs_env(qs_env, nkind=nkind, atomic_kind_set=atomic_kind_set, particle_set=particle_set, cell=cell)
2026 :
2027 : ALLOCATE (exp_ao_kind(nkind), exp_RI_kind(nkind), radius_ao_kind(nkind), &
2028 56 : radius_ao_product_kind(nkind), radius_RI_kind(nkind))
2029 :
2030 24 : exp_min_RI = 10.0_dp
2031 24 : exp_min_ao = 10.0_dp
2032 24 : exp_RI_kind = 10.0_dp
2033 24 : exp_AO_kind = 10.0_dp
2034 :
2035 8 : eps = bs_env%eps_filter*bs_env%heuristic_filter_factor
2036 :
2037 24 : DO ikind = 1, nkind
2038 :
2039 16 : CALL get_gto_basis_set(bs_env%basis_set_RI(ikind)%gto_basis_set, zet=exp_RI)
2040 16 : CALL get_gto_basis_set(bs_env%basis_set_ao(ikind)%gto_basis_set, zet=exp_ao)
2041 :
2042 : ! we need to remove all exponents lower than a lower bound, e.g. 1E-3, because
2043 : ! for contracted basis sets, there might be exponents = 0 in zet
2044 32 : DO i = 1, SIZE(exp_RI, 1)
2045 56 : DO j = 1, SIZE(exp_RI, 2)
2046 24 : IF (exp_RI(i, j) < exp_min_RI .AND. exp_RI(i, j) > 1E-3_dp) exp_min_RI = exp_RI(i, j)
2047 24 : IF (exp_RI(i, j) < exp_RI_kind(ikind) .AND. exp_RI(i, j) > 1E-3_dp) &
2048 32 : exp_RI_kind(ikind) = exp_RI(i, j)
2049 : END DO
2050 : END DO
2051 80 : DO i = 1, SIZE(exp_ao, 1)
2052 192 : DO j = 1, SIZE(exp_ao, 2)
2053 112 : IF (exp_ao(i, j) < exp_min_ao .AND. exp_ao(i, j) > 1E-3_dp) exp_min_ao = exp_ao(i, j)
2054 112 : IF (exp_ao(i, j) < exp_ao_kind(ikind) .AND. exp_ao(i, j) > 1E-3_dp) &
2055 112 : exp_ao_kind(ikind) = exp_ao(i, j)
2056 : END DO
2057 : END DO
2058 16 : radius_ao_kind(ikind) = SQRT(-LOG(eps)/exp_ao_kind(ikind))
2059 16 : radius_ao_product_kind(ikind) = SQRT(-LOG(eps)/(2.0_dp*exp_ao_kind(ikind)))
2060 24 : radius_RI_kind(ikind) = SQRT(-LOG(eps)/exp_RI_kind(ikind))
2061 : END DO
2062 :
2063 8 : radius_ao = SQRT(-LOG(eps)/exp_min_ao)
2064 8 : radius_ao_product = SQRT(-LOG(eps)/(2.0_dp*exp_min_ao))
2065 8 : radius_RI = SQRT(-LOG(eps)/exp_min_RI)
2066 :
2067 8 : CALL get_atomic_kind_set(atomic_kind_set=atomic_kind_set, kind_of=kind_of)
2068 :
2069 : ! For a 3c integral (μR υS | P0) we have that cell R and cell S need to be within radius_3c
2070 8 : cell_radius_3c = radius_ao_product + radius_RI + bs_env%ri_metric%cutoff_radius
2071 :
2072 32 : n_max(1:3) = bs_env%periodic(1:3)*30
2073 :
2074 8 : nimages_3c_max = 0
2075 :
2076 8 : i_cell_x_min = 0
2077 8 : i_cell_x_max = 0
2078 8 : j_cell_y_min = 0
2079 8 : j_cell_y_max = 0
2080 8 : k_cell_z_min = 0
2081 8 : k_cell_z_max = 0
2082 :
2083 136 : DO i_cell_x = -n_max(1), n_max(1)
2084 7944 : DO j_cell_y = -n_max(2), n_max(2)
2085 37704 : DO k_cell_z = -n_max(3), n_max(3)
2086 :
2087 119072 : cell_index(1:3) = [i_cell_x, j_cell_y, k_cell_z]
2088 :
2089 29768 : CALL get_cell_dist(cell_index, bs_env%hmat, cell_dist)
2090 :
2091 37576 : IF (cell_dist < cell_radius_3c) THEN
2092 192 : nimages_3c_max = nimages_3c_max + 1
2093 192 : i_cell_x_min = MIN(i_cell_x_min, i_cell_x)
2094 192 : i_cell_x_max = MAX(i_cell_x_max, i_cell_x)
2095 192 : j_cell_y_min = MIN(j_cell_y_min, j_cell_y)
2096 192 : j_cell_y_max = MAX(j_cell_y_max, j_cell_y)
2097 192 : k_cell_z_min = MIN(k_cell_z_min, k_cell_z)
2098 192 : k_cell_z_max = MAX(k_cell_z_max, k_cell_z)
2099 : END IF
2100 :
2101 : END DO
2102 : END DO
2103 : END DO
2104 :
2105 : ! get index_to_cell_3c_max for the maximum possible cell range;
2106 : ! compute 3c integrals later in this routine and check really which cell is needed
2107 24 : ALLOCATE (index_to_cell_3c_max(3, nimages_3c_max))
2108 :
2109 8 : img = 0
2110 136 : DO i_cell_x = -n_max(1), n_max(1)
2111 7944 : DO j_cell_y = -n_max(2), n_max(2)
2112 37704 : DO k_cell_z = -n_max(3), n_max(3)
2113 :
2114 119072 : cell_index(1:3) = [i_cell_x, j_cell_y, k_cell_z]
2115 :
2116 29768 : CALL get_cell_dist(cell_index, bs_env%hmat, cell_dist)
2117 :
2118 37576 : IF (cell_dist < cell_radius_3c) THEN
2119 192 : img = img + 1
2120 768 : index_to_cell_3c_max(1:3, img) = cell_index(1:3)
2121 : END IF
2122 :
2123 : END DO
2124 : END DO
2125 : END DO
2126 :
2127 : ! get pairs of R and S which have non-zero 3c integral (μR υS | P0)
2128 32 : ALLOCATE (nblocks_3c_max(nimages_3c_max, nimages_3c_max))
2129 8 : nblocks_3c_max(:, :) = 0
2130 :
2131 8 : block_count = 0
2132 200 : DO j_cell = 1, nimages_3c_max
2133 4832 : DO k_cell = 1, nimages_3c_max
2134 :
2135 17838 : DO atom_j = 1, bs_env%n_atom
2136 54924 : DO atom_k = 1, bs_env%n_atom
2137 158598 : DO atom_i = 1, bs_env%n_atom
2138 :
2139 108306 : block_count = block_count + 1
2140 108306 : IF (MODULO(block_count, bs_env%para_env%num_pe) /= bs_env%para_env%mepos) CYCLE
2141 :
2142 216612 : CALL scaled_to_real(vec_cell_j, REAL(index_to_cell_3c_max(1:3, j_cell), kind=dp), cell)
2143 216612 : CALL scaled_to_real(vec_cell_k, REAL(index_to_cell_3c_max(1:3, k_cell), kind=dp), cell)
2144 :
2145 216612 : rij = pbc(particle_set(atom_j)%r(:), cell) - pbc(particle_set(atom_i)%r(:), cell) + vec_cell_j(:)
2146 : rjk = pbc(particle_set(atom_k)%r(:), cell) - pbc(particle_set(atom_j)%r(:), cell) &
2147 216612 : + vec_cell_k(:) - vec_cell_j(:)
2148 216612 : rik(:) = rij(:) + rjk(:)
2149 216612 : dij = NORM2(rij)
2150 216612 : dik = NORM2(rik)
2151 216612 : djk = NORM2(rjk)
2152 54153 : IF (djk > radius_ao_kind(kind_of(atom_j)) + radius_ao_kind(kind_of(atom_k))) CYCLE
2153 16971 : IF (dij > radius_ao_kind(kind_of(atom_j)) + radius_RI_kind(kind_of(atom_i)) &
2154 : + bs_env%ri_metric%cutoff_radius) CYCLE
2155 8645 : IF (dik > radius_RI_kind(kind_of(atom_i)) + radius_ao_kind(kind_of(atom_k)) &
2156 : + bs_env%ri_metric%cutoff_radius) CYCLE
2157 :
2158 5658 : j_size = bs_env%i_ao_end_from_atom(atom_j) - bs_env%i_ao_start_from_atom(atom_j) + 1
2159 5658 : k_size = bs_env%i_ao_end_from_atom(atom_k) - bs_env%i_ao_start_from_atom(atom_k) + 1
2160 5658 : i_size = bs_env%i_RI_end_from_atom(atom_i) - bs_env%i_RI_start_from_atom(atom_i) + 1
2161 :
2162 28290 : ALLOCATE (int_3c(j_size, k_size, i_size))
2163 :
2164 : ! compute 3-c int. ( μ(atom j) R , ν (atom k) S | P (atom i) 0 )
2165 : ! ("|": truncated Coulomb operator), inside build_3c_integrals: (j k | i)
2166 : CALL build_3c_integral_block(int_3c, qs_env, bs_env%ri_metric, &
2167 : basis_j=bs_env%basis_set_AO, &
2168 : basis_k=bs_env%basis_set_AO, &
2169 : basis_i=bs_env%basis_set_RI, &
2170 : cell_j=index_to_cell_3c_max(1:3, j_cell), &
2171 : cell_k=index_to_cell_3c_max(1:3, k_cell), &
2172 5658 : atom_k=atom_k, atom_j=atom_j, atom_i=atom_i)
2173 :
2174 300671 : frobenius_norm = SQRT(SUM(int_3c(:, :, :)**2))
2175 :
2176 5658 : DEALLOCATE (int_3c)
2177 :
2178 : ! we use a higher threshold here to safe memory when storing the 3c integrals
2179 : ! in every tensor group
2180 42936 : IF (frobenius_norm > eps) THEN
2181 1204 : nblocks_3c_max(j_cell, k_cell) = nblocks_3c_max(j_cell, k_cell) + 1
2182 : END IF
2183 :
2184 : END DO
2185 : END DO
2186 : END DO
2187 :
2188 : END DO
2189 : END DO
2190 :
2191 8 : CALL bs_env%para_env%sum(nblocks_3c_max)
2192 :
2193 24 : ALLOCATE (n_other_3c_images_max(nimages_3c_max))
2194 8 : n_other_3c_images_max(:) = 0
2195 :
2196 8 : nimages_3c = 0
2197 8 : nimage_pairs_3c = 0
2198 :
2199 200 : DO j_cell = 1, nimages_3c_max
2200 4824 : DO k_cell = 1, nimages_3c_max
2201 4824 : IF (nblocks_3c_max(j_cell, k_cell) > 0) THEN
2202 424 : n_other_3c_images_max(j_cell) = n_other_3c_images_max(j_cell) + 1
2203 424 : nimage_pairs_3c = nimage_pairs_3c + 1
2204 : END IF
2205 : END DO
2206 :
2207 200 : IF (n_other_3c_images_max(j_cell) > 0) nimages_3c = nimages_3c + 1
2208 :
2209 : END DO
2210 :
2211 8 : bs_env%nimages_3c = nimages_3c
2212 24 : ALLOCATE (bs_env%index_to_cell_3c(3, nimages_3c))
2213 : ALLOCATE (bs_env%cell_to_index_3c(i_cell_x_min:i_cell_x_max, &
2214 : j_cell_y_min:j_cell_y_max, &
2215 40 : k_cell_z_min:k_cell_z_max))
2216 400 : bs_env%cell_to_index_3c(:, :, :) = -1
2217 :
2218 32 : ALLOCATE (bs_env%nblocks_3c(nimages_3c, nimages_3c))
2219 8 : bs_env%nblocks_3c(nimages_3c, nimages_3c) = 0
2220 :
2221 8 : j_cell = 0
2222 200 : DO j_cell_max = 1, nimages_3c_max
2223 192 : IF (n_other_3c_images_max(j_cell_max) == 0) CYCLE
2224 82 : j_cell = j_cell + 1
2225 328 : cell_index(1:3) = index_to_cell_3c_max(1:3, j_cell_max)
2226 328 : bs_env%index_to_cell_3c(1:3, j_cell) = cell_index(1:3)
2227 82 : bs_env%cell_to_index_3c(cell_index(1), cell_index(2), cell_index(3)) = j_cell
2228 :
2229 82 : k_cell = 0
2230 2100 : DO k_cell_max = 1, nimages_3c_max
2231 2010 : IF (n_other_3c_images_max(k_cell_max) == 0) CYCLE
2232 914 : k_cell = k_cell + 1
2233 :
2234 2202 : bs_env%nblocks_3c(j_cell, k_cell) = nblocks_3c_max(j_cell_max, k_cell_max)
2235 : END DO
2236 :
2237 : END DO
2238 :
2239 : ! we use: 8*10^-9 GB / double precision number
2240 : mem_3c_GB = REAL(bs_env%n_RI, KIND=dp)*REAL(bs_env%n_ao, KIND=dp)**2 &
2241 8 : *REAL(nimage_pairs_3c, KIND=dp)*8E-9_dp
2242 :
2243 8 : CALL m_memory(mem_occ_per_proc)
2244 8 : CALL bs_env%para_env%max(mem_occ_per_proc)
2245 :
2246 8 : mem_occ_per_proc_GB = REAL(mem_occ_per_proc, KIND=dp)/1.0E9_dp
2247 :
2248 : ! number of processors per group that entirely stores the 3c integrals and does tensor ops
2249 8 : avail_mem_per_proc_GB = bs_env%input_memory_per_proc_GB - mem_occ_per_proc_GB
2250 :
2251 : ! careful: downconvering real to integer, 1.9 -> 1; thus add 1.0 for upconversion, 1.9 -> 2
2252 8 : bs_env%group_size_tensor = MAX(INT(mem_3c_GB/avail_mem_per_proc_GB + 1.0_dp), 1)
2253 :
2254 8 : u = bs_env%unit_nr
2255 :
2256 8 : IF (u > 0) THEN
2257 4 : WRITE (u, FMT="(T2,A,F52.1,A)") "Radius of atomic orbitals", radius_ao*angstrom, " Å"
2258 4 : WRITE (u, FMT="(T2,A,F55.1,A)") "Radius of RI functions", radius_RI*angstrom, " Å"
2259 4 : WRITE (u, FMT="(T2,A,I47)") "Number of cells for 3c integrals", nimages_3c
2260 4 : WRITE (u, FMT="(T2,A,I42)") "Number of cell pairs for 3c integrals", nimage_pairs_3c
2261 4 : WRITE (u, '(T2,A)') ''
2262 4 : WRITE (u, '(T2,A,F37.1,A)') 'Input: Available memory per MPI process', &
2263 8 : bs_env%input_memory_per_proc_GB, ' GB'
2264 4 : WRITE (u, '(T2,A,F35.1,A)') 'Used memory per MPI process before GW run', &
2265 8 : mem_occ_per_proc_GB, ' GB'
2266 4 : WRITE (u, '(T2,A,F44.1,A)') 'Memory of three-center integrals', mem_3c_GB, ' GB'
2267 : END IF
2268 :
2269 8 : CALL timestop(handle)
2270 :
2271 24 : END SUBROUTINE setup_cells_3c
2272 :
2273 : ! **************************************************************************************************
2274 : !> \brief ...
2275 : !> \param index_to_cell_1 ...
2276 : !> \param index_to_cell_2 ...
2277 : !> \param nimages_1 ...
2278 : !> \param nimages_2 ...
2279 : !> \param index_to_cell ...
2280 : !> \param cell_to_index ...
2281 : !> \param nimages ...
2282 : ! **************************************************************************************************
2283 8 : SUBROUTINE sum_two_R_grids(index_to_cell_1, index_to_cell_2, nimages_1, nimages_2, &
2284 : index_to_cell, cell_to_index, nimages)
2285 :
2286 : INTEGER, DIMENSION(:, :) :: index_to_cell_1, index_to_cell_2
2287 : INTEGER :: nimages_1, nimages_2
2288 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: index_to_cell
2289 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index
2290 : INTEGER :: nimages
2291 :
2292 : CHARACTER(LEN=*), PARAMETER :: routineN = 'sum_two_R_grids'
2293 :
2294 : INTEGER :: handle, i_dim, img_1, img_2, nimages_max
2295 8 : INTEGER, ALLOCATABLE, DIMENSION(:, :) :: index_to_cell_tmp
2296 : INTEGER, DIMENSION(3) :: cell_1, cell_2, R, R_max, R_min
2297 :
2298 8 : CALL timeset(routineN, handle)
2299 :
2300 32 : DO i_dim = 1, 3
2301 516 : R_min(i_dim) = MINVAL(index_to_cell_1(i_dim, :)) + MINVAL(index_to_cell_2(i_dim, :))
2302 548 : R_max(i_dim) = MAXVAL(index_to_cell_1(i_dim, :)) + MAXVAL(index_to_cell_2(i_dim, :))
2303 : END DO
2304 :
2305 8 : nimages_max = (R_max(1) - R_min(1) + 1)*(R_max(2) - R_min(2) + 1)*(R_max(3) - R_min(3) + 1)
2306 :
2307 24 : ALLOCATE (index_to_cell_tmp(3, nimages_max))
2308 1048 : index_to_cell_tmp(:, :) = -1
2309 :
2310 40 : ALLOCATE (cell_to_index(R_min(1):R_max(1), R_min(2):R_max(2), R_min(3):R_max(3)))
2311 532 : cell_to_index(:, :, :) = -1
2312 :
2313 8 : nimages = 0
2314 :
2315 90 : DO img_1 = 1, nimages_1
2316 :
2317 1004 : DO img_2 = 1, nimages_2
2318 :
2319 3656 : cell_1(1:3) = index_to_cell_1(1:3, img_1)
2320 3656 : cell_2(1:3) = index_to_cell_2(1:3, img_2)
2321 :
2322 3656 : R(1:3) = cell_1(1:3) + cell_2(1:3)
2323 :
2324 : ! check whether we have found a new cell
2325 996 : IF (cell_to_index(R(1), R(2), R(3)) == -1) THEN
2326 :
2327 236 : nimages = nimages + 1
2328 236 : cell_to_index(R(1), R(2), R(3)) = nimages
2329 944 : index_to_cell_tmp(1:3, nimages) = R(1:3)
2330 :
2331 : END IF
2332 :
2333 : END DO
2334 :
2335 : END DO
2336 :
2337 24 : ALLOCATE (index_to_cell(3, nimages))
2338 952 : index_to_cell(:, :) = index_to_cell_tmp(1:3, 1:nimages)
2339 :
2340 8 : CALL timestop(handle)
2341 :
2342 16 : END SUBROUTINE sum_two_R_grids
2343 :
2344 : ! **************************************************************************************************
2345 : !> \brief ...
2346 : !> \param qs_env ...
2347 : !> \param bs_env ...
2348 : ! **************************************************************************************************
2349 8 : SUBROUTINE compute_3c_integrals(qs_env, bs_env)
2350 :
2351 : TYPE(qs_environment_type), POINTER :: qs_env
2352 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2353 :
2354 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_3c_integrals'
2355 :
2356 : INTEGER :: handle, j_cell, k_cell, nimages_3c
2357 :
2358 8 : CALL timeset(routineN, handle)
2359 :
2360 8 : nimages_3c = bs_env%nimages_3c
2361 1092 : ALLOCATE (bs_env%t_3c_int(nimages_3c, nimages_3c))
2362 90 : DO j_cell = 1, nimages_3c
2363 1004 : DO k_cell = 1, nimages_3c
2364 996 : CALL dbt_create(bs_env%t_RI_AO__AO, bs_env%t_3c_int(j_cell, k_cell))
2365 : END DO
2366 : END DO
2367 :
2368 : CALL build_3c_integrals(bs_env%t_3c_int, &
2369 : bs_env%eps_filter, &
2370 : qs_env, &
2371 : bs_env%nl_3c, &
2372 : int_eps=bs_env%eps_filter*0.05_dp, &
2373 : basis_i=bs_env%basis_set_RI, &
2374 : basis_j=bs_env%basis_set_AO, &
2375 : basis_k=bs_env%basis_set_AO, &
2376 : potential_parameter=bs_env%ri_metric, &
2377 : desymmetrize=.FALSE., do_kpoints=.TRUE., cell_sym=.TRUE., &
2378 8 : cell_to_index_ext=bs_env%cell_to_index_3c)
2379 :
2380 8 : CALL bs_env%para_env%sync()
2381 :
2382 8 : CALL timestop(handle)
2383 :
2384 8 : END SUBROUTINE compute_3c_integrals
2385 :
2386 : ! **************************************************************************************************
2387 : !> \brief ...
2388 : !> \param cell_index ...
2389 : !> \param hmat ...
2390 : !> \param cell_dist ...
2391 : ! **************************************************************************************************
2392 59536 : SUBROUTINE get_cell_dist(cell_index, hmat, cell_dist)
2393 :
2394 : INTEGER, DIMENSION(3) :: cell_index
2395 : REAL(KIND=dp) :: hmat(3, 3), cell_dist
2396 :
2397 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_cell_dist'
2398 :
2399 : INTEGER :: handle, i_dim
2400 : INTEGER, DIMENSION(3) :: cell_index_adj
2401 : REAL(KIND=dp) :: cell_dist_3(3)
2402 :
2403 59536 : CALL timeset(routineN, handle)
2404 :
2405 : ! the distance of cells needs to be taken to adjacent neighbors, not
2406 : ! between the center of the cells. We thus need to rescale the cell index
2407 238144 : DO i_dim = 1, 3
2408 178608 : IF (cell_index(i_dim) > 0) cell_index_adj(i_dim) = cell_index(i_dim) - 1
2409 178608 : IF (cell_index(i_dim) < 0) cell_index_adj(i_dim) = cell_index(i_dim) + 1
2410 238144 : IF (cell_index(i_dim) == 0) cell_index_adj(i_dim) = cell_index(i_dim)
2411 : END DO
2412 :
2413 952576 : cell_dist_3(1:3) = MATMUL(hmat, REAL(cell_index_adj, KIND=dp))
2414 :
2415 238144 : cell_dist = SQRT(ABS(SUM(cell_dist_3(1:3)**2)))
2416 :
2417 59536 : CALL timestop(handle)
2418 :
2419 59536 : END SUBROUTINE get_cell_dist
2420 :
2421 : ! **************************************************************************************************
2422 : !> \brief ...
2423 : !> \param qs_env ...
2424 : !> \param bs_env ...
2425 : !> \param kpoints ...
2426 : !> \param do_print ...
2427 : ! **************************************************************************************************
2428 0 : SUBROUTINE setup_kpoints_scf_desymm(qs_env, bs_env, kpoints, do_print)
2429 : TYPE(qs_environment_type), POINTER :: qs_env
2430 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2431 : TYPE(kpoint_type), POINTER :: kpoints
2432 :
2433 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_kpoints_scf_desymm'
2434 :
2435 : INTEGER :: handle, i_cell_x, i_dim, img, j_cell_y, &
2436 : k_cell_z, nimages, nkp, u
2437 : INTEGER, DIMENSION(3) :: cell_grid, cixd, nkp_grid
2438 : TYPE(kpoint_type), POINTER :: kpoints_scf
2439 :
2440 : LOGICAL:: do_print
2441 :
2442 0 : CALL timeset(routineN, handle)
2443 :
2444 0 : NULLIFY (kpoints)
2445 0 : CALL kpoint_create(kpoints)
2446 :
2447 0 : CALL get_qs_env(qs_env=qs_env, kpoints=kpoints_scf)
2448 :
2449 0 : nkp_grid(1:3) = kpoints_scf%nkp_grid(1:3)
2450 0 : nkp = nkp_grid(1)*nkp_grid(2)*nkp_grid(3)
2451 :
2452 : ! we need in periodic directions at least 2 k-points in the SCF
2453 0 : DO i_dim = 1, 3
2454 0 : IF (bs_env%periodic(i_dim) == 1) THEN
2455 0 : CPASSERT(nkp_grid(i_dim) > 1)
2456 : END IF
2457 : END DO
2458 :
2459 0 : kpoints%kp_scheme = "GENERAL"
2460 0 : kpoints%nkp_grid(1:3) = nkp_grid(1:3)
2461 0 : kpoints%nkp = nkp
2462 0 : bs_env%nkp_scf_desymm = nkp
2463 :
2464 0 : ALLOCATE (kpoints%xkp(1:3, nkp))
2465 0 : CALL compute_xkp(kpoints%xkp, 1, nkp, nkp_grid)
2466 :
2467 0 : ALLOCATE (kpoints%wkp(nkp))
2468 0 : kpoints%wkp(:) = 1.0_dp/REAL(nkp, KIND=dp)
2469 :
2470 : ! for example 4x3x6 kpoint grid -> 3x3x5 cell grid because we need the same number of
2471 : ! neighbor cells on both sides of the unit cell
2472 0 : cell_grid(1:3) = nkp_grid(1:3) - MODULO(nkp_grid(1:3) + 1, 2)
2473 : ! cell index: for example for x: from -n_x/2 to +n_x/2, n_x: number of cells in x direction
2474 0 : cixd(1:3) = cell_grid(1:3)/2
2475 :
2476 0 : nimages = cell_grid(1)*cell_grid(2)*cell_grid(3)
2477 :
2478 0 : bs_env%nimages_scf_desymm = nimages
2479 :
2480 0 : ALLOCATE (kpoints%cell_to_index(-cixd(1):cixd(1), -cixd(2):cixd(2), -cixd(3):cixd(3)))
2481 0 : ALLOCATE (kpoints%index_to_cell(3, nimages))
2482 :
2483 0 : img = 0
2484 0 : DO i_cell_x = -cixd(1), cixd(1)
2485 0 : DO j_cell_y = -cixd(2), cixd(2)
2486 0 : DO k_cell_z = -cixd(3), cixd(3)
2487 0 : img = img + 1
2488 0 : kpoints%cell_to_index(i_cell_x, j_cell_y, k_cell_z) = img
2489 0 : kpoints%index_to_cell(1:3, img) = [i_cell_x, j_cell_y, k_cell_z]
2490 : END DO
2491 : END DO
2492 : END DO
2493 :
2494 0 : u = bs_env%unit_nr
2495 0 : IF (u > 0 .AND. do_print) THEN
2496 0 : WRITE (u, FMT="(T2,A,I49)") "Number of cells for G, χ, W, Σ", nimages
2497 : END IF
2498 :
2499 0 : CALL timestop(handle)
2500 :
2501 0 : END SUBROUTINE setup_kpoints_scf_desymm
2502 :
2503 : ! **************************************************************************************************
2504 : !> \brief ...
2505 : !> \param bs_env ...
2506 : ! **************************************************************************************************
2507 8 : SUBROUTINE setup_cells_Delta_R(bs_env)
2508 :
2509 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2510 :
2511 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_cells_Delta_R'
2512 :
2513 : INTEGER :: handle
2514 :
2515 8 : CALL timeset(routineN, handle)
2516 :
2517 : ! cell sums batch wise for fixed ΔR = S_1 - R_1; for example:
2518 : ! Σ_λσ^R = sum_PR1νS1 M^G_λ0,νS1,PR1 M^W_σR,νS1,PR1
2519 :
2520 : CALL sum_two_R_grids(bs_env%index_to_cell_3c, &
2521 : bs_env%index_to_cell_3c, &
2522 : bs_env%nimages_3c, bs_env%nimages_3c, &
2523 : bs_env%index_to_cell_Delta_R, &
2524 : bs_env%cell_to_index_Delta_R, &
2525 8 : bs_env%nimages_Delta_R)
2526 :
2527 8 : IF (bs_env%unit_nr > 0) THEN
2528 4 : WRITE (bs_env%unit_nr, FMT="(T2,A,I61)") "Number of cells ΔR", bs_env%nimages_Delta_R
2529 : END IF
2530 :
2531 8 : CALL timestop(handle)
2532 :
2533 8 : END SUBROUTINE setup_cells_Delta_R
2534 :
2535 : ! **************************************************************************************************
2536 : !> \brief ...
2537 : !> \param bs_env ...
2538 : ! **************************************************************************************************
2539 8 : SUBROUTINE setup_parallelization_Delta_R(bs_env)
2540 :
2541 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2542 :
2543 : CHARACTER(LEN=*), PARAMETER :: routineN = 'setup_parallelization_Delta_R'
2544 :
2545 : INTEGER :: handle, i_cell_Delta_R, i_task_local, &
2546 : n_tasks_local
2547 8 : INTEGER, ALLOCATABLE, DIMENSION(:) :: i_cell_Delta_R_group, &
2548 8 : n_tensor_ops_Delta_R
2549 :
2550 8 : CALL timeset(routineN, handle)
2551 :
2552 8 : CALL compute_n_tensor_ops_Delta_R(bs_env, n_tensor_ops_Delta_R)
2553 :
2554 8 : CALL compute_Delta_R_dist(bs_env, n_tensor_ops_Delta_R, i_cell_Delta_R_group, n_tasks_local)
2555 :
2556 8 : bs_env%n_tasks_Delta_R_local = n_tasks_local
2557 :
2558 24 : ALLOCATE (bs_env%task_Delta_R(n_tasks_local))
2559 :
2560 8 : i_task_local = 0
2561 244 : DO i_cell_Delta_R = 1, bs_env%nimages_Delta_R
2562 :
2563 236 : IF (i_cell_Delta_R_group(i_cell_Delta_R) /= bs_env%tensor_group_color) CYCLE
2564 :
2565 103 : i_task_local = i_task_local + 1
2566 :
2567 244 : bs_env%task_Delta_R(i_task_local) = i_cell_Delta_R
2568 :
2569 : END DO
2570 :
2571 16 : ALLOCATE (bs_env%skip_DR_chi(n_tasks_local))
2572 111 : bs_env%skip_DR_chi(:) = .FALSE.
2573 16 : ALLOCATE (bs_env%skip_DR_Sigma(n_tasks_local))
2574 111 : bs_env%skip_DR_Sigma(:) = .FALSE.
2575 :
2576 8 : CALL allocate_skip_3xR(bs_env%skip_DR_R12_S_Goccx3c_chi, bs_env)
2577 8 : CALL allocate_skip_3xR(bs_env%skip_DR_R12_S_Gvirx3c_chi, bs_env)
2578 8 : CALL allocate_skip_3xR(bs_env%skip_DR_R_R2_MxM_chi, bs_env)
2579 :
2580 8 : CALL allocate_skip_3xR(bs_env%skip_DR_R1_S2_Gx3c_Sigma, bs_env)
2581 8 : CALL allocate_skip_3xR(bs_env%skip_DR_R1_R_MxM_Sigma, bs_env)
2582 :
2583 8 : CALL timestop(handle)
2584 :
2585 16 : END SUBROUTINE setup_parallelization_Delta_R
2586 :
2587 : ! **************************************************************************************************
2588 : !> \brief ...
2589 : !> \param skip ...
2590 : !> \param bs_env ...
2591 : ! **************************************************************************************************
2592 40 : SUBROUTINE allocate_skip_3xR(skip, bs_env)
2593 : LOGICAL, ALLOCATABLE, DIMENSION(:, :, :) :: skip
2594 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2595 :
2596 : CHARACTER(LEN=*), PARAMETER :: routineN = 'allocate_skip_3xR'
2597 :
2598 : INTEGER :: handle
2599 :
2600 40 : CALL timeset(routineN, handle)
2601 :
2602 200 : ALLOCATE (skip(bs_env%n_tasks_Delta_R_local, bs_env%nimages_3c, bs_env%nimages_scf_desymm))
2603 40 : skip(:, :, :) = .FALSE.
2604 :
2605 40 : CALL timestop(handle)
2606 :
2607 40 : END SUBROUTINE allocate_skip_3xR
2608 :
2609 : ! **************************************************************************************************
2610 : !> \brief ...
2611 : !> \param bs_env ...
2612 : !> \param n_tensor_ops_Delta_R ...
2613 : !> \param i_cell_Delta_R_group ...
2614 : !> \param n_tasks_local ...
2615 : ! **************************************************************************************************
2616 8 : SUBROUTINE compute_Delta_R_dist(bs_env, n_tensor_ops_Delta_R, i_cell_Delta_R_group, n_tasks_local)
2617 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2618 : INTEGER, ALLOCATABLE, DIMENSION(:) :: n_tensor_ops_Delta_R, &
2619 : i_cell_Delta_R_group
2620 : INTEGER :: n_tasks_local
2621 :
2622 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_Delta_R_dist'
2623 :
2624 : INTEGER :: handle, i_Delta_R_max_op, i_group_min, &
2625 : nimages_Delta_R, u
2626 8 : INTEGER, ALLOCATABLE, DIMENSION(:) :: n_tensor_ops_Delta_R_in_group
2627 :
2628 8 : CALL timeset(routineN, handle)
2629 :
2630 8 : nimages_Delta_R = bs_env%nimages_Delta_R
2631 :
2632 8 : u = bs_env%unit_nr
2633 :
2634 8 : IF (u > 0 .AND. nimages_Delta_R < bs_env%num_tensor_groups) THEN
2635 0 : WRITE (u, FMT="(T2,A,I5,A,I5,A)") "There are only ", nimages_Delta_R, &
2636 0 : " tasks to work on but there are ", bs_env%num_tensor_groups, " groups."
2637 0 : WRITE (u, FMT="(T2,A)") "Please reduce the number of MPI processes."
2638 0 : WRITE (u, '(T2,A)') ''
2639 : END IF
2640 :
2641 24 : ALLOCATE (n_tensor_ops_Delta_R_in_group(bs_env%num_tensor_groups))
2642 8 : n_tensor_ops_Delta_R_in_group(:) = 0
2643 24 : ALLOCATE (i_cell_Delta_R_group(nimages_Delta_R))
2644 244 : i_cell_Delta_R_group(:) = -1
2645 :
2646 8 : n_tasks_local = 0
2647 :
2648 882 : DO WHILE (ANY(n_tensor_ops_Delta_R(:) /= 0))
2649 :
2650 : ! get largest element of n_tensor_ops_Delta_R
2651 6844 : i_Delta_R_max_op = MAXLOC(n_tensor_ops_Delta_R, 1)
2652 :
2653 : ! distribute i_Delta_R_max_op to tensor group which has currently the smallest load
2654 824 : i_group_min = MINLOC(n_tensor_ops_Delta_R_in_group, 1)
2655 :
2656 : ! the tensor groups are 0-index based; but i_group_min is 1-index based
2657 206 : i_cell_Delta_R_group(i_Delta_R_max_op) = i_group_min - 1
2658 : n_tensor_ops_Delta_R_in_group(i_group_min) = n_tensor_ops_Delta_R_in_group(i_group_min) + &
2659 206 : n_tensor_ops_Delta_R(i_Delta_R_max_op)
2660 :
2661 : ! remove i_Delta_R_max_op from n_tensor_ops_Delta_R
2662 206 : n_tensor_ops_Delta_R(i_Delta_R_max_op) = 0
2663 :
2664 214 : IF (bs_env%tensor_group_color == i_group_min - 1) n_tasks_local = n_tasks_local + 1
2665 :
2666 : END DO
2667 :
2668 8 : CALL timestop(handle)
2669 :
2670 16 : END SUBROUTINE compute_Delta_R_dist
2671 :
2672 : ! **************************************************************************************************
2673 : !> \brief ...
2674 : !> \param bs_env ...
2675 : !> \param n_tensor_ops_Delta_R ...
2676 : ! **************************************************************************************************
2677 8 : SUBROUTINE compute_n_tensor_ops_Delta_R(bs_env, n_tensor_ops_Delta_R)
2678 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2679 : INTEGER, ALLOCATABLE, DIMENSION(:) :: n_tensor_ops_Delta_R
2680 :
2681 : CHARACTER(LEN=*), PARAMETER :: routineN = 'compute_n_tensor_ops_Delta_R'
2682 :
2683 : INTEGER :: handle, i_cell_Delta_R, i_cell_R, i_cell_R1, i_cell_R1_minus_R, i_cell_R2, &
2684 : i_cell_R2_m_R1, i_cell_S1, i_cell_S1_m_R1_p_R2, i_cell_S1_minus_R, i_cell_S2, &
2685 : nimages_Delta_R
2686 : INTEGER, DIMENSION(3) :: cell_DR, cell_m_R1, cell_R, cell_R1, cell_R1_minus_R, cell_R2, &
2687 : cell_R2_m_R1, cell_S1, cell_S1_m_R2_p_R1, cell_S1_minus_R, cell_S1_p_S2_m_R1, cell_S2
2688 : LOGICAL :: cell_found
2689 :
2690 8 : CALL timeset(routineN, handle)
2691 :
2692 8 : nimages_Delta_R = bs_env%nimages_Delta_R
2693 :
2694 24 : ALLOCATE (n_tensor_ops_Delta_R(nimages_Delta_R))
2695 8 : n_tensor_ops_Delta_R(:) = 0
2696 :
2697 : ! compute number of tensor operations for specific Delta_R
2698 244 : DO i_cell_Delta_R = 1, nimages_Delta_R
2699 :
2700 236 : IF (MODULO(i_cell_Delta_R, bs_env%num_tensor_groups) /= bs_env%tensor_group_color) CYCLE
2701 :
2702 1451 : DO i_cell_R1 = 1, bs_env%nimages_3c
2703 :
2704 5300 : cell_R1(1:3) = bs_env%index_to_cell_3c(1:3, i_cell_R1)
2705 5300 : cell_DR(1:3) = bs_env%index_to_cell_Delta_R(1:3, i_cell_Delta_R)
2706 :
2707 : ! S_1 = R_1 + ΔR (from ΔR = S_1 - R_1)
2708 : CALL add_R(cell_R1, cell_DR, bs_env%index_to_cell_3c, cell_S1, &
2709 1325 : cell_found, bs_env%cell_to_index_3c, i_cell_S1)
2710 1325 : IF (.NOT. cell_found) CYCLE
2711 :
2712 4300 : DO i_cell_R2 = 1, bs_env%nimages_scf_desymm
2713 :
2714 15480 : cell_R2(1:3) = bs_env%kpoints_scf_desymm%index_to_cell(1:3, i_cell_R2)
2715 :
2716 : ! R_2 - R_1
2717 : CALL add_R(cell_R2, -cell_R1, bs_env%index_to_cell_3c, cell_R2_m_R1, &
2718 15480 : cell_found, bs_env%cell_to_index_3c, i_cell_R2_m_R1)
2719 3870 : IF (.NOT. cell_found) CYCLE
2720 :
2721 : ! S_1 - R_1 + R_2
2722 : CALL add_R(cell_S1, cell_R2_m_R1, bs_env%index_to_cell_3c, cell_S1_m_R2_p_R1, &
2723 2310 : cell_found, bs_env%cell_to_index_3c, i_cell_S1_m_R1_p_R2)
2724 2310 : IF (.NOT. cell_found) CYCLE
2725 :
2726 5836 : n_tensor_ops_Delta_R(i_cell_Delta_R) = n_tensor_ops_Delta_R(i_cell_Delta_R) + 1
2727 :
2728 : END DO ! i_cell_R2
2729 :
2730 4300 : DO i_cell_S2 = 1, bs_env%nimages_scf_desymm
2731 :
2732 15480 : cell_S2(1:3) = bs_env%kpoints_scf_desymm%index_to_cell(1:3, i_cell_S2)
2733 15480 : cell_m_R1(1:3) = -cell_R1(1:3)
2734 15480 : cell_S1_p_S2_m_R1(1:3) = cell_S1(1:3) + cell_S2(1:3) - cell_R1(1:3)
2735 :
2736 3870 : CALL is_cell_in_index_to_cell(cell_m_R1, bs_env%index_to_cell_3c, cell_found)
2737 3870 : IF (.NOT. cell_found) CYCLE
2738 :
2739 3141 : CALL is_cell_in_index_to_cell(cell_S1_p_S2_m_R1, bs_env%index_to_cell_3c, cell_found)
2740 430 : IF (.NOT. cell_found) CYCLE
2741 :
2742 : END DO ! i_cell_S2
2743 :
2744 5861 : DO i_cell_R = 1, bs_env%nimages_scf_desymm
2745 :
2746 15480 : cell_R = bs_env%kpoints_scf_desymm%index_to_cell(1:3, i_cell_R)
2747 :
2748 : ! R_1 - R
2749 : CALL add_R(cell_R1, -cell_R, bs_env%index_to_cell_3c, cell_R1_minus_R, &
2750 15480 : cell_found, bs_env%cell_to_index_3c, i_cell_R1_minus_R)
2751 3870 : IF (.NOT. cell_found) CYCLE
2752 :
2753 : ! S_1 - R
2754 : CALL add_R(cell_S1, -cell_R, bs_env%index_to_cell_3c, cell_S1_minus_R, &
2755 9996 : cell_found, bs_env%cell_to_index_3c, i_cell_S1_minus_R)
2756 1325 : IF (.NOT. cell_found) CYCLE
2757 :
2758 : END DO ! i_cell_R
2759 :
2760 : END DO ! i_cell_R1
2761 :
2762 : END DO ! i_cell_Delta_R
2763 :
2764 8 : CALL bs_env%para_env%sum(n_tensor_ops_Delta_R)
2765 :
2766 8 : CALL timestop(handle)
2767 :
2768 8 : END SUBROUTINE compute_n_tensor_ops_Delta_R
2769 :
2770 : ! **************************************************************************************************
2771 : !> \brief ...
2772 : !> \param cell_1 ...
2773 : !> \param cell_2 ...
2774 : !> \param index_to_cell ...
2775 : !> \param cell_1_plus_2 ...
2776 : !> \param cell_found ...
2777 : !> \param cell_to_index ...
2778 : !> \param i_cell_1_plus_2 ...
2779 : ! **************************************************************************************************
2780 123948 : SUBROUTINE add_R(cell_1, cell_2, index_to_cell, cell_1_plus_2, cell_found, &
2781 : cell_to_index, i_cell_1_plus_2)
2782 :
2783 : INTEGER, DIMENSION(3) :: cell_1, cell_2
2784 : INTEGER, DIMENSION(:, :) :: index_to_cell
2785 : INTEGER, DIMENSION(3) :: cell_1_plus_2
2786 : LOGICAL :: cell_found
2787 : INTEGER, DIMENSION(:, :, :), INTENT(IN), &
2788 : OPTIONAL, POINTER :: cell_to_index
2789 : INTEGER, INTENT(OUT), OPTIONAL :: i_cell_1_plus_2
2790 :
2791 : CHARACTER(LEN=*), PARAMETER :: routineN = 'add_R'
2792 :
2793 : INTEGER :: handle
2794 :
2795 123948 : CALL timeset(routineN, handle)
2796 :
2797 495792 : cell_1_plus_2(1:3) = cell_1(1:3) + cell_2(1:3)
2798 :
2799 123948 : CALL is_cell_in_index_to_cell(cell_1_plus_2, index_to_cell, cell_found)
2800 :
2801 123948 : IF (PRESENT(i_cell_1_plus_2)) THEN
2802 123948 : IF (cell_found) THEN
2803 70638 : CPASSERT(PRESENT(cell_to_index))
2804 70638 : i_cell_1_plus_2 = cell_to_index(cell_1_plus_2(1), cell_1_plus_2(2), cell_1_plus_2(3))
2805 : ELSE
2806 53310 : i_cell_1_plus_2 = -1000
2807 : END IF
2808 : END IF
2809 :
2810 123948 : CALL timestop(handle)
2811 :
2812 123948 : END SUBROUTINE add_R
2813 :
2814 : ! **************************************************************************************************
2815 : !> \brief ...
2816 : !> \param cell ...
2817 : !> \param index_to_cell ...
2818 : !> \param cell_found ...
2819 : ! **************************************************************************************************
2820 194559 : SUBROUTINE is_cell_in_index_to_cell(cell, index_to_cell, cell_found)
2821 : INTEGER, DIMENSION(3) :: cell
2822 : INTEGER, DIMENSION(:, :) :: index_to_cell
2823 : LOGICAL :: cell_found
2824 :
2825 : CHARACTER(LEN=*), PARAMETER :: routineN = 'is_cell_in_index_to_cell'
2826 :
2827 : INTEGER :: handle, i_cell, nimg
2828 : INTEGER, DIMENSION(3) :: cell_i
2829 :
2830 194559 : CALL timeset(routineN, handle)
2831 :
2832 194559 : nimg = SIZE(index_to_cell, 2)
2833 :
2834 194559 : cell_found = .FALSE.
2835 :
2836 2443734 : DO i_cell = 1, nimg
2837 :
2838 8996700 : cell_i(1:3) = index_to_cell(1:3, i_cell)
2839 :
2840 2443734 : IF (cell_i(1) == cell(1) .AND. cell_i(2) == cell(2) .AND. cell_i(3) == cell(3)) THEN
2841 116455 : cell_found = .TRUE.
2842 : END IF
2843 :
2844 : END DO
2845 :
2846 194559 : CALL timestop(handle)
2847 :
2848 194559 : END SUBROUTINE is_cell_in_index_to_cell
2849 :
2850 : ! **************************************************************************************************
2851 : !> \brief ...
2852 : !> \param qs_env ...
2853 : !> \param bs_env ...
2854 : ! **************************************************************************************************
2855 8 : SUBROUTINE allocate_matrices_small_cell_full_kp(qs_env, bs_env)
2856 : TYPE(qs_environment_type), POINTER :: qs_env
2857 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2858 :
2859 : CHARACTER(LEN=*), PARAMETER :: routineN = 'allocate_matrices_small_cell_full_kp'
2860 :
2861 : INTEGER :: handle, i_spin, i_t, img, n_spin, &
2862 : nimages_scf, num_time_freq_points
2863 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
2864 : TYPE(mp_para_env_type), POINTER :: para_env
2865 :
2866 8 : CALL timeset(routineN, handle)
2867 :
2868 8 : nimages_scf = bs_env%nimages_scf_desymm
2869 8 : num_time_freq_points = bs_env%num_time_freq_points
2870 8 : n_spin = bs_env%n_spin
2871 :
2872 8 : CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
2873 :
2874 96 : ALLOCATE (bs_env%fm_G_S(nimages_scf))
2875 88 : ALLOCATE (bs_env%fm_Sigma_x_R(nimages_scf))
2876 592 : ALLOCATE (bs_env%fm_chi_R_t(nimages_scf, num_time_freq_points))
2877 584 : ALLOCATE (bs_env%fm_MWM_R_t(nimages_scf, num_time_freq_points))
2878 608 : ALLOCATE (bs_env%fm_Sigma_c_R_neg_tau(nimages_scf, num_time_freq_points, n_spin))
2879 600 : ALLOCATE (bs_env%fm_Sigma_c_R_pos_tau(nimages_scf, num_time_freq_points, n_spin))
2880 80 : DO img = 1, nimages_scf
2881 72 : CALL cp_fm_create(bs_env%fm_G_S(img), bs_env%fm_work_mo(1)%matrix_struct)
2882 72 : CALL cp_fm_create(bs_env%fm_Sigma_x_R(img), bs_env%fm_work_mo(1)%matrix_struct)
2883 584 : DO i_t = 1, num_time_freq_points
2884 504 : CALL cp_fm_create(bs_env%fm_chi_R_t(img, i_t), bs_env%fm_RI_RI%matrix_struct)
2885 504 : CALL cp_fm_create(bs_env%fm_MWM_R_t(img, i_t), bs_env%fm_RI_RI%matrix_struct)
2886 504 : CALL cp_fm_set_all(bs_env%fm_MWM_R_t(img, i_t), 0.0_dp)
2887 1080 : DO i_spin = 1, n_spin
2888 : CALL cp_fm_create(bs_env%fm_Sigma_c_R_neg_tau(img, i_t, i_spin), &
2889 504 : bs_env%fm_work_mo(1)%matrix_struct)
2890 : CALL cp_fm_create(bs_env%fm_Sigma_c_R_pos_tau(img, i_t, i_spin), &
2891 504 : bs_env%fm_work_mo(1)%matrix_struct)
2892 504 : CALL cp_fm_set_all(bs_env%fm_Sigma_c_R_neg_tau(img, i_t, i_spin), 0.0_dp)
2893 1008 : CALL cp_fm_set_all(bs_env%fm_Sigma_c_R_pos_tau(img, i_t, i_spin), 0.0_dp)
2894 : END DO
2895 : END DO
2896 : END DO
2897 :
2898 8 : CALL timestop(handle)
2899 :
2900 8 : END SUBROUTINE allocate_matrices_small_cell_full_kp
2901 :
2902 : ! **************************************************************************************************
2903 : !> \brief ...
2904 : !> \param qs_env ...
2905 : !> \param bs_env ...
2906 : ! **************************************************************************************************
2907 8 : SUBROUTINE trafo_V_xc_R_to_kp(qs_env, bs_env)
2908 : TYPE(qs_environment_type), POINTER :: qs_env
2909 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2910 :
2911 : CHARACTER(LEN=*), PARAMETER :: routineN = 'trafo_V_xc_R_to_kp'
2912 :
2913 : INTEGER :: handle, ikp, img, ispin, n_ao
2914 8 : INTEGER, DIMENSION(:, :, :), POINTER :: cell_to_index_scf
2915 : TYPE(cp_cfm_type) :: cfm_mo_coeff, cfm_tmp, cfm_V_xc
2916 : TYPE(cp_fm_type) :: fm_V_xc_re
2917 8 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks
2918 : TYPE(kpoint_type), POINTER :: kpoints_scf
2919 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
2920 8 : POINTER :: sab_nl
2921 :
2922 8 : CALL timeset(routineN, handle)
2923 :
2924 8 : n_ao = bs_env%n_ao
2925 :
2926 8 : CALL get_qs_env(qs_env, matrix_ks_kp=matrix_ks, kpoints=kpoints_scf)
2927 :
2928 8 : NULLIFY (sab_nl)
2929 8 : CALL get_kpoint_info(kpoints_scf, sab_nl=sab_nl, cell_to_index=cell_to_index_scf)
2930 :
2931 8 : CALL cp_cfm_create(cfm_V_xc, bs_env%cfm_work_mo%matrix_struct)
2932 8 : CALL cp_cfm_create(cfm_mo_coeff, bs_env%cfm_work_mo%matrix_struct)
2933 8 : CALL cp_cfm_create(cfm_tmp, bs_env%cfm_work_mo%matrix_struct)
2934 8 : CALL cp_fm_create(fm_V_xc_re, bs_env%cfm_work_mo%matrix_struct)
2935 :
2936 256 : DO img = 1, bs_env%nimages_scf
2937 504 : DO ispin = 1, bs_env%n_spin
2938 : ! JW kind of hack because the format of matrix_ks remains dubious...
2939 248 : CALL dbcsr_set(matrix_ks(ispin, img)%matrix, 0.0_dp)
2940 496 : CALL copy_fm_to_dbcsr(bs_env%fm_V_xc_R(img, ispin), matrix_ks(ispin, img)%matrix)
2941 : END DO
2942 : END DO
2943 :
2944 40 : ALLOCATE (bs_env%v_xc_n(n_ao, bs_env%nkp_bs_and_DOS, bs_env%n_spin))
2945 :
2946 16 : DO ispin = 1, bs_env%n_spin
2947 206 : DO ikp = 1, bs_env%nkp_bs_and_DOS
2948 :
2949 : ! v^xc^R -> v^xc(k) (matrix_ks stores v^xc^R, see SUBROUTINE compute_V_xc)
2950 : CALL rsmat_to_kp(matrix_ks, ispin, bs_env%kpoints_DOS%xkp(1:3, ikp), &
2951 190 : cell_to_index_scf, sab_nl, bs_env, cfm_V_xc)
2952 :
2953 : ! get C_µn(k)
2954 190 : CALL cp_cfm_to_cfm(bs_env%cfm_mo_coeff_kp(ikp, ispin), cfm_mo_coeff)
2955 :
2956 : ! v^xc_nm(k_i) = sum_µν C^*_µn(k_i) v^xc_µν(k_i) C_νn(k_i)
2957 : CALL parallel_gemm('N', 'N', n_ao, n_ao, n_ao, z_one, cfm_V_xc, cfm_mo_coeff, &
2958 190 : z_zero, cfm_tmp)
2959 : CALL parallel_gemm('C', 'N', n_ao, n_ao, n_ao, z_one, cfm_mo_coeff, cfm_tmp, &
2960 190 : z_zero, cfm_V_xc)
2961 :
2962 : ! get v^xc_nn(k_i) which is a real quantity as v^xc is Hermitian
2963 190 : CALL cp_cfm_to_fm(cfm_V_xc, fm_V_xc_re)
2964 198 : CALL cp_fm_get_diag(fm_V_xc_re, bs_env%v_xc_n(:, ikp, ispin))
2965 :
2966 : END DO
2967 :
2968 : END DO
2969 :
2970 : ! just rebuild the overwritten KS matrix again
2971 8 : CALL qs_ks_build_kohn_sham_matrix(qs_env, calculate_forces=.FALSE., just_energy=.FALSE.)
2972 :
2973 8 : CALL cp_cfm_release(cfm_V_xc)
2974 8 : CALL cp_cfm_release(cfm_mo_coeff)
2975 8 : CALL cp_cfm_release(cfm_tmp)
2976 8 : CALL cp_fm_release(fm_V_xc_re)
2977 :
2978 8 : CALL timestop(handle)
2979 :
2980 16 : END SUBROUTINE trafo_V_xc_R_to_kp
2981 :
2982 : ! **************************************************************************************************
2983 : !> \brief ...
2984 : !> \param qs_env ...
2985 : !> \param bs_env ...
2986 : ! **************************************************************************************************
2987 8 : SUBROUTINE heuristic_RI_regularization(qs_env, bs_env)
2988 : TYPE(qs_environment_type), POINTER :: qs_env
2989 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
2990 :
2991 : CHARACTER(LEN=*), PARAMETER :: routineN = 'heuristic_RI_regularization'
2992 :
2993 8 : COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: M
2994 : INTEGER :: handle, ikp, ikp_local, n_RI, nkp, &
2995 : nkp_local, u
2996 : REAL(KIND=dp) :: cond_nr, cond_nr_max, max_ev, &
2997 : max_ev_ikp, min_ev, min_ev_ikp
2998 8 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: M_R
2999 :
3000 8 : CALL timeset(routineN, handle)
3001 :
3002 : ! compute M^R_PQ = <phi_P,0|V^tr(rc)|phi_Q,R> for RI metric
3003 8 : CALL get_V_tr_R(M_R, bs_env%ri_metric, 0.0_dp, bs_env, qs_env)
3004 :
3005 8 : nkp = bs_env%nkp_chi_eps_W_orig_plus_extra
3006 8 : n_RI = bs_env%n_RI
3007 :
3008 8 : nkp_local = 0
3009 5128 : DO ikp = 1, nkp
3010 : ! trivial parallelization over k-points
3011 5120 : IF (MODULO(ikp, bs_env%para_env%num_pe) /= bs_env%para_env%mepos) CYCLE
3012 5128 : nkp_local = nkp_local + 1
3013 : END DO
3014 :
3015 40 : ALLOCATE (M(n_RI, n_RI, nkp_local))
3016 :
3017 8 : ikp_local = 0
3018 8 : cond_nr_max = 0.0_dp
3019 8 : min_ev = 1000.0_dp
3020 8 : max_ev = -1000.0_dp
3021 :
3022 5128 : DO ikp = 1, nkp
3023 :
3024 : ! trivial parallelization
3025 5120 : IF (MODULO(ikp, bs_env%para_env%num_pe) /= bs_env%para_env%mepos) CYCLE
3026 :
3027 2560 : ikp_local = ikp_local + 1
3028 :
3029 : ! M(k) = sum_R e^ikR M^R
3030 : CALL rs_to_kp(M_R, M(:, :, ikp_local), &
3031 : bs_env%kpoints_scf_desymm%index_to_cell, &
3032 2560 : bs_env%kpoints_chi_eps_W%xkp(1:3, ikp))
3033 :
3034 : ! compute condition number of M_PQ(k)
3035 2560 : CALL power(M(:, :, ikp_local), 1.0_dp, 0.0_dp, cond_nr, min_ev_ikp, max_ev_ikp)
3036 :
3037 2560 : IF (cond_nr > cond_nr_max) cond_nr_max = cond_nr
3038 2560 : IF (max_ev_ikp > max_ev) max_ev = max_ev_ikp
3039 2568 : IF (min_ev_ikp < min_ev) min_ev = min_ev_ikp
3040 :
3041 : END DO ! ikp
3042 :
3043 8 : CALL bs_env%para_env%max(cond_nr_max)
3044 8 : CALL bs_env%para_env%min(min_ev)
3045 8 : CALL bs_env%para_env%max(max_ev)
3046 :
3047 8 : u = bs_env%unit_nr
3048 8 : IF (u > 0) THEN
3049 4 : WRITE (u, FMT="(T2,A,ES34.1)") "Min. abs. eigenvalue of RI metric matrix M(k)", min_ev
3050 4 : WRITE (u, FMT="(T2,A,ES34.1)") "Max. abs. eigenvalue of RI metric matrix M(k)", max_ev
3051 4 : WRITE (u, FMT="(T2,A,ES50.1)") "Max. condition number of M(k)", cond_nr_max
3052 : END IF
3053 :
3054 8 : CALL timestop(handle)
3055 :
3056 16 : END SUBROUTINE heuristic_RI_regularization
3057 :
3058 : ! **************************************************************************************************
3059 : !> \brief ...
3060 : !> \param V_tr_R ...
3061 : !> \param pot_type ...
3062 : !> \param regularization_RI ...
3063 : !> \param bs_env ...
3064 : !> \param qs_env ...
3065 : ! **************************************************************************************************
3066 88 : SUBROUTINE get_V_tr_R(V_tr_R, pot_type, regularization_RI, bs_env, qs_env)
3067 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :, :) :: V_tr_R
3068 : TYPE(libint_potential_type) :: pot_type
3069 : REAL(KIND=dp) :: regularization_RI
3070 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
3071 : TYPE(qs_environment_type), POINTER :: qs_env
3072 :
3073 : CHARACTER(LEN=*), PARAMETER :: routineN = 'get_V_tr_R'
3074 :
3075 : INTEGER :: handle, img, nimages_scf_desymm
3076 : INTEGER, ALLOCATABLE, DIMENSION(:) :: sizes_RI
3077 88 : INTEGER, DIMENSION(:), POINTER :: col_bsize, row_bsize
3078 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
3079 88 : TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: fm_V_tr_R
3080 : TYPE(dbcsr_distribution_type) :: dbcsr_dist
3081 88 : TYPE(dbcsr_type), ALLOCATABLE, DIMENSION(:) :: mat_V_tr_R
3082 : TYPE(distribution_2d_type), POINTER :: dist_2d
3083 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
3084 88 : POINTER :: sab_RI
3085 88 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
3086 88 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
3087 :
3088 88 : CALL timeset(routineN, handle)
3089 :
3090 88 : NULLIFY (sab_RI, dist_2d)
3091 :
3092 : CALL get_qs_env(qs_env=qs_env, &
3093 : blacs_env=blacs_env, &
3094 : distribution_2d=dist_2d, &
3095 : qs_kind_set=qs_kind_set, &
3096 88 : particle_set=particle_set)
3097 :
3098 264 : ALLOCATE (sizes_RI(bs_env%n_atom))
3099 88 : CALL get_particle_set(particle_set, qs_kind_set, nsgf=sizes_RI, basis=bs_env%basis_set_RI)
3100 : CALL build_2c_neighbor_lists(sab_RI, bs_env%basis_set_RI, bs_env%basis_set_RI, &
3101 : pot_type, "2c_nl_RI", qs_env, sym_ij=.FALSE., &
3102 88 : dist_2d=dist_2d)
3103 88 : CALL cp_dbcsr_dist2d_to_dist(dist_2d, dbcsr_dist)
3104 264 : ALLOCATE (row_bsize(SIZE(sizes_RI)))
3105 176 : ALLOCATE (col_bsize(SIZE(sizes_RI)))
3106 324 : row_bsize(:) = sizes_RI
3107 324 : col_bsize(:) = sizes_RI
3108 :
3109 88 : nimages_scf_desymm = bs_env%nimages_scf_desymm
3110 1056 : ALLOCATE (mat_V_tr_R(nimages_scf_desymm))
3111 : CALL dbcsr_create(mat_V_tr_R(1), "(RI|RI)", dbcsr_dist, dbcsr_type_no_symmetry, &
3112 88 : row_bsize, col_bsize)
3113 88 : DEALLOCATE (row_bsize, col_bsize)
3114 :
3115 792 : DO img = 2, nimages_scf_desymm
3116 792 : CALL dbcsr_create(mat_V_tr_R(img), template=mat_V_tr_R(1))
3117 : END DO
3118 :
3119 : CALL build_2c_integrals(mat_V_tr_R, 0.0_dp, qs_env, sab_RI, bs_env%basis_set_RI, &
3120 : bs_env%basis_set_RI, pot_type, do_kpoints=.TRUE., &
3121 : ext_kpoints=bs_env%kpoints_scf_desymm, &
3122 88 : regularization_RI=regularization_RI)
3123 :
3124 1056 : ALLOCATE (fm_V_tr_R(nimages_scf_desymm))
3125 880 : DO img = 1, nimages_scf_desymm
3126 792 : CALL cp_fm_create(fm_V_tr_R(img), bs_env%fm_RI_RI%matrix_struct)
3127 792 : CALL copy_dbcsr_to_fm(mat_V_tr_R(img), fm_V_tr_R(img))
3128 880 : CALL dbcsr_release(mat_V_tr_R(img))
3129 : END DO
3130 :
3131 88 : IF (.NOT. ALLOCATED(V_tr_R)) THEN
3132 440 : ALLOCATE (V_tr_R(bs_env%n_RI, bs_env%n_RI, nimages_scf_desymm))
3133 : END IF
3134 :
3135 88 : CALL fm_to_local_array(fm_V_tr_R, V_tr_R)
3136 :
3137 88 : CALL cp_fm_release(fm_V_tr_R)
3138 88 : CALL dbcsr_distribution_release(dbcsr_dist)
3139 88 : CALL release_neighbor_list_sets(sab_RI)
3140 :
3141 88 : CALL timestop(handle)
3142 :
3143 264 : END SUBROUTINE get_V_tr_R
3144 :
3145 : ! **************************************************************************************************
3146 : !> \brief ...
3147 : !> \param matrix ...
3148 : !> \param exponent ...
3149 : !> \param eps ...
3150 : !> \param cond_nr ...
3151 : !> \param min_ev ...
3152 : !> \param max_ev ...
3153 : ! **************************************************************************************************
3154 44032 : SUBROUTINE power(matrix, exponent, eps, cond_nr, min_ev, max_ev)
3155 : COMPLEX(KIND=dp), DIMENSION(:, :) :: matrix
3156 : REAL(KIND=dp) :: exponent, eps
3157 : REAL(KIND=dp), OPTIONAL :: cond_nr, min_ev, max_ev
3158 :
3159 : CHARACTER(len=*), PARAMETER :: routineN = 'power'
3160 :
3161 44032 : COMPLEX(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: eigenvectors
3162 : INTEGER :: handle, i, n
3163 : REAL(KIND=dp) :: pos_eval
3164 44032 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: eigenvalues
3165 :
3166 44032 : CALL timeset(routineN, handle)
3167 :
3168 : ! make matrix perfectly Hermitian
3169 3385216 : matrix(:, :) = 0.5_dp*(matrix(:, :) + CONJG(TRANSPOSE(matrix(:, :))))
3170 :
3171 44032 : n = SIZE(matrix, 1)
3172 264192 : ALLOCATE (eigenvalues(n), eigenvectors(n, n))
3173 44032 : CALL diag_complex(matrix, eigenvectors, eigenvalues)
3174 :
3175 73472 : IF (PRESENT(cond_nr)) cond_nr = MAXVAL(ABS(eigenvalues))/MINVAL(ABS(eigenvalues))
3176 58752 : IF (PRESENT(min_ev)) min_ev = MINVAL(ABS(eigenvalues))
3177 58752 : IF (PRESENT(max_ev)) max_ev = MAXVAL(ABS(eigenvalues))
3178 :
3179 293328 : DO i = 1, n
3180 249296 : IF (eps < eigenvalues(i)) THEN
3181 249296 : pos_eval = (eigenvalues(i))**(0.5_dp*exponent)
3182 : ELSE
3183 : pos_eval = 0.0_dp
3184 : END IF
3185 1714624 : eigenvectors(:, i) = eigenvectors(:, i)*pos_eval
3186 : END DO
3187 :
3188 44032 : CALL ZGEMM("N", "C", n, n, n, z_one, eigenvectors, n, eigenvectors, n, z_zero, matrix, n)
3189 :
3190 44032 : DEALLOCATE (eigenvalues, eigenvectors)
3191 :
3192 44032 : CALL timestop(handle)
3193 :
3194 44032 : END SUBROUTINE power
3195 :
3196 : ! **************************************************************************************************
3197 : !> \brief ...
3198 : !> \param bs_env ...
3199 : !> \param Sigma_c_n_time ...
3200 : !> \param Sigma_c_n_freq ...
3201 : !> \param ispin ...
3202 : ! **************************************************************************************************
3203 242 : SUBROUTINE time_to_freq(bs_env, Sigma_c_n_time, Sigma_c_n_freq, ispin)
3204 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
3205 : REAL(KIND=dp), DIMENSION(:, :, :) :: Sigma_c_n_time, Sigma_c_n_freq
3206 : INTEGER :: ispin
3207 :
3208 : CHARACTER(LEN=*), PARAMETER :: routineN = 'time_to_freq'
3209 :
3210 : INTEGER :: handle, i_t, j_w, n_occ
3211 : REAL(KIND=dp) :: freq_j, time_i, w_cos_ij, w_sin_ij
3212 242 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: Sigma_c_n_cos_time, Sigma_c_n_sin_time
3213 :
3214 242 : CALL timeset(routineN, handle)
3215 :
3216 968 : ALLOCATE (Sigma_c_n_cos_time(bs_env%n_ao, bs_env%num_time_freq_points))
3217 726 : ALLOCATE (Sigma_c_n_sin_time(bs_env%n_ao, bs_env%num_time_freq_points))
3218 :
3219 23234 : Sigma_c_n_cos_time(:, :) = 0.5_dp*(Sigma_c_n_time(:, :, 1) + Sigma_c_n_time(:, :, 2))
3220 23234 : Sigma_c_n_sin_time(:, :) = 0.5_dp*(Sigma_c_n_time(:, :, 1) - Sigma_c_n_time(:, :, 2))
3221 :
3222 46710 : Sigma_c_n_freq(:, :, :) = 0.0_dp
3223 :
3224 2396 : DO i_t = 1, bs_env%num_time_freq_points
3225 :
3226 24390 : DO j_w = 1, bs_env%num_time_freq_points
3227 :
3228 21994 : freq_j = bs_env%imag_freq_points(j_w)
3229 21994 : time_i = bs_env%imag_time_points(i_t)
3230 : ! integration weights for cosine and sine transform
3231 21994 : w_cos_ij = bs_env%weights_cos_t_to_w(j_w, i_t)*COS(freq_j*time_i)
3232 21994 : w_sin_ij = bs_env%weights_sin_t_to_w(j_w, i_t)*SIN(freq_j*time_i)
3233 :
3234 : ! 1. Re(Σ^c_nn(k_i,iω)) from cosine transform
3235 : Sigma_c_n_freq(:, j_w, 1) = Sigma_c_n_freq(:, j_w, 1) + &
3236 200352 : w_cos_ij*Sigma_c_n_cos_time(:, i_t)
3237 :
3238 : ! 2. Im(Σ^c_nn(k_i,iω)) from sine transform
3239 : Sigma_c_n_freq(:, j_w, 2) = Sigma_c_n_freq(:, j_w, 2) + &
3240 202506 : w_sin_ij*Sigma_c_n_sin_time(:, i_t)
3241 :
3242 : END DO
3243 :
3244 : END DO
3245 :
3246 : ! for occupied levels, we need the correlation self-energy for negative omega.
3247 : ! Therefore, weight_sin should be computed with -omega, which results in an
3248 : ! additional minus for the imaginary part:
3249 242 : n_occ = bs_env%n_occ(ispin)
3250 10322 : Sigma_c_n_freq(1:n_occ, :, 2) = -Sigma_c_n_freq(1:n_occ, :, 2)
3251 :
3252 242 : CALL timestop(handle)
3253 :
3254 484 : END SUBROUTINE time_to_freq
3255 :
3256 : ! **************************************************************************************************
3257 : !> \brief ...
3258 : !> \param bs_env ...
3259 : !> \param Sigma_c_ikp_n_freq ...
3260 : !> \param Sigma_x_ikp_n ...
3261 : !> \param V_xc_ikp_n ...
3262 : !> \param eigenval_scf ...
3263 : !> \param ikp ...
3264 : !> \param ispin ...
3265 : ! **************************************************************************************************
3266 242 : SUBROUTINE analyt_conti_and_print(bs_env, Sigma_c_ikp_n_freq, Sigma_x_ikp_n, V_xc_ikp_n, &
3267 242 : eigenval_scf, ikp, ispin)
3268 :
3269 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
3270 : REAL(KIND=dp), DIMENSION(:, :, :) :: Sigma_c_ikp_n_freq
3271 : REAL(KIND=dp), DIMENSION(:) :: Sigma_x_ikp_n, V_xc_ikp_n, eigenval_scf
3272 : INTEGER :: ikp, ispin
3273 :
3274 : CHARACTER(LEN=*), PARAMETER :: routineN = 'analyt_conti_and_print'
3275 :
3276 : CHARACTER(len=3) :: occ_vir
3277 : CHARACTER(len=default_path_length) :: fname
3278 : INTEGER :: handle, i_mo, ikp_for_print, iunit, &
3279 : n_mo, nkp
3280 : LOGICAL :: is_bandstruc_kpoint, print_DOS_kpoints, &
3281 : print_ikp
3282 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: dummy, Sigma_c_ikp_n_qp
3283 :
3284 242 : CALL timeset(routineN, handle)
3285 :
3286 242 : n_mo = bs_env%n_ao
3287 968 : ALLOCATE (dummy(n_mo), Sigma_c_ikp_n_qp(n_mo))
3288 242 : Sigma_c_ikp_n_qp(:) = 0.0_dp
3289 :
3290 2902 : DO i_mo = 1, n_mo
3291 :
3292 : ! parallelization
3293 2660 : IF (MODULO(i_mo, bs_env%para_env%num_pe) /= bs_env%para_env%mepos) CYCLE
3294 :
3295 : CALL continuation_pade(Sigma_c_ikp_n_qp, &
3296 : bs_env%imag_freq_points_fit, dummy, dummy, &
3297 : Sigma_c_ikp_n_freq(:, 1:bs_env%num_freq_points_fit, 1)*z_one + &
3298 : Sigma_c_ikp_n_freq(:, 1:bs_env%num_freq_points_fit, 2)*gaussi, &
3299 : Sigma_x_ikp_n(:) - V_xc_ikp_n(:), &
3300 : eigenval_scf(:), eigenval_scf(:), &
3301 : bs_env%do_hedin_shift, &
3302 : i_mo, bs_env%n_occ(ispin), bs_env%n_vir(ispin), &
3303 : bs_env%nparam_pade, bs_env%num_freq_points_fit, &
3304 : ri_rpa_g0w0_crossing_newton, bs_env%n_occ(ispin), &
3305 86712 : 0.0_dp, .TRUE., .FALSE., 1, e_fermi_ext=bs_env%e_fermi(ispin))
3306 : END DO
3307 :
3308 242 : CALL bs_env%para_env%sum(Sigma_c_ikp_n_qp)
3309 :
3310 242 : CALL correct_obvious_fitting_fails(Sigma_c_ikp_n_qp, ispin, bs_env)
3311 :
3312 : bs_env%eigenval_G0W0(:, ikp, ispin) = eigenval_scf(:) + &
3313 : Sigma_c_ikp_n_qp(:) + &
3314 : Sigma_x_ikp_n(:) - &
3315 2902 : V_xc_ikp_n(:)
3316 :
3317 2902 : bs_env%eigenval_HF(:, ikp, ispin) = eigenval_scf(:) + Sigma_x_ikp_n(:) - V_xc_ikp_n(:)
3318 :
3319 : ! only print eigenvalues of DOS k-points in case no bandstructure path has been given
3320 242 : print_DOS_kpoints = (bs_env%nkp_only_bs <= 0)
3321 : ! in kpoints_DOS, the last nkp_only_bs are bandstructure k-points
3322 242 : is_bandstruc_kpoint = (ikp > bs_env%nkp_only_DOS)
3323 242 : print_ikp = print_DOS_kpoints .OR. is_bandstruc_kpoint
3324 :
3325 242 : IF (bs_env%para_env%is_source() .AND. print_ikp) THEN
3326 :
3327 105 : IF (print_DOS_kpoints) THEN
3328 74 : nkp = bs_env%nkp_only_DOS
3329 74 : ikp_for_print = ikp
3330 : ELSE
3331 31 : nkp = bs_env%nkp_only_bs
3332 31 : ikp_for_print = ikp - bs_env%nkp_only_DOS
3333 : END IF
3334 :
3335 105 : fname = "bandstructure_SCF_and_G0W0"
3336 :
3337 105 : IF (ikp_for_print == 1 .AND. ispin == 1) THEN
3338 : CALL open_file(TRIM(fname), unit_number=iunit, file_status="REPLACE", &
3339 21 : file_action="WRITE")
3340 : ELSE
3341 : CALL open_file(TRIM(fname), unit_number=iunit, file_status="OLD", &
3342 84 : file_action="WRITE", file_position="APPEND")
3343 : END IF
3344 :
3345 105 : WRITE (iunit, "(A)") " "
3346 105 : WRITE (iunit, "(A10,I7,A25,3F10.4,T90,A7,I2)") "kpoint: ", ikp_for_print, "coordinate: ", &
3347 525 : bs_env%kpoints_DOS%xkp(:, ikp), "spin: ", ispin
3348 105 : WRITE (iunit, "(A)") " "
3349 105 : WRITE (iunit, "(A5,A12,3A17,A16,A18)") "n", "k", "ϵ_nk^DFT (eV)", "Σ^c_nk (eV)", &
3350 210 : "Σ^x_nk (eV)", "v_nk^xc (eV)", "ϵ_nk^G0W0 (eV)"
3351 105 : WRITE (iunit, "(A)") " "
3352 :
3353 1291 : DO i_mo = 1, n_mo
3354 1186 : IF (i_mo <= bs_env%n_occ(ispin)) occ_vir = 'occ'
3355 1186 : IF (i_mo > bs_env%n_occ(ispin)) occ_vir = 'vir'
3356 1186 : WRITE (iunit, "(I5,3A,I5,4F16.3,F17.3)") i_mo, ' (', occ_vir, ') ', ikp_for_print, &
3357 1186 : eigenval_scf(i_mo)*evolt, &
3358 1186 : Sigma_c_ikp_n_qp(i_mo)*evolt, &
3359 1186 : Sigma_x_ikp_n(i_mo)*evolt, &
3360 1186 : V_xc_ikp_n(i_mo)*evolt, &
3361 2477 : bs_env%eigenval_G0W0(i_mo, ikp, ispin)*evolt
3362 : END DO
3363 :
3364 105 : WRITE (iunit, "(A)") " "
3365 :
3366 105 : CALL close_file(iunit)
3367 :
3368 : END IF
3369 :
3370 242 : CALL timestop(handle)
3371 :
3372 484 : END SUBROUTINE analyt_conti_and_print
3373 :
3374 : ! **************************************************************************************************
3375 : !> \brief ...
3376 : !> \param Sigma_c_ikp_n_qp ...
3377 : !> \param ispin ...
3378 : !> \param bs_env ...
3379 : ! **************************************************************************************************
3380 242 : SUBROUTINE correct_obvious_fitting_fails(Sigma_c_ikp_n_qp, ispin, bs_env)
3381 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: Sigma_c_ikp_n_qp
3382 : INTEGER :: ispin
3383 : TYPE(post_scf_bandstructure_type), POINTER :: bs_env
3384 :
3385 : CHARACTER(LEN=*), PARAMETER :: routineN = 'correct_obvious_fitting_fails'
3386 :
3387 : INTEGER :: handle, homo, i_mo, j_mo, &
3388 : n_levels_scissor, n_mo
3389 : LOGICAL :: is_occ, is_vir
3390 : REAL(KIND=dp) :: sum_Sigma_c
3391 :
3392 242 : CALL timeset(routineN, handle)
3393 :
3394 242 : n_mo = bs_env%n_ao
3395 242 : homo = bs_env%n_occ(ispin)
3396 :
3397 2902 : DO i_mo = 1, n_mo
3398 :
3399 : ! if |𝚺^c| > 13 eV, we use a scissors shift
3400 2902 : IF (ABS(Sigma_c_ikp_n_qp(i_mo)) > 13.0_dp/evolt) THEN
3401 :
3402 0 : is_occ = (i_mo <= homo)
3403 0 : is_vir = (i_mo > homo)
3404 :
3405 0 : n_levels_scissor = 0
3406 0 : sum_Sigma_c = 0.0_dp
3407 :
3408 : ! compute scissor
3409 0 : DO j_mo = 1, n_mo
3410 :
3411 : ! only compute scissor from other GW levels close in energy
3412 0 : IF (is_occ .AND. j_mo > homo) CYCLE
3413 0 : IF (is_vir .AND. j_mo <= homo) CYCLE
3414 0 : IF (ABS(i_mo - j_mo) > 10) CYCLE
3415 0 : IF (i_mo == j_mo) CYCLE
3416 :
3417 0 : n_levels_scissor = n_levels_scissor + 1
3418 0 : sum_Sigma_c = sum_Sigma_c + Sigma_c_ikp_n_qp(j_mo)
3419 :
3420 : END DO
3421 :
3422 : ! overwrite the self-energy with scissor shift
3423 0 : Sigma_c_ikp_n_qp(i_mo) = sum_Sigma_c/REAL(n_levels_scissor, KIND=dp)
3424 :
3425 : END IF
3426 :
3427 : END DO ! i_mo
3428 :
3429 242 : CALL timestop(handle)
3430 :
3431 242 : END SUBROUTINE correct_obvious_fitting_fails
3432 :
3433 : END MODULE gw_utils
|