(git:33f85d8)
Loading...
Searching...
No Matches
mp2_gpw.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Calls routines to get RI integrals and calculate total energies
10!> \par History
11!> 10.2011 created [Joost VandeVondele and Mauro Del Ben]
12!> 07.2019 split from mp2_gpw.F [Frederick Stein]
13! **************************************************************************************************
14MODULE mp2_gpw
19 USE cell_types, ONLY: cell_type,&
26 USE cp_dbcsr_api, ONLY: &
27 dbcsr_clear_mempools, dbcsr_copy, dbcsr_create, dbcsr_distribution_release, &
30 dbcsr_p_type, dbcsr_release, dbcsr_type, dbcsr_type_no_symmetry, dbcsr_type_symmetric
35 USE cp_fm_types, ONLY: cp_fm_get_info,&
38 USE cp_log_handling, ONLY: &
42 USE dbt_api, ONLY: dbt_type
53 USE hfx_types, ONLY: block_ind_type,&
55 USE input_constants, ONLY: &
59 USE kinds, ONLY: dp
60 USE kpoint_types, ONLY: kpoint_type
63 USE machine, ONLY: default_output_unit,&
74 USE mp2_types, ONLY: mp2_type,&
82 USE qs_kind_types, ONLY: get_qs_kind,&
85 USE qs_mo_types, ONLY: get_mo_set,&
95 USE rpa_rse, ONLY: rse_energy
96#include "./base/base_uses.f90"
97
98 IMPLICIT NONE
99
100 PRIVATE
101
102 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'mp2_gpw'
103
105
106CONTAINS
107
108! **************************************************************************************************
109!> \brief with a big bang to mp2
110!> \param qs_env ...
111!> \param mp2_env ...
112!> \param Emp2 ...
113!> \param Emp2_Cou ...
114!> \param Emp2_EX ...
115!> \param Emp2_S ...
116!> \param Emp2_T ...
117!> \param mos_mp2 ...
118!> \param para_env ...
119!> \param unit_nr ...
120!> \param calc_forces ...
121!> \param calc_ex ...
122!> \param do_ri_mp2 ...
123!> \param do_ri_rpa ...
124!> \param do_ri_sos_laplace_mp2 ...
125!> \author Mauro Del Ben and Joost VandeVondele
126! **************************************************************************************************
127 SUBROUTINE mp2_gpw_main(qs_env, mp2_env, Emp2, Emp2_Cou, Emp2_EX, Emp2_S, Emp2_T, &
128 mos_mp2, para_env, unit_nr, calc_forces, calc_ex, do_ri_mp2, do_ri_rpa, &
129 do_ri_sos_laplace_mp2)
130 TYPE(qs_environment_type), POINTER :: qs_env
131 TYPE(mp2_type) :: mp2_env
132 REAL(kind=dp), INTENT(OUT) :: emp2, emp2_cou, emp2_ex, emp2_s, emp2_t
133 TYPE(mo_set_type), DIMENSION(:), INTENT(IN) :: mos_mp2
134 TYPE(mp_para_env_type), POINTER :: para_env
135 INTEGER, INTENT(IN) :: unit_nr
136 LOGICAL, INTENT(IN) :: calc_forces, calc_ex
137 LOGICAL, INTENT(IN), OPTIONAL :: do_ri_mp2, do_ri_rpa, &
138 do_ri_sos_laplace_mp2
139
140 CHARACTER(LEN=*), PARAMETER :: routinen = 'mp2_gpw_main'
141
142 INTEGER :: blacs_grid_layout, bse_lev_virt, color_sub, dimen_ri, dimen_ri_red, eri_method, &
143 handle, ispin, local_unit_nr, my_group_l_end, my_group_l_size, my_group_l_start, nmo, &
144 nspins, potential_type, ri_metric_type
145 INTEGER, ALLOCATABLE, DIMENSION(:) :: ends_array_mc, ends_array_mc_block, gw_corr_lev_occ, &
146 gw_corr_lev_virt, homo, starts_array_mc, starts_array_mc_block
147 INTEGER, DIMENSION(3) :: periodic
148 LOGICAL :: blacs_repeatable, do_bse, do_im_time, do_kpoints_cubic_rpa, my_do_gw, &
149 my_do_ri_mp2, my_do_ri_rpa, my_do_ri_sos_laplace_mp2
150 REAL(kind=dp) :: emp2_ab, emp2_bb, emp2_cou_bb, &
151 emp2_ex_bb, eps_gvg_rspace_old, &
152 eps_pgf_orb_old, eps_rho_rspace_old
153 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: eigenval
154 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: bib_c_bse_ab, bib_c_bse_ij
155 REAL(kind=dp), DIMENSION(:), POINTER :: mo_eigenvalues
156 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
157 TYPE(block_ind_type), ALLOCATABLE, &
158 DIMENSION(:, :, :) :: t_3c_o_ind
159 TYPE(cell_type), POINTER :: cell
160 TYPE(cp_blacs_env_type), POINTER :: blacs_env_sub, blacs_env_sub_mat_munu
161 TYPE(cp_fm_type) :: fm_matrix_pq
162 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: mo_coeff
163 TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:, :) :: fm_matrix_l_kpoints, fm_matrix_minv, &
164 fm_matrix_minv_l_kpoints, &
165 fm_matrix_minv_vtrunc_minv
166 TYPE(cp_fm_type), POINTER :: mo_coeff_ptr
167 TYPE(cp_logger_type), POINTER :: logger, logger_sub
168 TYPE(dbcsr_p_type) :: mat_munu, mat_p_global
169 TYPE(dbcsr_p_type), ALLOCATABLE, DIMENSION(:) :: mo_coeff_all, mo_coeff_gw, mo_coeff_o, &
170 mo_coeff_o_bse, mo_coeff_v, &
171 mo_coeff_v_bse
172 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s
173 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s_kp
174 TYPE(dbt_type) :: t_3c_m
175 TYPE(dbt_type), ALLOCATABLE, DIMENSION(:, :) :: t_3c_o
176 TYPE(dft_control_type), POINTER :: dft_control
177 TYPE(group_dist_d1_type) :: gd_array, gd_b_all, gd_b_occ_bse, &
178 gd_b_virt_bse
179 TYPE(group_dist_d1_type), ALLOCATABLE, &
180 DIMENSION(:) :: gd_b_virtual
181 TYPE(hfx_compression_type), ALLOCATABLE, &
182 DIMENSION(:, :, :) :: t_3c_o_compressed
183 TYPE(kpoint_type), POINTER :: kpoints, kpoints_from_dft
184 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
185 TYPE(mp_para_env_type), POINTER :: para_env_sub
186 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
187 POINTER :: sab_orb_sub
188 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
189 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
190 TYPE(qs_ks_env_type), POINTER :: ks_env
191 TYPE(three_dim_real_array), ALLOCATABLE, &
192 DIMENSION(:) :: bib_c, bib_c_gw
193
194 CALL timeset(routinen, handle)
195
196 ! check if we want to do ri-mp2
197 my_do_ri_mp2 = .false.
198 IF (PRESENT(do_ri_mp2)) my_do_ri_mp2 = do_ri_mp2
199
200 ! check if we want to do ri-rpa
201 my_do_ri_rpa = .false.
202 IF (PRESENT(do_ri_rpa)) my_do_ri_rpa = do_ri_rpa
203
204 ! check if we want to do ri-sos-laplace-mp2
205 my_do_ri_sos_laplace_mp2 = .false.
206 IF (PRESENT(do_ri_sos_laplace_mp2)) my_do_ri_sos_laplace_mp2 = do_ri_sos_laplace_mp2
207
208 ! GW and SOS-MP2 cannot be used together
209 IF (my_do_ri_sos_laplace_mp2) THEN
210 cpassert(.NOT. mp2_env%ri_rpa%do_ri_g0w0)
211 END IF
212
213 ! check if we want to do imaginary time
214 do_im_time = mp2_env%do_im_time
215 do_bse = qs_env%mp2_env%bse%do_bse
216 do_kpoints_cubic_rpa = qs_env%mp2_env%ri_rpa_im_time%do_im_time_kpoints
217
218 IF (do_kpoints_cubic_rpa .AND. mp2_env%ri_rpa%do_ri_g0w0) THEN
219 cpabort("Full RPA k-points (DO_KPOINTS in LOW_SCALING section) not implemented with GW")
220 END IF
221
222 ! Get the number of spins
223 nspins = SIZE(mos_mp2)
224
225 ! ... setup needed to be able to qs_integrate in a subgroup.
226 IF (do_kpoints_cubic_rpa) THEN
227 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, kpoints=kpoints_from_dft)
228 mos(1:nspins) => kpoints_from_dft%kp_env(1)%kpoint_env%mos(1:nspins, 1)
229 ELSE
230 CALL get_qs_env(qs_env=qs_env, dft_control=dft_control, mos=mos)
231 END IF
232 CALL get_mo_set(mo_set=mos_mp2(1), nmo=nmo)
233 ALLOCATE (homo(nspins), eigenval(nmo, nspins), mo_coeff(nspins))
234 DO ispin = 1, nspins
235 CALL get_mo_set(mo_set=mos_mp2(ispin), &
236 eigenvalues=mo_eigenvalues, homo=homo(ispin), &
237 mo_coeff=mo_coeff_ptr)
238 mo_coeff(ispin) = mo_coeff_ptr
239 eigenval(:, ispin) = mo_eigenvalues(1:nmo)
240 END DO
241
242 ! a para_env
243 color_sub = para_env%mepos/mp2_env%mp2_num_proc
244 ALLOCATE (para_env_sub)
245 CALL para_env_sub%from_split(para_env, color_sub)
246
247 ! each of the sub groups might need to generate output
248 logger => cp_get_default_logger()
249 IF (para_env%is_source()) THEN
250 local_unit_nr = cp_logger_get_default_unit_nr(logger, local=.false.)
251 ELSE
252 local_unit_nr = default_output_unit
253 END IF
254
255 ! get stuff
256 CALL get_qs_env(qs_env, &
257 ks_env=ks_env, &
258 qs_kind_set=qs_kind_set, &
259 cell=cell, &
260 particle_set=particle_set, &
261 atomic_kind_set=atomic_kind_set, &
262 dft_control=dft_control, &
263 matrix_s_kp=matrix_s_kp)
264
265 CALL get_cell(cell=cell, periodic=periodic)
266
267 IF (do_im_time) THEN
268 IF (mp2_env%ri_metric%potential_type == ri_default) THEN
269 IF (sum(periodic) == 1 .OR. sum(periodic) == 3) THEN
270 mp2_env%ri_metric%potential_type = do_potential_id
271 ELSE
272 mp2_env%ri_metric%potential_type = do_potential_truncated
273 END IF
274 END IF
275
276 ! statically initialize libint
278
279 END IF
280
281 IF (mp2_env%ri_metric%potential_type == ri_default) THEN
282 mp2_env%ri_metric%potential_type = do_potential_coulomb
283 END IF
284
285 IF (mp2_env%eri_method == eri_default) THEN
286 IF (sum(periodic) > 0) mp2_env%eri_method = do_eri_gpw
287 IF (sum(periodic) == 0) mp2_env%eri_method = do_eri_os
288 IF (sum(mp2_env%ri_rpa_im_time%kp_grid) > 0) mp2_env%eri_method = do_eri_os
289 IF (mp2_env%method == mp2_method_gpw) mp2_env%eri_method = do_eri_gpw
290 IF (mp2_env%method == ri_mp2_method_gpw) mp2_env%eri_method = do_eri_gpw
291 IF (mp2_env%ri_rpa_im_time%do_im_time_kpoints) mp2_env%eri_method = do_eri_os
292 IF (calc_forces .AND. mp2_env%eri_method == do_eri_os) mp2_env%eri_method = do_eri_gpw
293 END IF
294 eri_method = mp2_env%eri_method
295
296 IF (unit_nr > 0 .AND. mp2_env%eri_method == do_eri_gpw) THEN
297 WRITE (unit=unit_nr, fmt="(T3,A,T71,F10.1)") &
298 "GPW_INFO| Density cutoff [a.u.]:", mp2_env%mp2_gpw%cutoff*0.5_dp
299 WRITE (unit=unit_nr, fmt="(T3,A,T71,F10.1)") &
300 "GPW_INFO| Relative density cutoff [a.u.]:", mp2_env%mp2_gpw%relative_cutoff*0.5_dp
301 CALL m_flush(unit_nr)
302 END IF
303
304 ! MG: Disable logger layer for BSE, misses some key information to print cube files properly
305 IF (.NOT. (mp2_env%ri_g0w0%print_local_bandgap .OR. mp2_env%bse%do_nto_analysis)) THEN
306 ! a logger
307 NULLIFY (logger_sub)
308 CALL cp_logger_create(logger_sub, para_env=para_env_sub, &
309 default_global_unit_nr=local_unit_nr, &
310 close_global_unit_on_dealloc=.false.)
311 CALL cp_logger_set(logger_sub, local_filename="MP2_localLog")
312 ! set to a custom print level (we could also have a different print level for para_env%source)
313 logger_sub%iter_info%print_level = mp2_env%mp2_gpw%print_level
314 CALL cp_add_default_logger(logger_sub)
315 END IF
316
317 ! a blacs_env (ignore the globenv stored defaults for now)
318 blacs_grid_layout = blacs_grid_square
319 blacs_repeatable = .true.
320 NULLIFY (blacs_env_sub)
321 CALL cp_blacs_env_create(blacs_env_sub, para_env_sub, &
322 blacs_grid_layout, &
323 blacs_repeatable)
324
325 blacs_env_sub_mat_munu => blacs_env_sub
326
327 matrix_s(1:1) => matrix_s_kp(1:1, 1)
328
329 CALL get_eps_old(dft_control, eps_pgf_orb_old, eps_rho_rspace_old, eps_gvg_rspace_old)
330
331 CALL create_mat_munu(mat_munu, qs_env, mp2_env%mp2_gpw%eps_grid, &
332 blacs_env_sub_mat_munu, do_alloc_blocks_from_nbl=.NOT. do_im_time, sab_orb_sub=sab_orb_sub, &
333 do_kpoints=mp2_env%ri_rpa_im_time%do_im_time_kpoints, &
334 dbcsr_sym_type=dbcsr_type_symmetric)
335
336 ! which RI metric we want to have
337 ri_metric_type = mp2_env%ri_metric%potential_type
338
339 ! which interaction potential
340 potential_type = mp2_env%potential_parameter%potential_type
341
342 ! check if we want to do ri-g0w0 on top of ri-rpa
343 my_do_gw = mp2_env%ri_rpa%do_ri_g0w0
344 ALLOCATE (gw_corr_lev_occ(nspins), gw_corr_lev_virt(nspins))
345 gw_corr_lev_occ(1) = mp2_env%ri_g0w0%corr_mos_occ
346 gw_corr_lev_virt(1) = mp2_env%ri_g0w0%corr_mos_virt
347 IF (nspins == 2) THEN
348 gw_corr_lev_occ(2) = mp2_env%ri_g0w0%corr_mos_occ_beta
349 gw_corr_lev_virt(2) = mp2_env%ri_g0w0%corr_mos_virt_beta
350 END IF
351
352 IF (do_bse) THEN
353 IF (nspins > 1) THEN
354 cpabort("BSE not implemented for open shell calculations")
355 END IF
356 !Keep default behavior for occupied
357 ! We do not implement an explicit bse_lev_occ here, because the small number of occupied levels
358 ! does not critically influence the memory
359 bse_lev_virt = gw_corr_lev_virt(1)
360 END IF
361
362 ! After the components are inside of the routines, we can move this line insight the branch
363 ALLOCATE (mo_coeff_o(nspins), mo_coeff_v(nspins), mo_coeff_all(nspins), mo_coeff_gw(nspins))
364
365 ! Always allocate for usage in call of replicate_mat_to_subgroup
366 ALLOCATE (mo_coeff_o_bse(1), mo_coeff_v_bse(1))
367
368 ! for imag. time, we do not need this
369 IF (.NOT. do_im_time) THEN
370
371 ! new routine: replicate a full matrix from one para_env to a smaller one
372 ! keeping the memory usage as small as possible in this case the
373 ! output the two part of the C matrix (virtual, occupied)
374 DO ispin = 1, nspins
375
376 CALL replicate_mat_to_subgroup(para_env, para_env_sub, mo_coeff(ispin), nmo, homo(ispin), mat_munu%matrix, &
377 mo_coeff_o(ispin)%matrix, mo_coeff_v(ispin)%matrix, &
378 mo_coeff_all(ispin)%matrix, mo_coeff_gw(ispin)%matrix, &
379 my_do_gw, gw_corr_lev_occ(ispin), gw_corr_lev_virt(ispin), do_bse, &
380 bse_lev_virt, mo_coeff_o_bse(1)%matrix, mo_coeff_v_bse(1)%matrix, &
381 mp2_env%mp2_gpw%eps_filter)
382
383 END DO
384
385 END IF
386
387 ! now we're kind of ready to go....
388 emp2_s = 0.0_dp
389 emp2_t = 0.0_dp
390 IF (my_do_ri_mp2 .OR. my_do_ri_rpa .OR. my_do_ri_sos_laplace_mp2) THEN
391 ! RI-GPW integrals (same stuff for both RPA and MP2)
392 IF (nspins == 2) THEN
393 ! open shell case (RI) here the (ia|K) integrals are computed for both the alpha and beta components
395 bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, gd_array, gd_b_virtual, dimen_ri, dimen_ri_red, qs_env, &
396 para_env, para_env_sub, color_sub, cell, particle_set, &
397 atomic_kind_set, qs_kind_set, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, &
398 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, nmo, homo, mat_munu, sab_orb_sub, &
399 mo_coeff_o, mo_coeff_v, mo_coeff_all, mo_coeff_gw, mo_coeff_o_bse, mo_coeff_v_bse, &
400 mp2_env%mp2_gpw%eps_filter, unit_nr, &
401 mp2_env%mp2_memory, mp2_env%calc_PQ_cond_num, calc_forces, blacs_env_sub, my_do_gw .AND. .NOT. do_im_time, &
402 do_bse, gd_b_all, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, &
403 gw_corr_lev_occ(1), gw_corr_lev_virt(1), &
404 bse_lev_virt, &
405 do_im_time, do_kpoints_cubic_rpa, kpoints, &
406 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
407 mp2_env%ri_metric, &
408 gd_b_occ_bse, gd_b_virt_bse)
409 ELSE
410 ! closed shell case (RI)
411 CALL mp2_ri_gpw_compute_in(bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, gd_array, gd_b_virtual, &
412 dimen_ri, dimen_ri_red, qs_env, para_env, para_env_sub, &
413 color_sub, cell, particle_set, &
414 atomic_kind_set, qs_kind_set, fm_matrix_pq, &
415 fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, &
416 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, nmo, homo, &
417 mat_munu, sab_orb_sub, &
418 mo_coeff_o, mo_coeff_v, mo_coeff_all, mo_coeff_gw, mo_coeff_o_bse, mo_coeff_v_bse, &
419 mp2_env%mp2_gpw%eps_filter, unit_nr, &
420 mp2_env%mp2_memory, mp2_env%calc_PQ_cond_num, calc_forces, &
421 blacs_env_sub, my_do_gw .AND. .NOT. do_im_time, do_bse, gd_b_all, &
422 starts_array_mc, ends_array_mc, &
423 starts_array_mc_block, ends_array_mc_block, &
424 gw_corr_lev_occ(1), gw_corr_lev_virt(1), &
425 bse_lev_virt, &
426 do_im_time, do_kpoints_cubic_rpa, kpoints, &
427 t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
428 mp2_env%ri_metric, gd_b_occ_bse, gd_b_virt_bse)
429 END IF
430
431 ELSE
432 ! Canonical MP2-GPW
433 IF (nspins == 2) THEN
434 ! alpha-alpha and alpha-beta components
435 IF (unit_nr > 0) WRITE (unit_nr, *)
436 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Alpha (ia|'
437 CALL mp2_gpw_compute( &
438 emp2, emp2_cou, emp2_ex, qs_env, para_env, para_env_sub, color_sub, &
439 cell, particle_set, &
440 atomic_kind_set, qs_kind_set, eigenval, nmo, homo, mat_munu, &
441 sab_orb_sub, mo_coeff_o, mo_coeff_v, mp2_env%mp2_gpw%eps_filter, unit_nr, &
442 mp2_env%mp2_memory, calc_ex, blacs_env_sub, emp2_ab)
443
444 ! beta-beta component
445 IF (unit_nr > 0) WRITE (unit_nr, *)
446 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Beta (ia|'
447 CALL mp2_gpw_compute( &
448 emp2_bb, emp2_cou_bb, emp2_ex_bb, qs_env, para_env, para_env_sub, color_sub, cell, particle_set, &
449 atomic_kind_set, qs_kind_set, eigenval(:, 2:2), nmo, homo(2:2), mat_munu, &
450 sab_orb_sub, mo_coeff_o(2:2), mo_coeff_v(2:2), mp2_env%mp2_gpw%eps_filter, unit_nr, &
451 mp2_env%mp2_memory, calc_ex, blacs_env_sub)
452
453 ! make order on the MP2 energy contributions
454 emp2_cou = emp2_cou*0.25_dp
455 emp2_ex = emp2_ex*0.5_dp
456
457 emp2_cou_bb = emp2_cou_bb*0.25_dp
458 emp2_ex_bb = emp2_ex_bb*0.5_dp
459
460 emp2_s = emp2_ab
461 emp2_t = emp2_cou + emp2_cou_bb + emp2_ex + emp2_ex_bb
462
463 emp2_cou = emp2_cou + emp2_cou_bb + emp2_ab
464 emp2_ex = emp2_ex + emp2_ex_bb
465 emp2 = emp2_ex + emp2_cou
466
467 ELSE
468 ! closed shell case
469 CALL mp2_gpw_compute( &
470 emp2, emp2_cou, emp2_ex, qs_env, para_env, para_env_sub, color_sub, cell, particle_set, &
471 atomic_kind_set, qs_kind_set, eigenval(:, 1:1), nmo, homo(1:1), mat_munu, &
472 sab_orb_sub, mo_coeff_o(1:1), mo_coeff_v(1:1), mp2_env%mp2_gpw%eps_filter, unit_nr, &
473 mp2_env%mp2_memory, calc_ex, blacs_env_sub)
474 END IF
475 END IF
476
477 ! Free possibly large buffers allocated by dbcsr on the GPU,
478 ! large hybrid dgemm/pdgemm's coming later will need the space.
479 CALL dbcsr_clear_mempools()
480
481 IF (calc_forces .AND. .NOT. do_im_time) THEN
482 ! make a copy of mo_coeff_o and mo_coeff_v
483 ALLOCATE (mp2_env%ri_grad%mo_coeff_o(nspins), mp2_env%ri_grad%mo_coeff_v(nspins))
484 DO ispin = 1, nspins
485 NULLIFY (mp2_env%ri_grad%mo_coeff_o(ispin)%matrix)
486 CALL dbcsr_init_p(mp2_env%ri_grad%mo_coeff_o(ispin)%matrix)
487 CALL dbcsr_copy(mp2_env%ri_grad%mo_coeff_o(ispin)%matrix, mo_coeff_o(ispin)%matrix, &
488 name="mo_coeff_o"//cp_to_string(ispin))
489 NULLIFY (mp2_env%ri_grad%mo_coeff_v(ispin)%matrix)
490 CALL dbcsr_init_p(mp2_env%ri_grad%mo_coeff_v(ispin)%matrix)
491 CALL dbcsr_copy(mp2_env%ri_grad%mo_coeff_v(ispin)%matrix, mo_coeff_v(ispin)%matrix, &
492 name="mo_coeff_v"//cp_to_string(ispin))
493 END DO
494 CALL get_group_dist(gd_array, color_sub, my_group_l_start, my_group_l_end, my_group_l_size)
495 END IF
496 ! Copy mo coeffs for RPA exchange correction
497 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) THEN
498 ALLOCATE (mp2_env%ri_rpa%mo_coeff_o(nspins), mp2_env%ri_rpa%mo_coeff_v(nspins))
499 DO ispin = 1, nspins
500 CALL dbcsr_copy(mp2_env%ri_rpa%mo_coeff_o(ispin), mo_coeff_o(ispin)%matrix, name="mo_coeff_o")
501 CALL dbcsr_copy(mp2_env%ri_rpa%mo_coeff_v(ispin), mo_coeff_v(ispin)%matrix, name="mo_coeff_v")
502 END DO
503 END IF
504
505 IF (.NOT. do_im_time) THEN
506
507 DO ispin = 1, nspins
508 CALL dbcsr_release(mo_coeff_o(ispin)%matrix)
509 DEALLOCATE (mo_coeff_o(ispin)%matrix)
510 CALL dbcsr_release(mo_coeff_v(ispin)%matrix)
511 DEALLOCATE (mo_coeff_v(ispin)%matrix)
512 IF (my_do_gw) THEN
513 CALL dbcsr_release(mo_coeff_all(ispin)%matrix)
514 DEALLOCATE (mo_coeff_all(ispin)%matrix)
515 END IF
516 END DO
517 DEALLOCATE (mo_coeff_o, mo_coeff_v)
518 IF (my_do_gw) DEALLOCATE (mo_coeff_all)
519
520 END IF
521 IF (do_bse) THEN
522 CALL dbcsr_release(mo_coeff_o_bse(1)%matrix)
523 CALL dbcsr_release(mo_coeff_v_bse(1)%matrix)
524 DEALLOCATE (mo_coeff_o_bse(1)%matrix)
525 DEALLOCATE (mo_coeff_v_bse(1)%matrix)
526 END IF
527 DEALLOCATE (mo_coeff_o_bse, mo_coeff_v_bse)
528
529 ! Release some memory for RPA exchange correction
530 IF (calc_forces .AND. do_im_time .OR. &
531 (.NOT. calc_forces .AND. mp2_env%ri_rpa%exchange_correction == rpa_exchange_none)) THEN
532
533 CALL dbcsr_release(mat_munu%matrix)
534 DEALLOCATE (mat_munu%matrix)
535
536 CALL release_neighbor_list_sets(sab_orb_sub)
537
538 END IF
539
540 ! decide if to do RI-RPA or RI-MP2
541 IF (my_do_ri_rpa .OR. my_do_ri_sos_laplace_mp2) THEN
542
543 IF (do_im_time) CALL create_matrix_p(mat_p_global, qs_env, mp2_env, para_env)
544
545 IF (.NOT. ALLOCATED(bib_c)) ALLOCATE (bib_c(nspins))
546 IF (.NOT. ALLOCATED(bib_c_gw)) ALLOCATE (bib_c_gw(nspins))
547 IF (.NOT. ALLOCATED(gd_b_virtual)) ALLOCATE (gd_b_virtual(nspins))
548
549 ! RI-RPA
550 CALL rpa_ri_compute_en(qs_env, emp2, mp2_env, bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, &
551 para_env, para_env_sub, color_sub, &
552 gd_array, gd_b_virtual, gd_b_all, gd_b_occ_bse, gd_b_virt_bse, &
553 mo_coeff, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, &
554 fm_matrix_minv, fm_matrix_minv_vtrunc_minv, kpoints, &
555 eigenval, nmo, homo, dimen_ri, dimen_ri_red, gw_corr_lev_occ, gw_corr_lev_virt, &
556 bse_lev_virt, &
557 unit_nr, my_do_ri_sos_laplace_mp2, my_do_gw, do_im_time, do_bse, matrix_s, &
558 mat_munu, mat_p_global, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, &
559 starts_array_mc, ends_array_mc, &
560 starts_array_mc_block, ends_array_mc_block, calc_forces)
561
562 IF (mp2_env%ri_rpa%do_rse) &
563 CALL rse_energy(qs_env, mp2_env, para_env, dft_control, mo_coeff, nmo, homo, eigenval)
564
565 IF (do_im_time) THEN
566 IF (ASSOCIATED(mat_p_global%matrix)) THEN
567 CALL dbcsr_release(mat_p_global%matrix)
568 DEALLOCATE (mat_p_global%matrix)
569 END IF
570
572 IF (calc_forces) CALL cp_fm_release(fm_matrix_pq)
573 END IF
574
575 ! Release some memory for RPA exchange correction
576 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) THEN
577
578 CALL dbcsr_release(mat_munu%matrix)
579 DEALLOCATE (mat_munu%matrix)
580
581 CALL release_neighbor_list_sets(sab_orb_sub)
582
583 END IF
584
585 ELSE
586 IF (my_do_ri_mp2) THEN
587 emp2 = 0.0_dp
588 emp2_cou = 0.0_dp
589 emp2_ex = 0.0_dp
590
591 ! RI-MP2-GPW compute energy
593 emp2_cou, emp2_ex, emp2_s, emp2_t, bib_c, mp2_env, para_env, para_env_sub, color_sub, &
594 gd_array, gd_b_virtual, &
595 eigenval, nmo, homo, dimen_ri_red, unit_nr, calc_forces, calc_ex)
596
597 END IF
598 END IF
599
600 ! if we need forces time to calculate the MP2 non-separable contribution
601 ! and start computing the Lagrangian
602 IF (calc_forces .AND. .NOT. do_im_time) THEN
603
604 CALL calc_ri_mp2_nonsep(qs_env, mp2_env, para_env, para_env_sub, cell, &
605 particle_set, atomic_kind_set, qs_kind_set, &
606 mo_coeff, nmo, homo, dimen_ri, eigenval, &
607 my_group_l_start, my_group_l_end, my_group_l_size, &
608 sab_orb_sub, mat_munu, blacs_env_sub)
609
610 DO ispin = 1, nspins
611 CALL dbcsr_release(mp2_env%ri_grad%mo_coeff_o(ispin)%matrix)
612 DEALLOCATE (mp2_env%ri_grad%mo_coeff_o(ispin)%matrix)
613
614 CALL dbcsr_release(mp2_env%ri_grad%mo_coeff_v(ispin)%matrix)
615 DEALLOCATE (mp2_env%ri_grad%mo_coeff_v(ispin)%matrix)
616 END DO
617 DEALLOCATE (mp2_env%ri_grad%mo_coeff_o, mp2_env%ri_grad%mo_coeff_v)
618
619 CALL dbcsr_release(mat_munu%matrix)
620 DEALLOCATE (mat_munu%matrix)
621
622 CALL release_neighbor_list_sets(sab_orb_sub)
623
624 END IF
625
626 !XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXx
627 ! moved from above
628 IF (my_do_gw .AND. .NOT. do_im_time) THEN
629 DO ispin = 1, nspins
630 CALL dbcsr_release(mo_coeff_gw(ispin)%matrix)
631 DEALLOCATE (mo_coeff_gw(ispin)%matrix)
632 END DO
633 DEALLOCATE (mo_coeff_gw)
634 END IF
635
636 ! re-init the radii to be able to generate pair lists with MP2-appropriate screening
637 dft_control%qs_control%eps_pgf_orb = eps_pgf_orb_old
638 dft_control%qs_control%eps_rho_rspace = eps_rho_rspace_old
639 dft_control%qs_control%eps_gvg_rspace = eps_gvg_rspace_old
640 CALL init_interaction_radii(dft_control%qs_control, qs_kind_set)
641
642 CALL cp_blacs_env_release(blacs_env_sub)
643
644 IF (.NOT. (mp2_env%ri_g0w0%print_local_bandgap .OR. mp2_env%bse%do_nto_analysis)) THEN
646 CALL cp_logger_release(logger_sub)
647 END IF
648
649 CALL mp_para_env_release(para_env_sub)
650
651 ! finally solve the z-vector equation if forces are required
652 IF (calc_forces .AND. .NOT. do_im_time) THEN
653 CALL solve_z_vector_eq(qs_env, mp2_env, para_env, dft_control, &
654 mo_coeff, nmo, homo, eigenval, unit_nr)
655 END IF
656
657 DEALLOCATE (eigenval, mo_coeff)
658
659 CALL timestop(handle)
660
661 END SUBROUTINE mp2_gpw_main
662
663! **************************************************************************************************
664!> \brief ...
665!> \param para_env ...
666!> \param para_env_sub ...
667!> \param mo_coeff ...
668!> \param dimen ...
669!> \param homo ...
670!> \param mat_munu ...
671!> \param mo_coeff_o ...
672!> \param mo_coeff_v ...
673!> \param mo_coeff_all ...
674!> \param mo_coeff_gw ...
675!> \param my_do_gw ...
676!> \param gw_corr_lev_occ ...
677!> \param gw_corr_lev_virt ...
678!> \param my_do_bse ...
679!> \param bse_lev_virt ...
680!> \param mo_coeff_o_bse ...
681!> \param mo_coeff_v_bse ...
682!> \param eps_filter ...
683! **************************************************************************************************
684 SUBROUTINE replicate_mat_to_subgroup(para_env, para_env_sub, mo_coeff, dimen, homo, mat_munu, &
685 mo_coeff_o, mo_coeff_v, mo_coeff_all, mo_coeff_gw, my_do_gw, &
686 gw_corr_lev_occ, gw_corr_lev_virt, my_do_bse, &
687 bse_lev_virt, mo_coeff_o_bse, mo_coeff_v_bse, eps_filter)
688 TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_sub
689 TYPE(cp_fm_type), INTENT(IN) :: mo_coeff
690 INTEGER, INTENT(IN) :: dimen, homo
691 TYPE(dbcsr_type), INTENT(INOUT) :: mat_munu
692 TYPE(dbcsr_type), POINTER :: mo_coeff_o, mo_coeff_v, mo_coeff_all, &
693 mo_coeff_gw
694 LOGICAL, INTENT(IN) :: my_do_gw
695 INTEGER, INTENT(IN) :: gw_corr_lev_occ, gw_corr_lev_virt
696 LOGICAL, INTENT(IN) :: my_do_bse
697 INTEGER, INTENT(IN) :: bse_lev_virt
698 TYPE(dbcsr_type), POINTER :: mo_coeff_o_bse, mo_coeff_v_bse
699 REAL(kind=dp), INTENT(IN) :: eps_filter
700
701 CHARACTER(LEN=*), PARAMETER :: routinen = 'replicate_mat_to_subgroup'
702
703 INTEGER :: handle
704 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: c
705 TYPE(group_dist_d1_type) :: gd_array
706
707 CALL timeset(routinen, handle)
708
709 CALL grep_rows_in_subgroups(para_env, para_env_sub, mo_coeff, gd_array, c)
710
711 ! create and fill mo_coeff_o, mo_coeff_v and mo_coeff_all
712 ALLOCATE (mo_coeff_o)
713 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_o, c(:, 1:homo), &
714 mat_munu, gd_array, eps_filter)
715
716 ALLOCATE (mo_coeff_v)
717 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_v, c(:, homo + 1:dimen), &
718 mat_munu, gd_array, eps_filter)
719
720 IF (my_do_gw) THEN
721 ALLOCATE (mo_coeff_gw)
722 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_gw, c(:, homo - gw_corr_lev_occ + 1:homo + gw_corr_lev_virt), &
723 mat_munu, gd_array, eps_filter)
724
725 ! all levels
726 ALLOCATE (mo_coeff_all)
727 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_all, c, &
728 mat_munu, gd_array, eps_filter)
729
730 END IF
731
732 IF (my_do_bse) THEN
733
734 ALLOCATE (mo_coeff_o_bse)
735 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_o_bse, c(:, 1:homo), &
736 mat_munu, gd_array, eps_filter)
737
738 ALLOCATE (mo_coeff_v_bse)
739 CALL build_dbcsr_from_rows(para_env_sub, mo_coeff_v_bse, c(:, homo + 1:homo + bse_lev_virt), &
740 mat_munu, gd_array, eps_filter)
741
742 END IF
743 DEALLOCATE (c)
744 CALL release_group_dist(gd_array)
745
746 CALL timestop(handle)
747
748 END SUBROUTINE replicate_mat_to_subgroup
749
750! **************************************************************************************************
751!> \brief ...
752!> \param para_env ...
753!> \param para_env_sub ...
754!> \param mo_coeff ...
755!> \param gd_array ...
756!> \param C ...
757! **************************************************************************************************
758 SUBROUTINE grep_rows_in_subgroups(para_env, para_env_sub, mo_coeff, gd_array, C)
759 TYPE(mp_para_env_type), INTENT(IN) :: para_env, para_env_sub
760 TYPE(cp_fm_type), INTENT(IN) :: mo_coeff
761 TYPE(group_dist_d1_type), INTENT(OUT) :: gd_array
762 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :), &
763 INTENT(OUT) :: c
764
765 CHARACTER(LEN=*), PARAMETER :: routinen = 'grep_rows_in_subgroups'
766
767 INTEGER :: handle, i_global, iib, j_global, jjb, max_row_col_local, my_mu_end, my_mu_size, &
768 my_mu_start, ncol_global, ncol_local, ncol_rec, nrow_global, nrow_local, nrow_rec, &
769 proc_receive_static, proc_send_static, proc_shift
770 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: local_col_row_info, rec_col_row_info
771 INTEGER, DIMENSION(:), POINTER :: col_indices, col_indices_rec, &
772 row_indices, row_indices_rec
773 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: local_c, rec_c
774 REAL(kind=dp), CONTIGUOUS, DIMENSION(:, :), &
775 POINTER :: local_c_internal
776
777 CALL timeset(routinen, handle)
778
779 CALL cp_fm_get_info(matrix=mo_coeff, &
780 ncol_global=ncol_global, &
781 nrow_global=nrow_global, &
782 nrow_local=nrow_local, &
783 ncol_local=ncol_local, &
784 row_indices=row_indices, &
785 col_indices=col_indices, &
786 local_data=local_c_internal)
787
788 CALL create_group_dist(gd_array, para_env_sub%num_pe, nrow_global)
789 CALL get_group_dist(gd_array, para_env_sub%mepos, my_mu_start, my_mu_end, my_mu_size)
790
791 ! local storage for the C matrix
792 ALLOCATE (c(my_mu_size, ncol_global))
793 c = 0.0_dp
794
795 ALLOCATE (local_c(nrow_local, ncol_local))
796 local_c(:, :) = local_c_internal(1:nrow_local, 1:ncol_local)
797 NULLIFY (local_c_internal)
798
799 max_row_col_local = max(nrow_local, ncol_local)
800 CALL para_env%max(max_row_col_local)
801
802 ALLOCATE (local_col_row_info(0:max_row_col_local, 2))
803 local_col_row_info = 0
804 ! 0,1 nrows
805 local_col_row_info(0, 1) = nrow_local
806 local_col_row_info(1:nrow_local, 1) = row_indices(1:nrow_local)
807 ! 0,2 ncols
808 local_col_row_info(0, 2) = ncol_local
809 local_col_row_info(1:ncol_local, 2) = col_indices(1:ncol_local)
810
811 ALLOCATE (rec_col_row_info(0:max_row_col_local, 2))
812
813 ! accumulate data on C buffer starting from myself
814 DO iib = 1, nrow_local
815 i_global = row_indices(iib)
816 IF (i_global >= my_mu_start .AND. i_global <= my_mu_end) THEN
817 DO jjb = 1, ncol_local
818 j_global = col_indices(jjb)
819 c(i_global - my_mu_start + 1, j_global) = local_c(iib, jjb)
820 END DO
821 END IF
822 END DO
823
824 ! start ring communication for collecting the data from the other
825 proc_send_static = modulo(para_env%mepos + 1, para_env%num_pe)
826 proc_receive_static = modulo(para_env%mepos - 1, para_env%num_pe)
827 DO proc_shift = 1, para_env%num_pe - 1
828 ! first exchange information on the local data
829 rec_col_row_info = 0
830 CALL para_env%sendrecv(local_col_row_info, proc_send_static, rec_col_row_info, proc_receive_static)
831 nrow_rec = rec_col_row_info(0, 1)
832 ncol_rec = rec_col_row_info(0, 2)
833
834 ALLOCATE (row_indices_rec(nrow_rec))
835 row_indices_rec = rec_col_row_info(1:nrow_rec, 1)
836
837 ALLOCATE (col_indices_rec(ncol_rec))
838 col_indices_rec = rec_col_row_info(1:ncol_rec, 2)
839
840 ALLOCATE (rec_c(nrow_rec, ncol_rec))
841 rec_c = 0.0_dp
842
843 ! then send and receive the real data
844 CALL para_env%sendrecv(local_c, proc_send_static, rec_c, proc_receive_static)
845
846 ! accumulate the received data on C buffer
847 DO iib = 1, nrow_rec
848 i_global = row_indices_rec(iib)
849 IF (i_global >= my_mu_start .AND. i_global <= my_mu_end) THEN
850 DO jjb = 1, ncol_rec
851 j_global = col_indices_rec(jjb)
852 c(i_global - my_mu_start + 1, j_global) = rec_c(iib, jjb)
853 END DO
854 END IF
855 END DO
856
857 local_col_row_info(:, :) = rec_col_row_info
858 DEALLOCATE (local_c)
859 ALLOCATE (local_c(nrow_rec, ncol_rec))
860 local_c(:, :) = rec_c
861
862 DEALLOCATE (col_indices_rec)
863 DEALLOCATE (row_indices_rec)
864 DEALLOCATE (rec_c)
865 END DO
866
867 DEALLOCATE (local_c)
868 DEALLOCATE (local_col_row_info)
869 DEALLOCATE (rec_col_row_info)
870
871 CALL timestop(handle)
872
873 END SUBROUTINE grep_rows_in_subgroups
874
875! **************************************************************************************************
876!> \brief Encapsulate the building of dbcsr_matrices mo_coeff_(v,o,all)
877!> \param para_env_sub ...
878!> \param mo_coeff_to_build ...
879!> \param Cread ...
880!> \param mat_munu ...
881!> \param gd_array ...
882!> \param eps_filter ...
883!> \author Jan Wilhelm, Code by Mauro Del Ben
884! **************************************************************************************************
885 SUBROUTINE build_dbcsr_from_rows(para_env_sub, mo_coeff_to_build, Cread, &
886 mat_munu, gd_array, eps_filter)
887 TYPE(mp_para_env_type), INTENT(IN) :: para_env_sub
888 TYPE(dbcsr_type) :: mo_coeff_to_build
889 REAL(kind=dp), DIMENSION(:, :), INTENT(IN) :: cread
890 TYPE(dbcsr_type), INTENT(INOUT) :: mat_munu
891 TYPE(group_dist_d1_type), INTENT(IN) :: gd_array
892 REAL(kind=dp), INTENT(IN) :: eps_filter
893
894 CHARACTER(LEN=*), PARAMETER :: routinen = 'build_dbcsr_from_rows'
895
896 INTEGER :: col, col_offset, col_size, handle, i, i_global, j, j_global, my_mu_end, &
897 my_mu_start, ncol_global, proc_receive, proc_send, proc_shift, rec_mu_end, rec_mu_size, &
898 rec_mu_start, row, row_offset, row_size
899 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: rec_c
900 REAL(kind=dp), DIMENSION(:, :), POINTER :: data_block
901 TYPE(dbcsr_iterator_type) :: iter
902
903 CALL timeset(routinen, handle)
904
905 ncol_global = SIZE(cread, 2)
906
907 CALL get_group_dist(gd_array, para_env_sub%mepos, my_mu_start, my_mu_end)
908
909 CALL cp_dbcsr_m_by_n_from_row_template(mo_coeff_to_build, template=mat_munu, n=ncol_global, &
910 sym=dbcsr_type_no_symmetry)
911 CALL dbcsr_reserve_all_blocks(mo_coeff_to_build)
912
913 ! accumulate data on mo_coeff_to_build starting from myself
914 CALL dbcsr_iterator_start(iter, mo_coeff_to_build)
915 DO WHILE (dbcsr_iterator_blocks_left(iter))
916 CALL dbcsr_iterator_next_block(iter, row, col, data_block, &
917 row_size=row_size, col_size=col_size, &
918 row_offset=row_offset, col_offset=col_offset)
919 DO i = 1, row_size
920 i_global = row_offset + i - 1
921 IF (i_global >= my_mu_start .AND. i_global <= my_mu_end) THEN
922 DO j = 1, col_size
923 j_global = col_offset + j - 1
924 data_block(i, j) = cread(i_global - my_mu_start + 1, col_offset + j - 1)
925 END DO
926 END IF
927 END DO
928 END DO
929 CALL dbcsr_iterator_stop(iter)
930
931 ! start ring communication in the subgroup for collecting the data from the other
932 ! proc (occupied)
933 DO proc_shift = 1, para_env_sub%num_pe - 1
934 proc_send = modulo(para_env_sub%mepos + proc_shift, para_env_sub%num_pe)
935 proc_receive = modulo(para_env_sub%mepos - proc_shift, para_env_sub%num_pe)
936
937 CALL get_group_dist(gd_array, proc_receive, rec_mu_start, rec_mu_end, rec_mu_size)
938
939 ALLOCATE (rec_c(rec_mu_size, ncol_global))
940 rec_c = 0.0_dp
941
942 ! then send and receive the real data
943 CALL para_env_sub%sendrecv(cread, proc_send, rec_c, proc_receive)
944
945 ! accumulate data on mo_coeff_to_build the data received from proc_rec
946 CALL dbcsr_iterator_start(iter, mo_coeff_to_build)
947 DO WHILE (dbcsr_iterator_blocks_left(iter))
948 CALL dbcsr_iterator_next_block(iter, row, col, data_block, &
949 row_size=row_size, col_size=col_size, &
950 row_offset=row_offset, col_offset=col_offset)
951 DO i = 1, row_size
952 i_global = row_offset + i - 1
953 IF (i_global >= rec_mu_start .AND. i_global <= rec_mu_end) THEN
954 DO j = 1, col_size
955 j_global = col_offset + j - 1
956 data_block(i, j) = rec_c(i_global - rec_mu_start + 1, col_offset + j - 1)
957 END DO
958 END IF
959 END DO
960 END DO
961 CALL dbcsr_iterator_stop(iter)
962
963 DEALLOCATE (rec_c)
964
965 END DO
966 CALL dbcsr_filter(mo_coeff_to_build, eps_filter)
967
968 CALL timestop(handle)
969
970 END SUBROUTINE build_dbcsr_from_rows
971
972! **************************************************************************************************
973!> \brief Encapsulate the building of dbcsr_matrix mat_munu
974!> \param mat_munu ...
975!> \param qs_env ...
976!> \param eps_grid ...
977!> \param blacs_env_sub ...
978!> \param do_ri_aux_basis ...
979!> \param do_mixed_basis ...
980!> \param group_size_prim ...
981!> \param do_alloc_blocks_from_nbl ...
982!> \param do_kpoints ...
983!> \param sab_orb_sub ...
984!> \param dbcsr_sym_type ...
985!> \author Jan Wilhelm, code by Mauro Del Ben
986! **************************************************************************************************
987 SUBROUTINE create_mat_munu(mat_munu, qs_env, eps_grid, blacs_env_sub, &
988 do_ri_aux_basis, do_mixed_basis, group_size_prim, &
989 do_alloc_blocks_from_nbl, do_kpoints, sab_orb_sub, dbcsr_sym_type)
990
991 TYPE(dbcsr_p_type), INTENT(OUT) :: mat_munu
992 TYPE(qs_environment_type), POINTER :: qs_env
993 REAL(kind=dp) :: eps_grid
994 TYPE(cp_blacs_env_type), POINTER :: blacs_env_sub
995 LOGICAL, INTENT(IN), OPTIONAL :: do_ri_aux_basis, do_mixed_basis
996 INTEGER, INTENT(IN), OPTIONAL :: group_size_prim
997 LOGICAL, INTENT(IN), OPTIONAL :: do_alloc_blocks_from_nbl, do_kpoints
998 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
999 OPTIONAL, POINTER :: sab_orb_sub
1000 CHARACTER, OPTIONAL :: dbcsr_sym_type
1001
1002 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_mat_munu'
1003
1004 CHARACTER :: my_dbcsr_sym_type
1005 INTEGER :: handle, ikind, natom, nkind
1006 INTEGER, DIMENSION(:), POINTER :: col_blk_sizes, row_blk_sizes
1007 LOGICAL :: my_do_alloc_blocks_from_nbl, &
1008 my_do_kpoints, my_do_mixed_basis, &
1009 my_do_ri_aux_basis
1010 LOGICAL, ALLOCATABLE, DIMENSION(:) :: orb_present
1011 REAL(dp), ALLOCATABLE, DIMENSION(:) :: orb_radius
1012 REAL(dp), ALLOCATABLE, DIMENSION(:, :) :: pair_radius
1013 REAL(kind=dp) :: subcells
1014 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
1015 TYPE(cell_type), POINTER :: cell
1016 TYPE(dbcsr_distribution_type), POINTER :: dbcsr_dist_sub
1017 TYPE(dft_control_type), POINTER :: dft_control
1018 TYPE(distribution_1d_type), POINTER :: local_molecules_sub, local_particles_sub
1019 TYPE(distribution_2d_type), POINTER :: distribution_2d_sub
1020 TYPE(gto_basis_set_p_type), DIMENSION(:), POINTER :: basis_set_ri_aux
1021 TYPE(gto_basis_set_type), POINTER :: orb_basis_set
1022 TYPE(local_atoms_type), ALLOCATABLE, DIMENSION(:) :: atom2d
1023 TYPE(molecule_kind_type), DIMENSION(:), POINTER :: molecule_kind_set
1024 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
1025 TYPE(neighbor_list_set_p_type), DIMENSION(:), &
1026 POINTER :: my_sab_orb_sub
1027 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
1028 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
1029
1030 CALL timeset(routinen, handle)
1031
1032 NULLIFY (basis_set_ri_aux)
1033
1034 my_do_ri_aux_basis = .false.
1035 IF (PRESENT(do_ri_aux_basis)) THEN
1036 my_do_ri_aux_basis = do_ri_aux_basis
1037 END IF
1038
1039 my_do_mixed_basis = .false.
1040 IF (PRESENT(do_mixed_basis)) THEN
1041 my_do_mixed_basis = do_mixed_basis
1042 END IF
1043
1044 my_do_alloc_blocks_from_nbl = .false.
1045 IF (PRESENT(do_alloc_blocks_from_nbl)) THEN
1046 my_do_alloc_blocks_from_nbl = do_alloc_blocks_from_nbl
1047 END IF
1048
1049 my_do_kpoints = .false.
1050 IF (PRESENT(do_kpoints)) THEN
1051 my_do_kpoints = do_kpoints
1052 END IF
1053
1054 my_dbcsr_sym_type = dbcsr_type_no_symmetry
1055 IF (PRESENT(dbcsr_sym_type)) THEN
1056 my_dbcsr_sym_type = dbcsr_sym_type
1057 END IF
1058
1059 CALL get_qs_env(qs_env, &
1060 qs_kind_set=qs_kind_set, &
1061 cell=cell, &
1062 particle_set=particle_set, &
1063 atomic_kind_set=atomic_kind_set, &
1064 molecule_set=molecule_set, &
1065 molecule_kind_set=molecule_kind_set, &
1066 dft_control=dft_control)
1067
1068 IF (my_do_kpoints) THEN
1069 ! please choose EPS_PGF_ORB in QS section smaller than EPS_GRID in WFC_GPW section
1070 IF (eps_grid < dft_control%qs_control%eps_pgf_orb) THEN
1071 eps_grid = dft_control%qs_control%eps_pgf_orb
1072 cpwarn("WFC_GPW%EPS_GRID has been set to QS%EPS_PGF_ORB")
1073 END IF
1074 END IF
1075
1076 ! hack hack hack XXXXXXXXXXXXXXX ... to be fixed
1077 dft_control%qs_control%eps_pgf_orb = eps_grid
1078 dft_control%qs_control%eps_rho_rspace = eps_grid
1079 dft_control%qs_control%eps_gvg_rspace = eps_grid
1080 CALL init_interaction_radii(dft_control%qs_control, qs_kind_set)
1081
1082 ! get a distribution_1d
1083 NULLIFY (local_particles_sub, local_molecules_sub)
1084 CALL distribute_molecules_1d(atomic_kind_set=atomic_kind_set, &
1085 particle_set=particle_set, &
1086 local_particles=local_particles_sub, &
1087 molecule_kind_set=molecule_kind_set, &
1088 molecule_set=molecule_set, &
1089 local_molecules=local_molecules_sub, &
1090 force_env_section=qs_env%input)
1091
1092 ! get a distribution_2d
1093 NULLIFY (distribution_2d_sub)
1094 CALL distribute_molecules_2d(cell=cell, &
1095 atomic_kind_set=atomic_kind_set, &
1096 qs_kind_set=qs_kind_set, &
1097 particle_set=particle_set, &
1098 molecule_kind_set=molecule_kind_set, &
1099 molecule_set=molecule_set, &
1100 distribution_2d=distribution_2d_sub, &
1101 blacs_env=blacs_env_sub, &
1102 force_env_section=qs_env%input)
1103
1104 ! Build the sub orbital-orbital overlap neighbor lists
1105 CALL section_vals_val_get(qs_env%input, "DFT%SUBCELLS", r_val=subcells)
1106 nkind = SIZE(atomic_kind_set)
1107 ALLOCATE (atom2d(nkind))
1108
1109 CALL atom2d_build(atom2d, local_particles_sub, distribution_2d_sub, atomic_kind_set, &
1110 molecule_set, molecule_only=.false., particle_set=particle_set)
1111
1112 ALLOCATE (orb_present(nkind))
1113 ALLOCATE (orb_radius(nkind))
1114 ALLOCATE (pair_radius(nkind, nkind))
1115
1116 DO ikind = 1, nkind
1117 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis_set)
1118 IF (ASSOCIATED(orb_basis_set)) THEN
1119 orb_present(ikind) = .true.
1120 CALL get_gto_basis_set(gto_basis_set=orb_basis_set, kind_radius=orb_radius(ikind))
1121 ELSE
1122 orb_present(ikind) = .false.
1123 orb_radius(ikind) = 0.0_dp
1124 END IF
1125 END DO
1126
1127 CALL pair_radius_setup(orb_present, orb_present, orb_radius, orb_radius, pair_radius)
1128
1129 IF (PRESENT(sab_orb_sub)) THEN
1130 NULLIFY (sab_orb_sub)
1131 ! for cubic RPA/GW with kpoints, we need all neighbors and not only the symmetric ones
1132 IF (my_do_kpoints) THEN
1133 CALL build_neighbor_lists(sab_orb_sub, particle_set, atom2d, cell, pair_radius, &
1134 mic=.false., subcells=subcells, molecular=.false., nlname="sab_orb_sub", &
1135 symmetric=.false.)
1136 ELSE
1137 CALL build_neighbor_lists(sab_orb_sub, particle_set, atom2d, cell, pair_radius, &
1138 mic=.false., subcells=subcells, molecular=.false., nlname="sab_orb_sub")
1139 END IF
1140 ELSE
1141 NULLIFY (my_sab_orb_sub)
1142 ! for cubic RPA/GW with kpoints, we need all neighbors and not only the symmetric ones
1143 IF (my_do_kpoints) THEN
1144 CALL build_neighbor_lists(my_sab_orb_sub, particle_set, atom2d, cell, pair_radius, &
1145 mic=.false., subcells=subcells, molecular=.false., nlname="sab_orb_sub", &
1146 symmetric=.false.)
1147 ELSE
1148 CALL build_neighbor_lists(my_sab_orb_sub, particle_set, atom2d, cell, pair_radius, &
1149 mic=.false., subcells=subcells, molecular=.false., nlname="sab_orb_sub")
1150 END IF
1151 END IF
1152 CALL atom2d_cleanup(atom2d)
1153 DEALLOCATE (atom2d)
1154 DEALLOCATE (orb_present, orb_radius, pair_radius)
1155
1156 ! a dbcsr_dist
1157 ALLOCATE (dbcsr_dist_sub)
1158 CALL cp_dbcsr_dist2d_to_dist(distribution_2d_sub, dbcsr_dist_sub)
1159
1160 ! build a dbcsr matrix the hard way
1161 natom = SIZE(particle_set)
1162 ALLOCATE (row_blk_sizes(natom))
1163 IF (my_do_ri_aux_basis) THEN
1164
1165 ALLOCATE (basis_set_ri_aux(nkind))
1166 CALL basis_set_list_setup(basis_set_ri_aux, "RI_AUX", qs_kind_set)
1167 CALL get_particle_set(particle_set, qs_kind_set, nsgf=row_blk_sizes, basis=basis_set_ri_aux)
1168 DEALLOCATE (basis_set_ri_aux)
1169
1170 ELSE IF (my_do_mixed_basis) THEN
1171
1172 ALLOCATE (basis_set_ri_aux(nkind))
1173 CALL basis_set_list_setup(basis_set_ri_aux, "RI_AUX", qs_kind_set)
1174 CALL get_particle_set(particle_set, qs_kind_set, nsgf=row_blk_sizes, basis=basis_set_ri_aux)
1175 DEALLOCATE (basis_set_ri_aux)
1176
1177 ALLOCATE (col_blk_sizes(natom))
1178
1179 CALL get_particle_set(particle_set, qs_kind_set, nsgf=col_blk_sizes)
1180 col_blk_sizes = col_blk_sizes*group_size_prim
1181
1182 ELSE
1183 CALL get_particle_set(particle_set, qs_kind_set, nsgf=row_blk_sizes)
1184 END IF
1185
1186 NULLIFY (mat_munu%matrix)
1187 ALLOCATE (mat_munu%matrix)
1188
1189 IF (my_do_ri_aux_basis) THEN
1190
1191 CALL dbcsr_create(matrix=mat_munu%matrix, &
1192 name="(ai|munu)", &
1193 dist=dbcsr_dist_sub, matrix_type=my_dbcsr_sym_type, &
1194 row_blk_size=row_blk_sizes, col_blk_size=row_blk_sizes)
1195
1196 ELSE IF (my_do_mixed_basis) THEN
1197
1198 CALL dbcsr_create(matrix=mat_munu%matrix, &
1199 name="(ai|munu)", &
1200 dist=dbcsr_dist_sub, matrix_type=my_dbcsr_sym_type, &
1201 row_blk_size=row_blk_sizes, col_blk_size=col_blk_sizes)
1202
1203 ELSE
1204
1205 CALL dbcsr_create(matrix=mat_munu%matrix, &
1206 name="(ai|munu)", &
1207 dist=dbcsr_dist_sub, matrix_type=my_dbcsr_sym_type, &
1208 row_blk_size=row_blk_sizes, col_blk_size=row_blk_sizes)
1209
1210 IF (my_do_alloc_blocks_from_nbl) THEN
1211
1212 IF (PRESENT(sab_orb_sub)) THEN
1213 CALL cp_dbcsr_alloc_block_from_nbl(mat_munu%matrix, sab_orb_sub)
1214 ELSE
1215 CALL cp_dbcsr_alloc_block_from_nbl(mat_munu%matrix, my_sab_orb_sub)
1216 END IF
1217
1218 END IF
1219
1220 END IF
1221
1222 DEALLOCATE (row_blk_sizes)
1223
1224 IF (my_do_mixed_basis) THEN
1225 DEALLOCATE (col_blk_sizes)
1226 END IF
1227
1228 CALL dbcsr_distribution_release(dbcsr_dist_sub)
1229 DEALLOCATE (dbcsr_dist_sub)
1230
1231 CALL distribution_2d_release(distribution_2d_sub)
1232
1233 CALL distribution_1d_release(local_particles_sub)
1234 CALL distribution_1d_release(local_molecules_sub)
1235
1236 IF (.NOT. PRESENT(sab_orb_sub)) THEN
1237 CALL release_neighbor_list_sets(my_sab_orb_sub)
1238 END IF
1239
1240 CALL timestop(handle)
1241
1242 END SUBROUTINE create_mat_munu
1243
1244! **************************************************************************************************
1245!> \brief ...
1246!> \param mat_P_global ...
1247!> \param qs_env ...
1248!> \param mp2_env ...
1249!> \param para_env ...
1250! **************************************************************************************************
1251 SUBROUTINE create_matrix_p(mat_P_global, qs_env, mp2_env, para_env)
1252
1253 TYPE(dbcsr_p_type), INTENT(OUT) :: mat_p_global
1254 TYPE(qs_environment_type), POINTER :: qs_env
1255 TYPE(mp2_type) :: mp2_env
1256 TYPE(mp_para_env_type), POINTER :: para_env
1257
1258 CHARACTER(LEN=*), PARAMETER :: routinen = 'create_matrix_P'
1259
1260 INTEGER :: blacs_grid_layout, handle
1261 LOGICAL :: blacs_repeatable
1262 TYPE(cp_blacs_env_type), POINTER :: blacs_env_global
1263
1264 CALL timeset(routinen, handle)
1265
1266 blacs_grid_layout = blacs_grid_square
1267 blacs_repeatable = .true.
1268 NULLIFY (blacs_env_global)
1269 CALL cp_blacs_env_create(blacs_env_global, para_env, &
1270 blacs_grid_layout, &
1271 blacs_repeatable)
1272
1273 CALL create_mat_munu(mat_p_global, qs_env, mp2_env%mp2_gpw%eps_grid, &
1274 blacs_env_global, do_ri_aux_basis=.true., &
1275 do_kpoints=mp2_env%ri_rpa_im_time%do_im_time_kpoints)
1276
1277 CALL dbcsr_reserve_all_blocks(mat_p_global%matrix)
1278 CALL cp_blacs_env_release(blacs_env_global)
1279
1280 CALL timestop(handle)
1281
1282 END SUBROUTINE
1283
1284! **************************************************************************************************
1285!> \brief ...
1286!> \param dft_control ...
1287!> \param eps_pgf_orb_old ...
1288!> \param eps_rho_rspace_old ...
1289!> \param eps_gvg_rspace_old ...
1290! **************************************************************************************************
1291 PURE SUBROUTINE get_eps_old(dft_control, eps_pgf_orb_old, eps_rho_rspace_old, eps_gvg_rspace_old)
1292
1293 TYPE(dft_control_type), INTENT(INOUT) :: dft_control
1294 REAL(kind=dp), INTENT(OUT) :: eps_pgf_orb_old, eps_rho_rspace_old, &
1295 eps_gvg_rspace_old
1296
1297 ! re-init the radii to be able to generate pair lists with MP2-appropriate screening
1298 eps_pgf_orb_old = dft_control%qs_control%eps_pgf_orb
1299 eps_rho_rspace_old = dft_control%qs_control%eps_rho_rspace
1300 eps_gvg_rspace_old = dft_control%qs_control%eps_gvg_rspace
1301
1302 END SUBROUTINE get_eps_old
1303
1304END MODULE mp2_gpw
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Define the atomic kind types and their sub types.
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
...
Handles all functions related to the CELL.
Definition cell_types.F:15
subroutine, public get_cell(cell, alpha, beta, gamma, deth, orthorhombic, abc, periodic, h, h_inv, symmetry_id, tag)
Get informations about a simulation cell.
Definition cell_types.F:195
methods related to the blacs parallel environment
integer, parameter, public blacs_grid_square
subroutine, public cp_blacs_env_release(blacs_env)
releases the given blacs_env
subroutine, public cp_blacs_env_create(blacs_env, para_env, blacs_grid_layout, blacs_repeatable, row_major, grid_2d)
allocates and initializes a type that represent a blacs context
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_distribution_release(dist)
...
subroutine, public dbcsr_iterator_next_block(iterator, row, column, block, block_number_argument_has_been_removed, row_size, col_size, row_offset, col_offset)
...
logical function, public dbcsr_iterator_blocks_left(iterator)
...
subroutine, public dbcsr_iterator_stop(iterator)
...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_init_p(matrix)
...
subroutine, public dbcsr_filter(matrix, eps)
...
subroutine, public dbcsr_iterator_start(iterator, matrix, shared, dynamic, dynamic_byrows)
...
subroutine, public dbcsr_release(matrix)
...
subroutine, public dbcsr_reserve_all_blocks(matrix)
Reserves all blocks.
DBCSR operations in CP2K.
subroutine, public cp_dbcsr_dist2d_to_dist(dist2d, dist)
Creates a DBCSR distribution from a distribution_2d.
subroutine, public cp_dbcsr_m_by_n_from_row_template(matrix, template, n, sym)
Utility function to create dbcsr matrix, m x n matrix (n arbitrary) with the same processor grid and ...
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_get_info(matrix, name, nrow_global, ncol_global, nrow_block, ncol_block, nrow_local, ncol_local, row_indices, col_indices, local_data, context, nrow_locals, ncol_locals, matrix_struct, para_env)
returns all kind of information about the full matrix
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
subroutine, public cp_logger_set(logger, local_filename, global_filename)
sets various attributes of the given logger
subroutine, public cp_rm_default_logger()
the cousin of cp_add_default_logger, decrements the stack, so that the default logger is what it has ...
subroutine, public cp_logger_release(logger)
releases this logger
subroutine, public cp_logger_create(logger, para_env, print_level, default_global_unit_nr, default_local_unit_nr, global_filename, local_filename, close_global_unit_on_dealloc, iter_info, close_local_unit_on_dealloc, suffix, template_logger)
initializes a logger
subroutine, public cp_add_default_logger(logger)
adds a default logger. MUST be called before logging occours
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
This is the start of a dbt_api, all publically needed functions are exported here....
Definition dbt_api.F:17
stores a lists of integer that are local to a processor. The idea is that these integers represent ob...
subroutine, public distribution_1d_release(distribution_1d)
releases the given distribution_1d
stores a mapping of 2D info (e.g. matrix) on a 2D processor distribution (i.e. blacs grid) where cpus...
subroutine, public distribution_2d_release(distribution_2d)
...
Distribution methods for atoms, particles, or molecules.
subroutine, public distribute_molecules_1d(atomic_kind_set, particle_set, local_particles, molecule_kind_set, molecule_set, local_molecules, force_env_section, prev_molecule_kind_set, prev_local_molecules)
Distribute molecules and particles.
subroutine, public distribute_molecules_2d(cell, atomic_kind_set, particle_set, qs_kind_set, molecule_kind_set, molecule_set, distribution_2d, blacs_env, force_env_section)
Distributes the particle pairs creating a 2d distribution optimally suited for quickstep.
Types to describe group distributions.
Types and set/get functions for HFX.
Definition hfx_types.F:15
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public rpa_exchange_none
integer, parameter, public do_eri_os
integer, parameter, public ri_mp2_method_gpw
integer, parameter, public do_potential_truncated
integer, parameter, public mp2_method_gpw
integer, parameter, public do_potential_id
integer, parameter, public eri_default
integer, parameter, public ri_default
integer, parameter, public do_potential_coulomb
integer, parameter, public do_eri_gpw
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Types and basic routines needed for a kpoint calculation.
Interface to the Libint-Library or a c++ wrapper.
subroutine, public cp_libint_static_cleanup()
subroutine, public cp_libint_static_init()
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
integer, parameter, public default_output_unit
Definition machine.F:53
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:130
Interface to the message passing library MPI.
subroutine, public mp_para_env_release(para_env)
releases the para object (to be called when you don't want anymore the shared copy of this object)
Define the molecule kind structure types and the corresponding functionality.
Define the data structure for the molecule information.
Routines to calculate CPHF like update and solve Z-vector equation for MP2 gradients (only GPW)
Definition mp2_cphf.F:14
subroutine, public solve_z_vector_eq(qs_env, mp2_env, para_env, dft_control, mo_coeff, nmo, homo, eigenval, unit_nr)
Solve Z-vector equations necessary for the calculation of the MP2 gradients, in order to be consisten...
Definition mp2_cphf.F:161
Routines to calculate MP2 energy using GPW method.
subroutine, public mp2_gpw_compute(emp2, emp2_cou, emp2_ex, qs_env, para_env, para_env_sub, color_sub, cell, particle_set, atomic_kind_set, qs_kind_set, eigenval, nmo, homo, mat_munu, sab_orb_sub, mo_coeff_o, mo_coeff_v, eps_filter, unit_nr, mp2_memory, calc_ex, blacs_env_sub, emp2_ab)
...
Calls routines to get RI integrals and calculate total energies.
Definition mp2_gpw.F:14
subroutine, public mp2_gpw_main(qs_env, mp2_env, emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, mos_mp2, para_env, unit_nr, calc_forces, calc_ex, do_ri_mp2, do_ri_rpa, do_ri_sos_laplace_mp2)
with a big bang to mp2
Definition mp2_gpw.F:130
subroutine, public grep_rows_in_subgroups(para_env, para_env_sub, mo_coeff, gd_array, c)
...
Definition mp2_gpw.F:759
subroutine, public build_dbcsr_from_rows(para_env_sub, mo_coeff_to_build, cread, mat_munu, gd_array, eps_filter)
Encapsulate the building of dbcsr_matrices mo_coeff_(v,o,all)
Definition mp2_gpw.F:887
subroutine, public create_mat_munu(mat_munu, qs_env, eps_grid, blacs_env_sub, do_ri_aux_basis, do_mixed_basis, group_size_prim, do_alloc_blocks_from_nbl, do_kpoints, sab_orb_sub, dbcsr_sym_type)
Encapsulate the building of dbcsr_matrix mat_munu.
Definition mp2_gpw.F:990
Routines to calculate and distribute 2c- and 3c- integrals for RI.
subroutine, public mp2_ri_gpw_compute_in(bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, gd_array, gd_b_virtual, dimen_ri, dimen_ri_red, qs_env, para_env, para_env_sub, color_sub, cell, particle_set, atomic_kind_set, qs_kind_set, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, nmo, homo, mat_munu, sab_orb_sub, mo_coeff_o, mo_coeff_v, mo_coeff_all, mo_coeff_gw, mo_coeff_o_bse, mo_coeff_v_bse, eps_filter, unit_nr, mp2_memory, calc_pq_cond_num, calc_forces, blacs_env_sub, my_do_gw, do_bse, gd_b_all, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, gw_corr_lev_occ, gw_corr_lev_virt, bse_lev_virt, do_im_time, do_kpoints_cubic_rpa, kpoints, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, ri_metric, gd_b_occ_bse, gd_b_virt_bse)
with ri mp2 gpw
Routines to calculate RI-GPW-MP2 energy using pw.
Definition mp2_ri_gpw.F:14
subroutine, public mp2_ri_gpw_compute_en(emp2_cou, emp2_ex, emp2_s, emp2_t, bib_c, mp2_env, para_env, para_env_sub, color_sub, gd_array, gd_b_virtual, eigenval, nmo, homo, dimen_ri, unit_nr, calc_forces, calc_ex)
...
Definition mp2_ri_gpw.F:75
Routines to calculate gradients of RI-GPW-MP2 energy using pw.
Definition mp2_ri_grad.F:13
subroutine, public calc_ri_mp2_nonsep(qs_env, mp2_env, para_env, para_env_sub, cell, particle_set, atomic_kind_set, qs_kind_set, mo_coeff, nmo, homo, dimen_ri, eigenval, my_group_l_start, my_group_l_end, my_group_l_size, sab_orb_sub, mat_munu, blacs_env_sub)
Calculate the non-separable part of the gradients and update the Lagrangian.
Types needed for MP2 calculations.
Definition mp2_types.F:14
Define methods related to particle_type.
subroutine, public get_particle_set(particle_set, qs_kind_set, first_sgf, last_sgf, nsgf, nmao, basis)
Get the components of a particle set.
Define the data structure for the particle information.
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Get the QUICKSTEP environment.
Some utility functions for the calculation of integrals.
subroutine, public basis_set_list_setup(basis_set_list, basis_type, qs_kind_set)
Set up an easy accessible list of the basis sets for all kinds.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public init_interaction_radii(qs_control, qs_kind_set)
Initialize all the atomic kind radii for a given threshold value.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public get_mo_set(mo_set, maxocc, homo, lfomo, nao, nelectron, n_el_f, nmo, eigenvalues, occupation_numbers, mo_coeff, mo_coeff_b, uniform_occupation, kts, mu, flexible_electron_count)
Get the components of a MO set data structure.
Define the neighbor list data types and the corresponding functionality.
subroutine, public release_neighbor_list_sets(nlists)
releases an array of neighbor_list_sets
Generate the atomic neighbor lists.
subroutine, public atom2d_cleanup(atom2d)
free the internals of atom2d
subroutine, public pair_radius_setup(present_a, present_b, radius_a, radius_b, pair_radius, prmin)
...
subroutine, public build_neighbor_lists(ab_list, particle_set, atom, cell, pair_radius, subcells, mic, symmetric, molecular, subset_of_mol, current_subset, operator_type, nlname, atomb_to_keep)
Build simple pair neighbor lists.
subroutine, public atom2d_build(atom2d, distribution_1d, distribution_2d, atomic_kind_set, molecule_set, molecule_only, particle_set)
Build some distribution structure of atoms, refactored from build_qs_neighbor_lists.
Routines to calculate RI-RPA energy.
Definition rpa_main.F:17
subroutine, public rpa_ri_compute_en(qs_env, erpa, mp2_env, bib_c, bib_c_gw, bib_c_bse_ij, bib_c_bse_ab, para_env, para_env_sub, color_sub, gd_array, gd_b_virtual, gd_b_all, gd_b_occ_bse, gd_b_virt_bse, mo_coeff, fm_matrix_pq, fm_matrix_l_kpoints, fm_matrix_minv_l_kpoints, fm_matrix_minv, fm_matrix_minv_vtrunc_minv, kpoints, eigenval, nmo, homo, dimen_ri, dimen_ri_red, gw_corr_lev_occ, gw_corr_lev_virt, bse_lev_virt, unit_nr, do_ri_sos_laplace_mp2, my_do_gw, do_im_time, do_bse, matrix_s, mat_munu, mat_p_global, t_3c_m, t_3c_o, t_3c_o_compressed, t_3c_o_ind, starts_array_mc, ends_array_mc, starts_array_mc_block, ends_array_mc_block, calc_forces)
...
Definition rpa_main.F:195
Routines to compute singles correction to RPA (RSE)
Definition rpa_rse.F:14
subroutine, public rse_energy(qs_env, mp2_env, para_env, dft_control, mo_coeff, nmo, homo, eigenval)
Single excitations energy corrections for RPA.
Definition rpa_rse.F:93
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
structure to store local (to a processor) ordered lists of integers.
distributes pairs on a 2d grid of processors
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
calculation environment to calculate the ks matrix, holds all the needed vars. assumes that the core ...