Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Routines to compute singles correction to RPA (RSE)
10 : !> \par History
11 : !> 08.2019 created [Vladimir Rybkin]
12 : !> \author Vladimir Rybkin
13 : ! **************************************************************************************************
14 : MODULE rpa_rse
15 :
16 : USE cp_blacs_env, ONLY: cp_blacs_env_type
17 : USE cp_control_types, ONLY: dft_control_type
18 : USE cp_dbcsr_api, ONLY: dbcsr_copy,&
19 : dbcsr_create,&
20 : dbcsr_init_p,&
21 : dbcsr_p_type,&
22 : dbcsr_release,&
23 : dbcsr_scale,&
24 : dbcsr_set,&
25 : dbcsr_type_symmetric
26 : USE cp_dbcsr_cp2k_link, ONLY: cp_dbcsr_alloc_block_from_nbl
27 : USE cp_dbcsr_operations, ONLY: copy_dbcsr_to_fm,&
28 : copy_fm_to_dbcsr,&
29 : dbcsr_allocate_matrix_set
30 : USE cp_fm_basic_linalg, ONLY: cp_fm_scale_and_add
31 : USE cp_fm_diag, ONLY: choose_eigv_solver
32 : USE cp_fm_struct, ONLY: cp_fm_struct_create,&
33 : cp_fm_struct_release,&
34 : cp_fm_struct_type
35 : USE cp_fm_types, ONLY: cp_fm_create,&
36 : cp_fm_get_diag,&
37 : cp_fm_get_info,&
38 : cp_fm_release,&
39 : cp_fm_set_all,&
40 : cp_fm_to_fm_submat,&
41 : cp_fm_type
42 : USE hfx_energy_potential, ONLY: integrate_four_center
43 : USE hfx_exx, ONLY: exx_post_hfx,&
44 : exx_pre_hfx
45 : USE hfx_ri, ONLY: hfx_ri_update_ks
46 : USE input_section_types, ONLY: section_vals_get,&
47 : section_vals_get_subs_vals,&
48 : section_vals_type,&
49 : section_vals_val_get
50 : USE kinds, ONLY: dp
51 : USE message_passing, ONLY: mp_para_env_type
52 : USE mp2_types, ONLY: mp2_type
53 : USE parallel_gemm_api, ONLY: parallel_gemm
54 : USE pw_types, ONLY: pw_r3d_rs_type
55 : USE qs_energy_types, ONLY: qs_energy_type
56 : USE qs_environment_types, ONLY: get_qs_env,&
57 : qs_environment_type
58 : USE qs_ks_types, ONLY: qs_ks_env_type
59 : USE qs_ks_utils, ONLY: compute_matrix_vxc
60 : USE qs_neighbor_list_types, ONLY: neighbor_list_set_p_type
61 : USE qs_rho_types, ONLY: qs_rho_get,&
62 : qs_rho_type
63 : USE qs_vxc, ONLY: qs_vxc_create
64 :
65 : !$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
66 :
67 : #include "./base/base_uses.f90"
68 :
69 : IMPLICIT NONE
70 :
71 : PRIVATE
72 :
73 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'rpa_rse'
74 :
75 : PUBLIC :: rse_energy
76 :
77 : CONTAINS
78 :
79 : ! **************************************************************************************************
80 : !> \brief Single excitations energy corrections for RPA
81 : !> \param qs_env ...
82 : !> \param mp2_env ...
83 : !> \param para_env ...
84 : !> \param dft_control ...
85 : !> \param mo_coeff ...
86 : !> \param nmo ...
87 : !> \param homo ...
88 : !> \param Eigenval ...
89 : !> \author Vladimir Rybkin, 08/2019
90 : ! **************************************************************************************************
91 4 : SUBROUTINE rse_energy(qs_env, mp2_env, para_env, dft_control, &
92 4 : mo_coeff, nmo, homo, Eigenval)
93 : TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
94 : TYPE(mp2_type), INTENT(INOUT) :: mp2_env
95 : TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
96 : TYPE(dft_control_type), INTENT(IN), POINTER :: dft_control
97 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
98 : INTEGER, INTENT(IN) :: nmo
99 : INTEGER, DIMENSION(:), INTENT(IN) :: homo
100 : REAL(KIND=dp), DIMENSION(:, :), INTENT(IN) :: Eigenval
101 :
102 : CHARACTER(LEN=*), PARAMETER :: routineN = 'rse_energy'
103 :
104 : INTEGER :: dimen, handle, i_global, iiB, ispin, &
105 : j_global, jjB, n_rep_hf, ncol_local, &
106 : nrow_local, nspins
107 4 : INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
108 : LOGICAL :: do_hfx, hfx_treat_lsd_in_core
109 : REAL(KIND=dp) :: coeff, corr, rse_corr
110 4 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: diag_diff
111 : TYPE(cp_blacs_env_type), POINTER :: blacs_env
112 : TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
113 : TYPE(cp_fm_type) :: fm_ao, fm_ao_mo
114 4 : TYPE(cp_fm_type), ALLOCATABLE, DIMENSION(:) :: fm_P_mu_nu, fm_X_mo, fm_XC_mo
115 4 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: mat_mu_nu, matrix_s, rho_ao
116 : TYPE(neighbor_list_set_p_type), DIMENSION(:), &
117 4 : POINTER :: sab_orb
118 : TYPE(qs_energy_type), POINTER :: energy
119 : TYPE(qs_rho_type), POINTER :: rho
120 : TYPE(section_vals_type), POINTER :: hfx_sections, input
121 :
122 4 : CALL timeset(routineN, handle)
123 :
124 4 : nspins = dft_control%nspins
125 :
126 : ! Pick the diagonal terms
127 : CALL cp_fm_get_info(matrix=mo_coeff(1), &
128 : nrow_local=nrow_local, &
129 : ncol_local=ncol_local, &
130 : row_indices=row_indices, &
131 4 : col_indices=col_indices)
132 :
133 : ! start collecting stuff
134 4 : dimen = nmo
135 4 : NULLIFY (input, matrix_s, blacs_env, rho, energy, sab_orb)
136 : CALL get_qs_env(qs_env, &
137 : input=input, &
138 : matrix_s=matrix_s, &
139 : blacs_env=blacs_env, &
140 : rho=rho, &
141 : energy=energy, &
142 4 : sab_orb=sab_orb)
143 :
144 4 : CALL qs_rho_get(rho, rho_ao=rho_ao)
145 :
146 : ! hfx section
147 4 : NULLIFY (hfx_sections)
148 4 : hfx_sections => section_vals_get_subs_vals(input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
149 4 : CALL section_vals_get(hfx_sections, explicit=do_hfx, n_repetition=n_rep_hf)
150 4 : IF (do_hfx) THEN
151 : CALL section_vals_val_get(hfx_sections, "TREAT_LSD_IN_CORE", l_val=hfx_treat_lsd_in_core, &
152 4 : i_rep_section=1)
153 : END IF
154 :
155 : ! create work array
156 4 : NULLIFY (mat_mu_nu)
157 4 : CALL dbcsr_allocate_matrix_set(mat_mu_nu, nspins)
158 10 : DO ispin = 1, nspins
159 6 : ALLOCATE (mat_mu_nu(ispin)%matrix)
160 : CALL dbcsr_create(matrix=mat_mu_nu(ispin)%matrix, template=matrix_s(1)%matrix, name="T_mu_nu", &
161 6 : matrix_type=dbcsr_type_symmetric, nze=0)
162 6 : CALL cp_dbcsr_alloc_block_from_nbl(mat_mu_nu(ispin)%matrix, sab_orb)
163 10 : CALL dbcsr_set(mat_mu_nu(ispin)%matrix, 0.0_dp)
164 : END DO
165 :
166 : ! Dense (full) matrices
167 18 : ALLOCATE (fm_P_mu_nu(nspins))
168 4 : NULLIFY (fm_struct_tmp)
169 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
170 4 : nrow_global=dimen, ncol_global=dimen)
171 10 : DO ispin = 1, nspins
172 6 : CALL cp_fm_create(fm_P_mu_nu(ispin), fm_struct_tmp, name="P_mu_nu")
173 10 : CALL cp_fm_set_all(fm_P_mu_nu(ispin), 0.0_dp)
174 : END DO
175 4 : CALL cp_fm_struct_release(fm_struct_tmp)
176 :
177 4 : NULLIFY (fm_struct_tmp)
178 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
179 4 : nrow_global=dimen, ncol_global=dimen)
180 32 : ALLOCATE (fm_X_mo(nspins), fm_XC_mo(nspins))
181 10 : DO ispin = 1, nspins
182 6 : CALL cp_fm_create(fm_X_mo(ispin), fm_struct_tmp, name="f_X_mo")
183 6 : CALL cp_fm_create(fm_XC_mo(ispin), fm_struct_tmp, name="f_XC_mo")
184 6 : CALL cp_fm_set_all(fm_X_mo(ispin), 0.0_dp)
185 10 : CALL cp_fm_set_all(fm_XC_mo(ispin), 0.0_dp)
186 : END DO
187 4 : CALL cp_fm_struct_release(fm_struct_tmp)
188 :
189 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
190 4 : nrow_global=dimen, ncol_global=dimen)
191 4 : CALL cp_fm_create(fm_ao, fm_struct_tmp, name="f_ao")
192 4 : CALL cp_fm_struct_release(fm_struct_tmp)
193 4 : CALL cp_fm_set_all(fm_ao, 0.0_dp)
194 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
195 4 : nrow_global=dimen, ncol_global=dimen)
196 4 : CALL cp_fm_create(fm_ao_mo, fm_struct_tmp, name="f_ao_mo")
197 4 : CALL cp_fm_struct_release(fm_struct_tmp)
198 4 : CALL cp_fm_set_all(fm_ao_mo, 0.0_dp)
199 :
200 : !
201 : ! Ready with preparations, do the real staff
202 : !
203 :
204 : ! Obtain density matrix like quantity
205 :
206 4 : coeff = 1.0_dp
207 4 : IF (nspins == 1) coeff = 2.0_dp
208 10 : DO ispin = 1, nspins
209 : CALL parallel_gemm(transa='N', transb='T', m=dimen, n=dimen, k=homo(ispin), alpha=coeff, &
210 : matrix_a=mo_coeff(ispin), matrix_b=mo_coeff(ispin), &
211 10 : beta=0.0_dp, matrix_c=fm_P_mu_nu(ispin))
212 : END DO
213 :
214 : ! Calculate exact exchange contribution
215 : CALL exchange_contribution(qs_env, para_env, dimen, mo_coeff, &
216 : hfx_sections, n_rep_hf, &
217 : rho, mat_mu_nu, fm_P_mu_nu, &
218 4 : fm_ao, fm_X_mo, fm_ao_mo)
219 :
220 : ! Calculate DFT exchange-correlation contribution
221 4 : CALL xc_contribution(qs_env, fm_ao, fm_ao_mo, fm_XC_mo, mo_coeff, dimen)
222 :
223 12 : ALLOCATE (diag_diff(dimen))
224 4 : rse_corr = 0.0_dp
225 :
226 10 : DO ispin = 1, nspins
227 : ! Compute the correction matrix: it is stored in fm_X_mo
228 6 : CALL cp_fm_scale_and_add(1.0_dp, fm_X_mo(ispin), -1.0_dp, fm_XC_mo(ispin))
229 :
230 : ! Pick the diagonal terms
231 6 : CALL cp_fm_get_diag(fm_X_mo(ispin), diag_diff)
232 :
233 : ! Compute the correction
234 : CALL cp_fm_get_info(matrix=fm_X_mo(ispin), &
235 : nrow_local=nrow_local, &
236 : ncol_local=ncol_local, &
237 : row_indices=row_indices, &
238 6 : col_indices=col_indices)
239 :
240 6 : corr = 0.0_dp
241 :
242 : !$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
243 : !$OMP REDUCTION(+: corr) &
244 6 : !$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,diag_diff,eigenval,fm_X_mo,homo,ispin)
245 : DO jjB = 1, ncol_local
246 : j_global = col_indices(jjB)
247 : DO iiB = 1, nrow_local
248 : i_global = row_indices(iiB)
249 : IF ((i_global .LE. homo(ispin)) .AND. (j_global .GT. homo(ispin))) THEN
250 : corr = corr + fm_X_mo(ispin)%local_data(iib, jjb)**2.0_dp/ &
251 : (eigenval(i_global, ispin) - eigenval(j_global, ispin) - diag_diff(i_global) + diag_diff(j_global))
252 : END IF
253 : END DO
254 : END DO
255 : !$OMP END PARALLEL DO
256 :
257 10 : rse_corr = rse_corr + corr
258 : END DO
259 :
260 4 : CALL para_env%sum(rse_corr)
261 :
262 4 : IF (nspins == 1) rse_corr = rse_corr*2.0_dp
263 :
264 4 : mp2_env%ri_rpa%rse_corr_diag = rse_corr
265 :
266 4 : CALL non_diag_rse(fm_X_mo, eigenval, dimen, homo, para_env, blacs_env, rse_corr)
267 :
268 4 : IF (nspins == 1) rse_corr = rse_corr*2.0_dp
269 :
270 4 : mp2_env%ri_rpa%rse_corr = rse_corr
271 :
272 : ! Release staff
273 4 : DEALLOCATE (diag_diff)
274 4 : CALL cp_fm_release(fm_ao)
275 4 : CALL cp_fm_release(fm_ao_mo)
276 4 : CALL cp_fm_release(fm_P_mu_nu)
277 4 : CALL cp_fm_release(fm_X_mo)
278 4 : CALL cp_fm_release(fm_XC_mo)
279 10 : DO ispin = 1, nspins
280 6 : CALL dbcsr_release(mat_mu_nu(ispin)%matrix)
281 10 : DEALLOCATE (mat_mu_nu(ispin)%matrix)
282 : END DO
283 4 : DEALLOCATE (mat_mu_nu)
284 :
285 4 : CALL timestop(handle)
286 :
287 20 : END SUBROUTINE rse_energy
288 :
289 : ! **************************************************************************************************
290 : !> \brief HF exchange occupied-virtual matrix
291 : !> \param qs_env ...
292 : !> \param para_env ...
293 : !> \param dimen ...
294 : !> \param mo_coeff ...
295 : !> \param hfx_sections ...
296 : !> \param n_rep_hf ...
297 : !> \param rho_work ...
298 : !> \param mat_mu_nu ...
299 : !> \param fm_P_mu_nu ...
300 : !> \param fm_X_ao ...
301 : !> \param fm_X_mo ...
302 : !> \param fm_X_ao_mo ...
303 : ! **************************************************************************************************
304 12 : SUBROUTINE exchange_contribution(qs_env, para_env, dimen, mo_coeff, &
305 : hfx_sections, n_rep_hf, &
306 4 : rho_work, mat_mu_nu, fm_P_mu_nu, &
307 4 : fm_X_ao, fm_X_mo, fm_X_ao_mo)
308 : TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
309 : TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
310 : INTEGER, INTENT(IN) :: dimen
311 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: mo_coeff
312 : TYPE(section_vals_type), INTENT(IN), POINTER :: hfx_sections
313 : INTEGER, INTENT(IN) :: n_rep_hf
314 : TYPE(qs_rho_type), INTENT(IN), POINTER :: rho_work
315 : TYPE(dbcsr_p_type), DIMENSION(:), INTENT(IN), &
316 : POINTER :: mat_mu_nu
317 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_P_mu_nu
318 : TYPE(cp_fm_type), INTENT(IN) :: fm_X_ao
319 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_X_mo
320 : TYPE(cp_fm_type), INTENT(IN) :: fm_X_ao_mo
321 :
322 : CHARACTER(LEN=*), PARAMETER :: routineN = 'exchange_contribution'
323 :
324 : INTEGER :: handle, irep, is, ns
325 : LOGICAL :: my_recalc_hfx_integrals
326 : REAL(KIND=dp) :: ehfx
327 4 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: P_mu_nu, rho_work_ao
328 4 : TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: mat_2d, rho_ao_2d
329 :
330 4 : CALL timeset(routineN, handle)
331 :
332 4 : CALL qs_rho_get(rho_work, rho_ao=rho_work_ao)
333 4 : ns = SIZE(rho_work_ao)
334 4 : NULLIFY (P_mu_nu)
335 4 : CALL dbcsr_allocate_matrix_set(P_mu_nu, ns)
336 10 : DO is = 1, ns
337 6 : CALL dbcsr_init_p(P_mu_nu(is)%matrix)
338 6 : CALL dbcsr_create(P_mu_nu(is)%matrix, template=rho_work_ao(1)%matrix)
339 6 : CALL dbcsr_copy(P_mu_nu(is)%matrix, rho_work_ao(1)%matrix)
340 10 : CALL dbcsr_set(P_mu_nu(is)%matrix, 0.0_dp)
341 : END DO
342 :
343 4 : my_recalc_hfx_integrals = .TRUE.
344 :
345 4 : CALL exx_pre_hfx(hfx_sections, qs_env%mp2_env%ri_rpa%x_data, qs_env%mp2_env%ri_rpa%reuse_hfx)
346 10 : DO is = 1, ns
347 6 : CALL copy_fm_to_dbcsr(fm_P_mu_nu(is), P_mu_nu(1)%matrix, keep_sparsity=.TRUE.)
348 :
349 6 : CALL dbcsr_set(mat_mu_nu(1)%matrix, 0.0_dp)
350 :
351 6 : IF (qs_env%mp2_env%ri_rpa%x_data(1, 1)%do_hfx_ri) THEN
352 :
353 0 : DO irep = 1, n_rep_hf
354 0 : rho_ao_2d(1:ns, 1:1) => P_mu_nu(1:ns)
355 0 : mat_2d(1:ns, 1:1) => mat_mu_nu(1:ns)
356 : CALL hfx_ri_update_ks(qs_env, qs_env%mp2_env%ri_rpa%x_data(irep, 1)%ri_data, mat_2d, ehfx, &
357 : rho_ao=rho_ao_2d, geometry_did_change=my_recalc_hfx_integrals, nspins=1, &
358 0 : hf_fraction=qs_env%mp2_env%ri_rpa%x_data(irep, 1)%general_parameter%fraction)
359 :
360 0 : IF (ns == 2) CALL dbcsr_scale(mat_mu_nu(1)%matrix, 2.0_dp)
361 0 : my_recalc_hfx_integrals = .FALSE.
362 : END DO
363 :
364 : ELSE
365 :
366 12 : DO irep = 1, n_rep_hf
367 6 : rho_ao_2d(1:ns, 1:1) => P_mu_nu(1:ns)
368 6 : mat_2d(1:ns, 1:1) => mat_mu_nu(1:ns)
369 : CALL integrate_four_center(qs_env, qs_env%mp2_env%ri_rpa%x_data, mat_2d, ehfx, rho_ao_2d, hfx_sections, &
370 : para_env, my_recalc_hfx_integrals, irep, .TRUE., &
371 6 : ispin=1)
372 :
373 12 : my_recalc_hfx_integrals = .FALSE.
374 : END DO
375 : END IF
376 :
377 : ! copy back to fm
378 6 : CALL cp_fm_set_all(fm_X_ao, 0.0_dp)
379 6 : CALL copy_dbcsr_to_fm(matrix=mat_mu_nu(1)%matrix, fm=fm_X_ao)
380 6 : CALL cp_fm_set_all(fm_X_mo(is), 0.0_dp)
381 :
382 : ! First index
383 : CALL parallel_gemm('T', 'N', dimen, dimen, dimen, 1.0_dp, &
384 6 : mo_coeff(is), fm_X_ao, 0.0_dp, fm_X_ao_mo)
385 :
386 : ! Second index
387 : CALL parallel_gemm('N', 'N', dimen, dimen, dimen, 1.0_dp, &
388 10 : fm_X_ao_mo, mo_coeff(is), 1.0_dp, fm_X_mo(is))
389 :
390 : END DO
391 4 : CALL exx_post_hfx(qs_env, qs_env%mp2_env%ri_rpa%x_data, qs_env%mp2_env%ri_rpa%reuse_hfx)
392 :
393 : ! Release dbcsr objects
394 10 : DO is = 1, SIZE(P_mu_nu)
395 6 : CALL dbcsr_release(P_mu_nu(is)%matrix)
396 10 : DEALLOCATE (P_mu_nu(is)%matrix)
397 : END DO
398 4 : DEALLOCATE (P_mu_nu)
399 :
400 4 : CALL timestop(handle)
401 :
402 4 : END SUBROUTINE exchange_contribution
403 :
404 : ! **************************************************************************************************
405 : !> \brief Exchange-correlation occupied-virtual matrix
406 : !> \param qs_env ...
407 : !> \param fm_XC_ao ...
408 : !> \param fm_XC_ao_mo ...
409 : !> \param fm_XC_mo ...
410 : !> \param mo_coeff ...
411 : !> \param dimen ...
412 : ! **************************************************************************************************
413 4 : SUBROUTINE xc_contribution(qs_env, fm_XC_ao, fm_XC_ao_mo, fm_XC_mo, mo_coeff, dimen)
414 : TYPE(qs_environment_type), INTENT(IN), POINTER :: qs_env
415 : TYPE(cp_fm_type), INTENT(IN) :: fm_XC_ao, fm_XC_ao_mo
416 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_XC_mo, mo_coeff
417 : INTEGER, INTENT(IN) :: dimen
418 :
419 : CHARACTER(LEN=*), PARAMETER :: routineN = 'xc_contribution'
420 :
421 : INTEGER :: handle, i
422 : REAL(KIND=dp) :: exc
423 4 : TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_vxc
424 4 : TYPE(pw_r3d_rs_type), DIMENSION(:), POINTER :: tau_rspace, v_rspace
425 : TYPE(qs_ks_env_type), POINTER :: ks_env
426 : TYPE(qs_rho_type), POINTER :: rho
427 : TYPE(section_vals_type), POINTER :: input, xc_section
428 :
429 4 : CALL timeset(routineN, handle)
430 :
431 4 : NULLIFY (matrix_vxc, v_rspace, tau_rspace, input, xc_section, ks_env, &
432 4 : rho)
433 4 : CALL get_qs_env(qs_env, matrix_vxc=matrix_vxc, input=input, ks_env=ks_env, rho=rho)
434 4 : xc_section => section_vals_get_subs_vals(input, "DFT%XC")
435 :
436 : ! Compute XC matrix in AO basis
437 : CALL qs_vxc_create(ks_env=ks_env, rho_struct=rho, xc_section=xc_section, &
438 4 : vxc_rho=v_rspace, vxc_tau=tau_rspace, exc=exc)
439 :
440 4 : IF (ASSOCIATED(v_rspace)) THEN
441 4 : CALL compute_matrix_vxc(qs_env=qs_env, v_rspace=v_rspace, matrix_vxc=matrix_vxc)
442 :
443 10 : DO i = 1, SIZE(v_rspace)
444 10 : CALL v_rspace(i)%release()
445 : END DO
446 4 : DEALLOCATE (v_rspace)
447 :
448 10 : DO i = 1, SIZE(matrix_vxc)
449 6 : CALL cp_fm_set_all(fm_XC_ao, 0.0_dp)
450 6 : CALL copy_dbcsr_to_fm(matrix=matrix_vxc(i)%matrix, fm=fm_XC_ao)
451 6 : CALL cp_fm_set_all(fm_XC_mo(i), 0.0_dp)
452 :
453 : ! First index
454 : CALL parallel_gemm('T', 'N', dimen, dimen, dimen, 1.0_dp, &
455 6 : mo_coeff(i), fm_XC_ao, 0.0_dp, fm_XC_ao_mo)
456 :
457 : ! Second index
458 : CALL parallel_gemm('N', 'N', dimen, dimen, dimen, 1.0_dp, &
459 10 : fm_XC_ao_mo, mo_coeff(i), 1.0_dp, fm_XC_mo(i))
460 :
461 : END DO
462 :
463 10 : DO i = 1, SIZE(matrix_vxc)
464 6 : CALL dbcsr_release(matrix_vxc(i)%matrix)
465 10 : DEALLOCATE (matrix_vxc(i)%matrix)
466 : END DO
467 4 : DEALLOCATE (matrix_vxc)
468 : END IF
469 :
470 4 : CALL timestop(handle)
471 :
472 4 : END SUBROUTINE xc_contribution
473 :
474 : ! **************************************************************************************************
475 : !> \brief ...
476 : !> \param fm_F_mo ...
477 : !> \param eigenval ...
478 : !> \param dimen ...
479 : !> \param homo ...
480 : !> \param para_env ...
481 : !> \param blacs_env ...
482 : !> \param rse_corr ...
483 : ! **************************************************************************************************
484 4 : SUBROUTINE non_diag_rse(fm_F_mo, eigenval, dimen, homo, para_env, &
485 : blacs_env, rse_corr)
486 : TYPE(cp_fm_type), DIMENSION(:), INTENT(IN) :: fm_F_mo
487 : REAL(KIND=dp), DIMENSION(:, :), INTENT(IN) :: Eigenval
488 : INTEGER, INTENT(IN) :: dimen
489 : INTEGER, DIMENSION(:), INTENT(IN) :: homo
490 : TYPE(mp_para_env_type), INTENT(IN), POINTER :: para_env
491 : TYPE(cp_blacs_env_type), INTENT(IN), POINTER :: blacs_env
492 : REAL(KIND=dp), INTENT(OUT) :: rse_corr
493 :
494 : CHARACTER(LEN=*), PARAMETER :: routineN = 'non_diag_rse'
495 :
496 : INTEGER :: handle, i_global, iiB, ispin, j_global, &
497 : jjB, ncol_local, nrow_local, nspins, &
498 : virtual
499 4 : INTEGER, DIMENSION(:), POINTER :: col_indices, row_indices
500 : REAL(KIND=dp) :: corr
501 4 : REAL(KIND=dp), ALLOCATABLE, DIMENSION(:) :: eig_o, eig_semi_can, eig_v
502 : TYPE(cp_fm_struct_type), POINTER :: fm_struct_tmp
503 : TYPE(cp_fm_type) :: fm_F_oo, fm_F_ov, fm_F_vv, fm_O, fm_tmp, &
504 : fm_U
505 :
506 4 : CALL timeset(routineN, handle)
507 :
508 4 : nspins = SIZE(fm_f_mo)
509 :
510 10 : DO ispin = 1, nspins
511 : ! Add eigenvalues on the diagonal
512 : CALL cp_fm_get_info(matrix=fm_F_mo(ispin), &
513 : nrow_local=nrow_local, &
514 : ncol_local=ncol_local, &
515 : row_indices=row_indices, &
516 6 : col_indices=col_indices)
517 :
518 : !$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
519 10 : !$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,fm_F_mo,eigenval,ispin)
520 : DO jjB = 1, ncol_local
521 : j_global = col_indices(jjB)
522 : DO iiB = 1, nrow_local
523 : i_global = row_indices(iiB)
524 : IF (i_global .EQ. j_global) fm_F_mo(ispin)%local_data(iib, jjb) = &
525 : fm_F_mo(ispin)%local_data(iib, jjb) + eigenval(i_global, ispin)
526 : END DO
527 : END DO
528 : !$OMP END PARALLEL DO
529 : END DO
530 :
531 4 : rse_corr = 0.0_dp
532 :
533 10 : DO ispin = 1, nspins
534 6 : IF (homo(ispin) <= 0 .OR. homo(ispin) >= dimen) CYCLE
535 : ! Create the occupied-occupied and virtual-virtual blocks, eigenvectors
536 6 : NULLIFY (fm_struct_tmp)
537 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
538 6 : nrow_global=homo(ispin), ncol_global=homo(ispin))
539 6 : CALL cp_fm_create(fm_F_oo, fm_struct_tmp, name="F_oo")
540 6 : CALL cp_fm_create(fm_O, fm_struct_tmp, name="O")
541 6 : CALL cp_fm_set_all(fm_F_oo, 0.0_dp)
542 6 : CALL cp_fm_set_all(fm_O, 0.0_dp)
543 6 : CALL cp_fm_struct_release(fm_struct_tmp)
544 :
545 : CALL cp_fm_to_fm_submat(msource=fm_F_mo(ispin), mtarget=fm_F_oo, &
546 : nrow=homo(ispin), ncol=homo(ispin), &
547 : s_firstrow=1, s_firstcol=1, &
548 6 : t_firstrow=1, t_firstcol=1)
549 6 : virtual = dimen - homo(ispin)
550 6 : NULLIFY (fm_struct_tmp)
551 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
552 6 : nrow_global=virtual, ncol_global=virtual)
553 6 : CALL cp_fm_create(fm_F_vv, fm_struct_tmp, name="F_vv")
554 6 : CALL cp_fm_create(fm_U, fm_struct_tmp, name="U")
555 6 : CALL cp_fm_set_all(fm_F_vv, 0.0_dp)
556 6 : CALL cp_fm_set_all(fm_U, 0.0_dp)
557 6 : CALL cp_fm_struct_release(fm_struct_tmp)
558 :
559 : CALL cp_fm_to_fm_submat(msource=fm_F_mo(ispin), mtarget=fm_F_vv, &
560 : nrow=virtual, ncol=virtual, &
561 : s_firstrow=homo(ispin) + 1, s_firstcol=homo(ispin) + 1, &
562 6 : t_firstrow=1, t_firstcol=1)
563 :
564 : ! Diagonalize occupied-occupied and virtual-virtual matrices
565 18 : ALLOCATE (eig_o(homo(ispin)))
566 18 : ALLOCATE (eig_v(virtual))
567 58 : eig_v = 0.0_dp
568 14 : eig_o = 0.0_dp
569 6 : CALL choose_eigv_solver(fm_F_oo, fm_O, eig_o)
570 6 : CALL choose_eigv_solver(fm_F_vv, fm_U, eig_v)
571 :
572 : ! Collect the eigenvalues to one array
573 18 : ALLOCATE (eig_semi_can(dimen))
574 66 : eig_semi_can = 0.0_dp
575 14 : eig_semi_can(1:homo(ispin)) = eig_o(:)
576 58 : eig_semi_can(homo(ispin) + 1:dimen) = eig_v(:)
577 :
578 : ! Create occupied-virtual block
579 6 : NULLIFY (fm_struct_tmp)
580 : CALL cp_fm_struct_create(fm_struct_tmp, para_env=para_env, context=blacs_env, &
581 6 : nrow_global=homo(ispin), ncol_global=virtual)
582 6 : CALL cp_fm_create(fm_F_ov, fm_struct_tmp, name="F_ov")
583 6 : CALL cp_fm_create(fm_tmp, fm_struct_tmp, name="tmp")
584 6 : CALL cp_fm_set_all(fm_F_ov, 0.0_dp)
585 6 : CALL cp_fm_set_all(fm_tmp, 0.0_dp)
586 6 : CALL cp_fm_struct_release(fm_struct_tmp)
587 :
588 : CALL cp_fm_to_fm_submat(msource=fm_F_mo(ispin), mtarget=fm_F_ov, &
589 : nrow=homo(ispin), ncol=virtual, &
590 : s_firstrow=1, s_firstcol=homo(ispin) + 1, &
591 6 : t_firstrow=1, t_firstcol=1)
592 :
593 : CALL parallel_gemm(transa='T', transb='N', m=homo(ispin), n=virtual, k=homo(ispin), alpha=1.0_dp, &
594 6 : matrix_a=fm_O, matrix_b=fm_F_ov, beta=0.0_dp, matrix_c=fm_tmp)
595 :
596 : CALL parallel_gemm(transa='N', transb='N', m=homo(ispin), n=virtual, k=virtual, alpha=1.0_dp, &
597 6 : matrix_a=fm_tmp, matrix_b=fm_U, beta=0.0_dp, matrix_c=fm_F_ov)
598 :
599 : ! Compute the correction
600 : CALL cp_fm_get_info(matrix=fm_F_ov, &
601 : nrow_local=nrow_local, &
602 : ncol_local=ncol_local, &
603 : row_indices=row_indices, &
604 6 : col_indices=col_indices)
605 6 : corr = 0.0_dp
606 : !$OMP PARALLEL DO DEFAULT(NONE) PRIVATE(jjB,iiB,i_global,j_global) &
607 : !$OMP REDUCTION(+:corr) &
608 6 : !$OMP SHARED(ncol_local,nrow_local,col_indices,row_indices,fm_F_ov,eig_semi_can,homo,ispin)
609 : DO jjB = 1, ncol_local
610 : j_global = col_indices(jjB)
611 : DO iiB = 1, nrow_local
612 : i_global = row_indices(iiB)
613 : corr = corr + fm_F_ov%local_data(iib, jjb)**2.0_dp/ &
614 : (eig_semi_can(i_global) - eig_semi_can(j_global + homo(ispin)))
615 : END DO
616 : END DO
617 : !$OMP END PARALLEL DO
618 :
619 6 : rse_corr = rse_corr + corr
620 :
621 : ! Release
622 6 : DEALLOCATE (eig_semi_can)
623 6 : DEALLOCATE (eig_o)
624 6 : DEALLOCATE (eig_v)
625 :
626 6 : CALL cp_fm_release(fm_F_ov)
627 6 : CALL cp_fm_release(fm_F_oo)
628 6 : CALL cp_fm_release(fm_F_vv)
629 6 : CALL cp_fm_release(fm_U)
630 6 : CALL cp_fm_release(fm_O)
631 34 : CALL cp_fm_release(fm_tmp)
632 :
633 : END DO
634 :
635 4 : CALL para_env%sum(rse_corr)
636 :
637 4 : CALL timestop(handle)
638 :
639 8 : END SUBROUTINE non_diag_rse
640 :
641 : END MODULE rpa_rse
|