(git:f664ec6)
Loading...
Searching...
No Matches
mp2.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines to calculate MP2 energy
10!> \par History
11!> 05.2011 created [Mauro Del Ben]
12!> \author Mauro Del Ben
13! **************************************************************************************************
14MODULE mp2
15 USE admm_types, ONLY: admm_type
20 USE bibliography, ONLY: bussy2023,&
24 stein2022,&
25 stein2024,&
26 cite_reference
29 USE cp_dbcsr_api, ONLY: dbcsr_copy,&
44 USE cp_fm_types, ONLY: cp_fm_create,&
55 USE hfx_exx, ONLY: calculate_exx
56 USE hfx_types, ONLY: &
71 USE kinds, ONLY: dp,&
72 int_8
73 USE kpoint_types, ONLY: kpoint_type
74 USE machine, ONLY: m_flush,&
75 m_memory,&
79 USE mp2_gpw, ONLY: mp2_gpw_main
81 USE mp2_types, ONLY: mp2_biel_type,&
85 mp2_type,&
95 USE qs_mo_types, ONLY: allocate_mo_set,&
100 USE qs_scf_methods, ONLY: eigensolver,&
105 USE virial_types, ONLY: virial_type
106
107!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
108
109#include "./base/base_uses.f90"
110
111 IMPLICIT NONE
112
113 PRIVATE
114
115 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'mp2'
116
117 PUBLIC :: mp2_main
118
119CONTAINS
120
121! **************************************************************************************************
122!> \brief the main entry point for MP2 calculations
123!> \param qs_env ...
124!> \param calc_forces ...
125!> \author Mauro Del Ben
126! **************************************************************************************************
127 SUBROUTINE mp2_main(qs_env, calc_forces)
128 TYPE(qs_environment_type), POINTER :: qs_env
129 LOGICAL, INTENT(IN) :: calc_forces
130
131 CHARACTER(len=*), PARAMETER :: routinen = 'mp2_main'
132
133 INTEGER :: bin, cholesky_method, dimen, handle, handle2, i, i_thread, iatom, ii, ikind, &
134 irep, ispin, max_nset, my_bin_size, n_rep_hf, n_threads, nao, natom, ndep, &
135 nfullcols_total, nfullrows_total, nkind, nmo, nspins, unit_nr
136 INTEGER(KIND=int_8) :: mem
137 INTEGER, ALLOCATABLE, DIMENSION(:) :: kind_of, nelec
138 LOGICAL :: calc_ex, do_admm, do_admm_rpa_exx, do_dynamic_load_balancing, do_exx, do_gw, &
139 do_im_time, do_kpoints_cubic_rpa, free_hfx_buffer, reuse_hfx, update_xc_energy
140 REAL(kind=dp) :: e_admm_from_gw(2), e_ex_from_gw, emp2, emp2_aa, emp2_aa_cou, emp2_aa_ex, &
141 emp2_ab, emp2_ab_cou, emp2_ab_ex, emp2_bb, emp2_bb_cou, emp2_bb_ex, emp2_cou, emp2_ex, &
142 emp2_s, emp2_t, maxocc, mem_real, t1, t2, t3
143 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: evals
144 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: auto
145 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :, :) :: c
146 REAL(kind=dp), DIMENSION(:), POINTER :: mo_eigenvalues
147 TYPE(admm_type), POINTER :: admm_env
148 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
149 TYPE(cp_blacs_env_type), POINTER :: blacs_env
150 TYPE(cp_fm_struct_type), POINTER :: fm_struct
151 TYPE(cp_fm_type) :: evecs, fm_matrix_ks, fm_matrix_s, &
152 fm_matrix_work
153 TYPE(cp_fm_type), POINTER :: fm_matrix_ks_red, fm_matrix_s_red, &
154 fm_work_red, mo_coeff
155 TYPE(cp_logger_type), POINTER :: logger
156 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_ks, matrix_s
157 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_ks_transl, matrix_s_kp
158 TYPE(dft_control_type), POINTER :: dft_control
159 TYPE(excited_energy_type), POINTER :: ex_env
160 TYPE(hfx_basis_info_type) :: basis_info
161 TYPE(hfx_basis_type), DIMENSION(:), POINTER :: basis_parameter
162 TYPE(hfx_container_type), DIMENSION(:), POINTER :: integral_containers
163 TYPE(hfx_container_type), POINTER :: maxval_container
164 TYPE(hfx_type), POINTER :: actual_x_data
165 TYPE(kpoint_type), POINTER :: kpoints
166 TYPE(mo_set_type), ALLOCATABLE, DIMENSION(:) :: mos_mp2
167 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
168 TYPE(mp2_biel_type) :: mp2_biel
169 TYPE(mp2_type), POINTER :: mp2_env
170 TYPE(mp_para_env_type), POINTER :: para_env
171 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
172 TYPE(qs_energy_type), POINTER :: energy
173 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
174 TYPE(qs_scf_env_type), POINTER :: scf_env
175 TYPE(scf_control_type), POINTER :: scf_control
176 TYPE(section_vals_type), POINTER :: hfx_sections, input
177 TYPE(virial_type), POINTER :: virial
178
179 ! If SCF has not converged we should abort MP2 calculation
180 IF (qs_env%mp2_env%hf_fail) THEN
181 CALL cp_abort(__location__, "SCF not converged: "// &
182 "not possible to run MP2")
183 END IF
184
185 NULLIFY (virial, dft_control, blacs_env, kpoints, fm_matrix_s_red, fm_matrix_ks_red, fm_work_red)
186 CALL timeset(routinen, handle)
187 logger => cp_get_default_logger()
188
189 CALL cite_reference(delben2012)
190
191 do_kpoints_cubic_rpa = qs_env%mp2_env%ri_rpa_im_time%do_im_time_kpoints
192
193 ! for cubic RPA and GW, we have kpoints and therefore, we get other matrices later
194 IF (do_kpoints_cubic_rpa) THEN
195
196 CALL get_qs_env(qs_env, &
197 input=input, &
198 atomic_kind_set=atomic_kind_set, &
199 qs_kind_set=qs_kind_set, &
200 dft_control=dft_control, &
201 particle_set=particle_set, &
202 para_env=para_env, &
203 blacs_env=blacs_env, &
204 energy=energy, &
205 kpoints=kpoints, &
206 scf_env=scf_env, &
207 scf_control=scf_control, &
208 matrix_ks_kp=matrix_ks_transl, &
209 matrix_s_kp=matrix_s_kp, &
210 mp2_env=mp2_env)
211
212 CALL get_gamma(matrix_s, matrix_ks, mos, &
213 matrix_s_kp, matrix_ks_transl, kpoints)
214
215 ELSE
216
217 CALL get_qs_env(qs_env, &
218 input=input, &
219 atomic_kind_set=atomic_kind_set, &
220 qs_kind_set=qs_kind_set, &
221 dft_control=dft_control, &
222 particle_set=particle_set, &
223 para_env=para_env, &
224 blacs_env=blacs_env, &
225 energy=energy, &
226 mos=mos, &
227 scf_env=scf_env, &
228 scf_control=scf_control, &
229 virial=virial, &
230 matrix_ks=matrix_ks, &
231 matrix_s=matrix_s, &
232 mp2_env=mp2_env, &
233 admm_env=admm_env)
234
235 END IF
236
237 ! IF DO_BSE In TDDFT, SAVE ks_matrix to ex_env
238 IF (dft_control%tddfpt2_control%do_bse .OR. dft_control%tddfpt2_control%do_bse_w_only .OR. &
239 dft_control%tddfpt2_control%do_bse_gw_only) THEN
240 NULLIFY (ex_env)
241 CALL get_qs_env(qs_env, exstate_env=ex_env)
242 nspins = 1 ! for now only open-shell
243 CALL dbcsr_allocate_matrix_set(ex_env%matrix_ks, nspins)
244 DO ispin = 1, nspins
245 ALLOCATE (ex_env%matrix_ks(ispin)%matrix)
246 CALL dbcsr_create(ex_env%matrix_ks(ispin)%matrix, template=matrix_s(1)%matrix)
247 CALL dbcsr_copy(ex_env%matrix_ks(ispin)%matrix, matrix_ks(ispin)%matrix)
248 END DO
249 END IF
250
251 unit_nr = cp_print_key_unit_nr(logger, input, "DFT%XC%WF_CORRELATION%PRINT", &
252 extension=".mp2Log")
253
254 IF (unit_nr > 0) THEN
255 IF (mp2_env%method /= ri_rpa_method_gpw) THEN
256 WRITE (unit_nr, *)
257 WRITE (unit_nr, *)
258 WRITE (unit_nr, '(T2,A)') 'MP2 section'
259 WRITE (unit_nr, '(T2,A)') '-----------'
260 WRITE (unit_nr, *)
261 ELSE
262 WRITE (unit_nr, *)
263 WRITE (unit_nr, *)
264 WRITE (unit_nr, '(T2,A)') 'RI-RPA section'
265 WRITE (unit_nr, '(T2,A)') '--------------'
266 WRITE (unit_nr, *)
267 END IF
268 END IF
269
270 IF (calc_forces) THEN
271 CALL cite_reference(delben2015b)
272 CALL cite_reference(rybkin2016)
273 CALL cite_reference(stein2022)
274 CALL cite_reference(bussy2023)
275 CALL cite_reference(stein2024)
276 !Gradients available for RI-MP2, and low-scaling Laplace MP2/RPA
277 IF (.NOT. (mp2_env%method == ri_mp2_method_gpw .OR. &
278 mp2_env%method == ri_rpa_method_gpw .OR. mp2_env%method == ri_mp2_laplace)) THEN
279 cpabort("No forces/gradients for the selected method.")
280 END IF
281 END IF
282
283 IF (.NOT. do_kpoints_cubic_rpa) THEN
284 IF (virial%pv_availability .AND. (.NOT. virial%pv_numer) .AND. mp2_env%eri_method == do_eri_mme) THEN
285 cpabort("analytical stress not implemented with ERI_METHOD = MME")
286 END IF
287 END IF
288
289 IF (mp2_env%do_im_time .AND. mp2_env%eri_method /= do_eri_gpw) THEN
290 mp2_env%mp2_num_proc = 1
291 END IF
292
293 IF (mp2_env%mp2_num_proc < 1 .OR. mp2_env%mp2_num_proc > para_env%num_pe) THEN
294 cpwarn("GROUP_SIZE is reset because of a too small or too large value.")
295 mp2_env%mp2_num_proc = max(min(para_env%num_pe, mp2_env%mp2_num_proc), 1)
296 END IF
297
298 IF (mod(para_env%num_pe, mp2_env%mp2_num_proc) /= 0) THEN
299 cpabort("GROUP_SIZE must be a divisor of the total number of MPI ranks!")
300 END IF
301
302 IF (.NOT. mp2_env%do_im_time) THEN
303 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T76,I5)') 'Used number of processes per group:', mp2_env%mp2_num_proc
304 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T68,F9.2,A4)') 'Maximum allowed memory usage per MPI process:', &
305 mp2_env%mp2_memory, ' MiB'
306 END IF
307
308 IF ((mp2_env%method /= mp2_method_gpw) .AND. &
309 (mp2_env%method /= ri_mp2_method_gpw) .AND. &
310 (mp2_env%method /= ri_rpa_method_gpw) .AND. &
311 (mp2_env%method /= ri_mp2_laplace)) THEN
312 CALL m_memory(mem)
313 mem_real = (mem + 1024*1024 - 1)/(1024*1024)
314 CALL para_env%max(mem_real)
315 mp2_env%mp2_memory = min(mp2_env%mp2_memory, mem_real)
316 IF (mp2_env%mp2_memory < 0.0_dp) mp2_env%mp2_memory = 1.0_dp
317
318 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T68,F9.2,A4)') 'Available memory per MPI process for MP2:', &
319 mp2_env%mp2_memory, ' MiB'
320 END IF
321
322 IF (unit_nr > 0) CALL m_flush(unit_nr)
323
324 nspins = dft_control%nspins
325 natom = SIZE(particle_set, 1)
326
327 CALL get_atomic_kind_set(atomic_kind_set, kind_of=kind_of)
328 nkind = SIZE(atomic_kind_set, 1)
329
330 do_admm_rpa_exx = mp2_env%ri_rpa%do_admm
331 IF (do_admm_rpa_exx .AND. .NOT. dft_control%do_admm) THEN
332 cpabort("Need an ADMM input section for ADMM RI_RPA EXX to work")
333 END IF
334 IF (do_admm_rpa_exx) THEN
335 CALL hfx_create_basis_types(basis_parameter, basis_info, qs_kind_set, "AUX_FIT")
336 ELSE
337 CALL hfx_create_basis_types(basis_parameter, basis_info, qs_kind_set, "ORB")
338 END IF
339
340 dimen = 0
341 max_nset = 0
342 DO iatom = 1, natom
343 ikind = kind_of(iatom)
344 dimen = dimen + sum(basis_parameter(ikind)%nsgf)
345 max_nset = max(max_nset, basis_parameter(ikind)%nset)
346 END DO
347
348 CALL get_mo_set(mo_set=mos(1), nao=nao)
349
350 ! diagonalize the KS matrix in order to have the full set of MO's
351 ! get S and KS matrices in fm_type (create also a working array)
352 NULLIFY (fm_struct)
353 CALL dbcsr_get_info(matrix_s(1)%matrix, nfullrows_total=nfullrows_total, nfullcols_total=nfullcols_total)
354 CALL cp_fm_struct_create(fm_struct, context=blacs_env, nrow_global=nfullrows_total, &
355 ncol_global=nfullcols_total, para_env=para_env)
356 CALL cp_fm_create(fm_matrix_s, fm_struct, name="fm_matrix_s")
357 CALL copy_dbcsr_to_fm(matrix_s(1)%matrix, fm_matrix_s)
358
359 CALL cp_fm_create(fm_matrix_ks, fm_struct, name="fm_matrix_ks")
360
361 CALL cp_fm_create(fm_matrix_work, fm_struct, name="fm_matrix_work")
362 CALL cp_fm_set_all(matrix=fm_matrix_work, alpha=0.0_dp)
363
364 CALL cp_fm_struct_release(fm_struct)
365
366 nmo = nao
367 ALLOCATE (nelec(nspins))
368 IF (scf_env%cholesky_method == cholesky_off) THEN
369 ALLOCATE (evals(nao))
370 evals = 0
371
372 CALL cp_fm_create(evecs, fm_matrix_s%matrix_struct)
373
374 ! Perform an EVD
375 CALL choose_eigv_solver(fm_matrix_s, evecs, evals)
376
377 ! Determine the number of neglectable eigenvalues assuming that the eigenvalues are in ascending order
378 ! (Required by Lapack)
379 ndep = 0
380 DO ii = 1, nao
381 IF (evals(ii) > scf_control%eps_eigval) THEN
382 ndep = ii - 1
383 EXIT
384 END IF
385 END DO
386 nmo = nao - ndep
387
388 DO ispin = 1, nspins
389 CALL get_mo_set(mo_set=mos(ispin), nelectron=nelec(ispin))
390 END DO
391 IF (maxval(nelec)/(3 - nspins) > nmo) THEN
392 ! Should not happen as the following MO calculation is the same as during the SCF steps
393 cpabort("Not enough MOs found!")
394 END IF
395
396 ! Set the eigenvalue of the eigenvectors belonging to the linear subspace to zero
397 evals(1:ndep) = 0.0_dp
398 ! Determine the eigenvalues of the inverse square root
399 evals(ndep + 1:nao) = 1.0_dp/sqrt(evals(ndep + 1:nao))
400
401 IF (ndep > 0) THEN
402 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T76,I5)') 'Number of removed MOs:', ndep
403 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T76,I5)') 'Number of available MOs:', nmo
404
405 ! Create reduced matrices
406 NULLIFY (fm_struct)
407 CALL cp_fm_struct_create(fm_struct, template_fmstruct=fm_matrix_s%matrix_struct, ncol_global=nmo)
408
409 ALLOCATE (fm_matrix_s_red, fm_work_red)
410 CALL cp_fm_create(fm_matrix_s_red, fm_struct)
411 CALL cp_fm_create(fm_work_red, fm_struct)
412 CALL cp_fm_struct_release(fm_struct)
413
414 ALLOCATE (fm_matrix_ks_red)
415 CALL cp_fm_struct_create(fm_struct, template_fmstruct=fm_matrix_s%matrix_struct, &
416 nrow_global=nmo, ncol_global=nmo)
417 CALL cp_fm_create(fm_matrix_ks_red, fm_struct)
418 CALL cp_fm_struct_release(fm_struct)
419
420 ! Scale the eigenvalues and copy them to
421 CALL cp_fm_to_fm(evecs, fm_matrix_s_red, nmo, ndep + 1)
422 CALL cp_fm_column_scale(fm_matrix_s_red, evals(ndep + 1:))
423
424 ! Obtain ortho from (P)DGEMM, skip the linear dependent columns
425 CALL parallel_gemm("N", "T", nao, nao, nmo, 1.0_dp, fm_matrix_s_red, evecs, &
426 0.0_dp, fm_matrix_s, b_first_col=ndep + 1)
427 ELSE
428 ! Take the square roots of the target values to allow application of SYRK
429 evals = sqrt(evals)
430 CALL cp_fm_column_scale(evecs, evals)
431 CALL cp_fm_syrk("U", "N", nao, 1.0_dp, evecs, 1, 1, 0.0_dp, fm_matrix_s)
432 CALL cp_fm_uplo_to_full(fm_matrix_s, fm_matrix_work)
433 END IF
434
435 CALL cp_fm_release(evecs)
436 cholesky_method = cholesky_off
437 ELSE
438 ! calculate S^(-1/2) (cholesky decomposition)
439 CALL cp_fm_cholesky_decompose(fm_matrix_s)
440 CALL cp_fm_triangular_invert(fm_matrix_s)
441 cholesky_method = cholesky_inverse
442 END IF
443
444 ALLOCATE (mos_mp2(nspins))
445 DO ispin = 1, nspins
446
447 CALL get_mo_set(mo_set=mos(ispin), maxocc=maxocc, nelectron=nelec(ispin))
448
449 CALL allocate_mo_set(mo_set=mos_mp2(ispin), &
450 nao=nao, &
451 nmo=nmo, &
452 nelectron=nelec(ispin), &
453 n_el_f=real(nelec(ispin), dp), &
454 maxocc=maxocc, &
455 flexible_electron_count=dft_control%relax_multiplicity)
456
457 CALL get_mo_set(mos_mp2(ispin), nao=nao)
458 CALL cp_fm_struct_create(fm_struct, nrow_global=nao, &
459 ncol_global=nmo, para_env=para_env, &
460 context=blacs_env)
461
462 CALL init_mo_set(mos_mp2(ispin), &
463 fm_struct=fm_struct, &
464 name="mp2_mos")
465 CALL cp_fm_struct_release(fm_struct)
466 END DO
467
468 DO ispin = 1, nspins
469
470 ! If ADMM we should make the ks matrix up-to-date
471 IF (dft_control%do_admm) THEN
472 CALL admm_correct_for_eigenvalues(ispin, admm_env, matrix_ks(ispin)%matrix)
473 END IF
474
475 CALL copy_dbcsr_to_fm(matrix_ks(ispin)%matrix, fm_matrix_ks)
476
477 IF (dft_control%do_admm) THEN
478 CALL admm_uncorrect_for_eigenvalues(ispin, admm_env, matrix_ks(ispin)%matrix)
479 END IF
480
481 IF (cholesky_method == cholesky_inverse) THEN
482
483 ! diagonalize KS matrix
484 CALL eigensolver(matrix_ks_fm=fm_matrix_ks, &
485 mo_set=mos_mp2(ispin), &
486 ortho=fm_matrix_s, &
487 work=fm_matrix_work, &
488 cholesky_method=cholesky_method, &
489 do_level_shift=.false., &
490 level_shift=0.0_dp, &
491 use_jacobi=.false.)
492
493 ELSE IF (cholesky_method == cholesky_off) THEN
494
495 IF (ASSOCIATED(fm_matrix_s_red)) THEN
496 CALL eigensolver_symm(matrix_ks_fm=fm_matrix_ks, &
497 mo_set=mos_mp2(ispin), &
498 ortho=fm_matrix_s, &
499 work=fm_matrix_work, &
500 do_level_shift=.false., &
501 level_shift=0.0_dp, &
502 use_jacobi=.false., &
503 jacobi_threshold=0.0_dp, &
504 ortho_red=fm_matrix_s_red, &
505 matrix_ks_fm_red=fm_matrix_ks_red, &
506 work_red=fm_work_red)
507 ELSE
508 CALL eigensolver_symm(matrix_ks_fm=fm_matrix_ks, &
509 mo_set=mos_mp2(ispin), &
510 ortho=fm_matrix_s, &
511 work=fm_matrix_work, &
512 do_level_shift=.false., &
513 level_shift=0.0_dp, &
514 use_jacobi=.false., &
515 jacobi_threshold=0.0_dp)
516 END IF
517 END IF
518
519 CALL get_mo_set(mos_mp2(ispin), mo_coeff=mo_coeff)
520 END DO
521
522 CALL cp_fm_release(fm_matrix_s)
523 CALL cp_fm_release(fm_matrix_ks)
524 CALL cp_fm_release(fm_matrix_work)
525 IF (ASSOCIATED(fm_matrix_s_red)) THEN
526 CALL cp_fm_release(fm_matrix_s_red)
527 DEALLOCATE (fm_matrix_s_red)
528 END IF
529 IF (ASSOCIATED(fm_matrix_ks_red)) THEN
530 CALL cp_fm_release(fm_matrix_ks_red)
531 DEALLOCATE (fm_matrix_ks_red)
532 END IF
533 IF (ASSOCIATED(fm_work_red)) THEN
534 CALL cp_fm_release(fm_work_red)
535 DEALLOCATE (fm_work_red)
536 END IF
537
538 hfx_sections => section_vals_get_subs_vals(input, "DFT%XC%HF")
539
540 ! build the table of index
541 t1 = m_walltime()
542 ALLOCATE (mp2_biel%index_table(natom, max_nset))
543
544 CALL build_index_table(natom, max_nset, mp2_biel%index_table, basis_parameter, kind_of)
545
546 ! free the hfx_container (for now if forces are required we don't release the HFX stuff)
547 free_hfx_buffer = .false.
548 IF (ASSOCIATED(qs_env%x_data)) THEN
549 free_hfx_buffer = .true.
550 IF (calc_forces .AND. (.NOT. mp2_env%ri_grad%free_hfx_buffer)) free_hfx_buffer = .false.
551 IF (qs_env%x_data(1, 1)%do_hfx_ri) free_hfx_buffer = .false.
552 IF (calc_forces .AND. mp2_env%do_im_time) free_hfx_buffer = .false.
553 IF (mp2_env%ri_rpa%reuse_hfx) free_hfx_buffer = .false.
554 END IF
555
556 IF (.NOT. do_kpoints_cubic_rpa) THEN
557 IF (virial%pv_numer) THEN
558 ! in the case of numerical stress we don't have to free the HFX integrals
559 free_hfx_buffer = .false.
560 mp2_env%ri_grad%free_hfx_buffer = free_hfx_buffer
561 END IF
562 END IF
563
564 ! calculate the matrix sigma_x - vxc for G0W0
565 t3 = 0
566 IF (mp2_env%ri_rpa%do_ri_g0w0) THEN
567 CALL compute_vec_sigma_x_minus_vxc_gw(qs_env, mp2_env, mos_mp2, e_ex_from_gw, e_admm_from_gw, t3, unit_nr)
568 END IF
569
570 IF (free_hfx_buffer) THEN
571 CALL timeset(routinen//"_free_hfx", handle2)
572 CALL section_vals_get(hfx_sections, n_repetition=n_rep_hf)
573 n_threads = 1
574!$ n_threads = omp_get_max_threads()
575
576 DO irep = 1, n_rep_hf
577 DO i_thread = 0, n_threads - 1
578 actual_x_data => qs_env%x_data(irep, i_thread + 1)
579
580 do_dynamic_load_balancing = .true.
581 IF (n_threads == 1 .OR. actual_x_data%memory_parameter%do_disk_storage) do_dynamic_load_balancing = .false.
582
583 IF (do_dynamic_load_balancing) THEN
584 my_bin_size = SIZE(actual_x_data%distribution_energy)
585 ELSE
586 my_bin_size = 1
587 END IF
588
589 IF (.NOT. actual_x_data%memory_parameter%do_all_on_the_fly) THEN
590 CALL dealloc_containers(actual_x_data%store_ints, actual_x_data%memory_parameter%actual_memory_usage)
591 END IF
592 END DO
593 END DO
594 CALL timestop(handle2)
595 END IF
596
597 emp2 = 0.d+00
598 emp2_cou = 0.d+00
599 emp2_ex = 0.d+00
600 calc_ex = .true.
601
602 t1 = m_walltime()
603 SELECT CASE (mp2_env%method)
604 CASE (mp2_method_direct)
605 IF (unit_nr > 0) WRITE (unit_nr, *)
606
607 ALLOCATE (auto(dimen, nspins))
608 ALLOCATE (c(dimen, dimen, nspins))
609
610 DO ispin = 1, nspins
611 ! get the alpha coeff and eigenvalues
612 CALL get_mo_set(mo_set=mos_mp2(ispin), &
613 eigenvalues=mo_eigenvalues, &
614 mo_coeff=mo_coeff)
615
616 CALL cp_fm_get_submatrix(mo_coeff, c(:, :, ispin), 1, 1, dimen, dimen, .false.)
617 auto(:, ispin) = mo_eigenvalues(:)
618 END DO
619
620 IF (nspins == 2) THEN
621 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Unrestricted Canonical Direct Methods:'
622 ! for now, require the mos to be always present
623
624 ! calculate the alpha-alpha MP2
625 emp2_aa = 0.0_dp
626 emp2_aa_cou = 0.0_dp
627 emp2_aa_ex = 0.0_dp
628 CALL mp2_direct_energy(dimen, nelec(1), nelec(1), mp2_biel, &
629 mp2_env, c(:, :, 1), auto(:, 1), emp2_aa, emp2_aa_cou, emp2_aa_ex, &
630 qs_env, para_env, unit_nr)
631 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy Alpha-Alpha = ', emp2_aa
632 IF (unit_nr > 0) WRITE (unit_nr, *)
633
634 emp2_bb = 0.0_dp
635 emp2_bb_cou = 0.0_dp
636 emp2_bb_ex = 0.0_dp
637 CALL mp2_direct_energy(dimen, nelec(2), nelec(2), mp2_biel, mp2_env, &
638 c(:, :, 2), auto(:, 2), emp2_bb, emp2_bb_cou, emp2_bb_ex, &
639 qs_env, para_env, unit_nr)
640 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy Beta-Beta= ', emp2_bb
641 IF (unit_nr > 0) WRITE (unit_nr, *)
642
643 emp2_ab = 0.0_dp
644 emp2_ab_cou = 0.0_dp
645 emp2_ab_ex = 0.0_dp
646 CALL mp2_direct_energy(dimen, nelec(1), nelec(2), mp2_biel, mp2_env, c(:, :, 1), &
647 auto(:, 1), emp2_ab, emp2_ab_cou, emp2_ab_ex, &
648 qs_env, para_env, unit_nr, c(:, :, 2), auto(:, 2))
649 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy Alpha-Beta= ', emp2_ab
650 IF (unit_nr > 0) WRITE (unit_nr, *)
651
652 emp2 = emp2_aa + emp2_bb + emp2_ab*2.0_dp !+Emp2_BA
653 emp2_cou = emp2_aa_cou + emp2_bb_cou + emp2_ab_cou*2.0_dp !+Emp2_BA
654 emp2_ex = emp2_aa_ex + emp2_bb_ex + emp2_ab_ex*2.0_dp !+Emp2_BA
655
656 emp2_s = emp2_ab*2.0_dp
657 emp2_t = emp2_aa + emp2_bb
658
659 ELSE
660
661 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A)') 'Canonical Direct Methods:'
662
663 CALL mp2_direct_energy(dimen, nelec(1)/2, nelec(1)/2, mp2_biel, mp2_env, &
664 c(:, :, 1), auto(:, 1), emp2, emp2_cou, emp2_ex, &
665 qs_env, para_env, unit_nr)
666
667 END IF
668
669 DEALLOCATE (c, auto)
670
672 ! optimize ri basis set or tests for RI-MP2 gradients
673 IF (unit_nr > 0) THEN
674 WRITE (unit_nr, *)
675 WRITE (unit_nr, '(T3,A)') 'Optimization of the auxiliary RI-MP2 basis'
676 WRITE (unit_nr, *)
677 END IF
678
679 ALLOCATE (auto(dimen, nspins))
680 ALLOCATE (c(dimen, dimen, nspins))
681
682 DO ispin = 1, nspins
683 ! get the alpha coeff and eigenvalues
684 CALL get_mo_set(mo_set=mos_mp2(ispin), &
685 eigenvalues=mo_eigenvalues, &
686 mo_coeff=mo_coeff)
687
688 CALL cp_fm_get_submatrix(mo_coeff, c(:, :, ispin), 1, 1, dimen, dimen, .false.)
689 auto(:, ispin) = mo_eigenvalues(:)
690 END DO
691
692 ! optimize basis
693 IF (nspins == 2) THEN
694 CALL optimize_ri_basis_main(emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, dimen, natom, nelec(1), &
695 mp2_biel, mp2_env, c(:, :, 1), auto(:, 1), &
696 kind_of, qs_env, para_env, unit_nr, &
697 nelec(2), c(:, :, 2), auto(:, 2))
698
699 ELSE
700 CALL optimize_ri_basis_main(emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, dimen, natom, nelec(1)/2, &
701 mp2_biel, mp2_env, c(:, :, 1), auto(:, 1), &
702 kind_of, qs_env, para_env, unit_nr)
703 END IF
704
705 DEALLOCATE (auto, c)
706
707 CASE (mp2_method_gpw)
708 ! check if calculate the exchange contribution
709 IF (mp2_env%scale_T == 0.0_dp .AND. (nspins == 2)) calc_ex = .false.
710
711 ! go with mp2_gpw
712 CALL mp2_gpw_main(qs_env, mp2_env, emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, &
713 mos_mp2, para_env, unit_nr, calc_forces, calc_ex)
714
715 CASE (ri_mp2_method_gpw)
716 ! check if calculate the exchange contribution
717 IF (mp2_env%scale_T == 0.0_dp .AND. (nspins == 2)) calc_ex = .false.
718
719 ! go with mp2_gpw
720 CALL mp2_gpw_main(qs_env, mp2_env, emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, &
721 mos_mp2, para_env, unit_nr, calc_forces, calc_ex, do_ri_mp2=.true.)
722
723 CASE (ri_rpa_method_gpw)
724 ! perform RI-RPA energy calculation (since most part of the calculation
725 ! is actually equal to the RI-MP2-GPW we decided to put RPA in the MP2
726 ! section to avoid code replication)
727
728 calc_ex = .false.
729
730 ! go with ri_rpa_gpw
731 CALL mp2_gpw_main(qs_env, mp2_env, emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, &
732 mos_mp2, para_env, unit_nr, calc_forces, calc_ex, do_ri_rpa=.true.)
733 ! Scale energy contributions
734 emp2 = emp2*mp2_env%ri_rpa%scale_rpa
735 mp2_env%ri_rpa%ener_exchange = mp2_env%ri_rpa%ener_exchange*mp2_env%ri_rpa%scale_rpa
736
737 CASE (ri_mp2_laplace)
738 ! perform RI-SOS-Laplace-MP2 energy calculation, most part of the code in common
739 ! with the RI-RPA part
740
741 ! In SOS-MP2 only the coulomb-like contribution of the MP2 energy is computed
742 calc_ex = .false.
743
744 ! go with sos_laplace_mp2_gpw
745 CALL mp2_gpw_main(qs_env, mp2_env, emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, &
746 mos_mp2, para_env, unit_nr, calc_forces, calc_ex, do_ri_sos_laplace_mp2=.true.)
747
748 CASE DEFAULT
749 cpabort("")
750 END SELECT
751
752 t2 = m_walltime()
753 IF (unit_nr > 0) WRITE (unit_nr, *)
754 IF (mp2_env%method /= ri_rpa_method_gpw) THEN
755 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.6)') 'Total MP2 Time=', t2 - t1
756 IF (mp2_env%method == ri_mp2_laplace) THEN
757 emp2_s = emp2
758 emp2_t = 0.0_dp
759 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy SO component (singlet) = ', emp2_s
760 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Scaling factor SO = ', mp2_env%scale_S
761 ELSE
762 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Coulomb Energy = ', emp2_cou/2.0_dp
763 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Exchange Energy = ', emp2_ex
764 IF (nspins == 1) THEN
765 ! valid only in the closed shell case
766 emp2_s = emp2_cou/2.0_dp
767 IF (calc_ex) THEN
768 emp2_t = emp2_ex + emp2_cou/2.0_dp
769 ELSE
770 ! unknown if Emp2_ex is not computed
771 emp2_t = 0.0_dp
772 END IF
773 END IF
774 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy SO component (singlet) = ', emp2_s
775 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'MP2 Energy SS component (triplet) = ', emp2_t
776 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Scaling factor SO = ', mp2_env%scale_S
777 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Scaling factor SS = ', mp2_env%scale_T
778 END IF
779 emp2_s = emp2_s*mp2_env%scale_S
780 emp2_t = emp2_t*mp2_env%scale_T
781 emp2 = emp2_s + emp2_t
782 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Second order perturbation energy = ', emp2
783 ELSE
784 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.6)') 'Total RI-RPA Time=', t2 - t1
785
786 update_xc_energy = .true.
787 IF (mp2_env%ri_rpa%do_ri_g0w0 .AND. .NOT. mp2_env%ri_g0w0%update_xc_energy) update_xc_energy = .false.
788 IF (.NOT. update_xc_energy) emp2 = 0.0_dp
789
790 IF (unit_nr > 0 .AND. update_xc_energy) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'RI-RPA energy = ', emp2
791 IF (unit_nr > 0 .AND. mp2_env%ri_rpa%sigma_param /= sigma_none) THEN
792 WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Sigma corr. to RI-RPA energy = ', &
793 mp2_env%ri_rpa%e_sigma_corr
794 END IF
795 IF (mp2_env%ri_rpa%exchange_correction == rpa_exchange_axk) THEN
796 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'RI-RPA-AXK energy=', mp2_env%ri_rpa%ener_exchange
797 ELSE IF (mp2_env%ri_rpa%exchange_correction == rpa_exchange_sosex) THEN
798 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'RI-RPA-SOSEX energy=', mp2_env%ri_rpa%ener_exchange
799 END IF
800 IF (mp2_env%ri_rpa%do_rse) THEN
801 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Diagonal singles correction (dRSE) = ', &
802 mp2_env%ri_rpa%rse_corr_diag
803 IF (unit_nr > 0) WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Full singles correction (RSE) =', &
804 mp2_env%ri_rpa%rse_corr
805 IF (dft_control%do_admm) cpabort("RPA RSE not implemented with RI_RPA%ADMM on")
806 END IF
807 END IF
808 IF (unit_nr > 0) WRITE (unit_nr, *)
809
810 ! we have it !!!!
811 IF (mp2_env%ri_rpa%exchange_correction /= rpa_exchange_none) THEN
812 emp2 = emp2 + mp2_env%ri_rpa%ener_exchange
813 END IF
814 IF (mp2_env%ri_rpa%do_rse) THEN
815 emp2 = emp2 + mp2_env%ri_rpa%rse_corr
816 END IF
817 IF (mp2_env%ri_rpa%sigma_param /= sigma_none) THEN
818 !WRITE (unit_nr, '(T3,A,T56,F25.14)') 'Sigma corr. to RI-RPA energy = ',&
819 emp2 = emp2 + mp2_env%ri_rpa%e_sigma_corr
820 END IF
821 energy%mp2 = emp2
822 energy%total = energy%total + emp2
823
824 DO ispin = 1, nspins
825 CALL deallocate_mo_set(mo_set=mos_mp2(ispin))
826 END DO
827 DEALLOCATE (mos_mp2)
828
829 ! if necessary reallocate hfx buffer
830 IF (free_hfx_buffer .AND. (.NOT. calc_forces) .AND. &
831 (mp2_env%ri_g0w0%do_ri_Sigma_x .OR. .NOT. mp2_env%ri_rpa_im_time%do_kpoints_from_Gamma)) THEN
832 CALL timeset(routinen//"_alloc_hfx", handle2)
833 DO irep = 1, n_rep_hf
834 DO i_thread = 0, n_threads - 1
835 actual_x_data => qs_env%x_data(irep, i_thread + 1)
836
837 do_dynamic_load_balancing = .true.
838 IF (n_threads == 1 .OR. actual_x_data%memory_parameter%do_disk_storage) do_dynamic_load_balancing = .false.
839
840 IF (do_dynamic_load_balancing) THEN
841 my_bin_size = SIZE(actual_x_data%distribution_energy)
842 ELSE
843 my_bin_size = 1
844 END IF
845
846 IF (.NOT. actual_x_data%memory_parameter%do_all_on_the_fly) THEN
847 CALL alloc_containers(actual_x_data%store_ints, my_bin_size)
848
849 DO bin = 1, my_bin_size
850 maxval_container => actual_x_data%store_ints%maxval_container(bin)
851 integral_containers => actual_x_data%store_ints%integral_containers(:, bin)
852 CALL hfx_init_container(maxval_container, actual_x_data%memory_parameter%actual_memory_usage, .false.)
853 DO i = 1, 64
854 CALL hfx_init_container(integral_containers(i), actual_x_data%memory_parameter%actual_memory_usage, .false.)
855 END DO
856 END DO
857 END IF
858 END DO
859 END DO
860 CALL timestop(handle2)
861 END IF
862
863 CALL hfx_release_basis_types(basis_parameter)
864
865 ! if required calculate the EXX contribution from the DFT density
866 IF (mp2_env%method == ri_rpa_method_gpw .AND. .NOT. calc_forces) THEN
867 do_exx = .false.
868 hfx_sections => section_vals_get_subs_vals(input, "DFT%XC%WF_CORRELATION%RI_RPA%HF")
869 CALL section_vals_get(hfx_sections, explicit=do_exx)
870 IF (do_exx) THEN
871 do_gw = mp2_env%ri_rpa%do_ri_g0w0
872 do_admm = mp2_env%ri_rpa%do_admm
873 reuse_hfx = qs_env%mp2_env%ri_rpa%reuse_hfx
874 do_im_time = qs_env%mp2_env%do_im_time
875
876 CALL calculate_exx(qs_env=qs_env, &
877 unit_nr=unit_nr, &
878 hfx_sections=hfx_sections, &
879 x_data=qs_env%mp2_env%ri_rpa%x_data, &
880 do_gw=do_gw, &
881 do_admm=do_admm, &
882 calc_forces=.false., &
883 reuse_hfx=reuse_hfx, &
884 do_im_time=do_im_time, &
885 e_ex_from_gw=e_ex_from_gw, &
886 e_admm_from_gw=e_admm_from_gw, &
887 t3=t3)
888
889 END IF
890 END IF
891
892 CALL cp_print_key_finished_output(unit_nr, logger, input, &
893 "DFT%XC%WF_CORRELATION%PRINT")
894
895 CALL timestop(handle)
896
897 END SUBROUTINE mp2_main
898
899! **************************************************************************************************
900!> \brief ...
901!> \param natom ...
902!> \param max_nset ...
903!> \param index_table ...
904!> \param basis_parameter ...
905!> \param kind_of ...
906! **************************************************************************************************
907 PURE SUBROUTINE build_index_table(natom, max_nset, index_table, basis_parameter, kind_of)
908 INTEGER, INTENT(IN) :: natom, max_nset
909 INTEGER, DIMENSION(natom, max_nset), INTENT(OUT) :: index_table
910 TYPE(hfx_basis_type), DIMENSION(:), POINTER :: basis_parameter
911 INTEGER, DIMENSION(natom), INTENT(IN) :: kind_of
912
913 INTEGER :: counter, iatom, ikind, iset, nset
914
915 index_table = -huge(0)
916 counter = 0
917 DO iatom = 1, natom
918 ikind = kind_of(iatom)
919 nset = basis_parameter(ikind)%nset
920 DO iset = 1, nset
921 index_table(iatom, iset) = counter + 1
922 counter = counter + basis_parameter(ikind)%nsgf(iset)
923 END DO
924 END DO
925
926 END SUBROUTINE build_index_table
927
928! **************************************************************************************************
929!> \brief ...
930!> \param matrix_s ...
931!> \param matrix_ks ...
932!> \param mos ...
933!> \param matrix_s_kp ...
934!> \param matrix_ks_transl ...
935!> \param kpoints ...
936! **************************************************************************************************
937 PURE SUBROUTINE get_gamma(matrix_s, matrix_ks, mos, matrix_s_kp, matrix_ks_transl, kpoints)
938
939 TYPE(dbcsr_p_type), DIMENSION(:), POINTER :: matrix_s, matrix_ks
940 TYPE(mo_set_type), DIMENSION(:), POINTER :: mos
941 TYPE(dbcsr_p_type), DIMENSION(:, :), POINTER :: matrix_s_kp, matrix_ks_transl
942 TYPE(kpoint_type), POINTER :: kpoints
943
944 INTEGER :: nspins
945
946 nspins = SIZE(matrix_ks_transl, 1)
947
948 matrix_ks(1:nspins) => matrix_ks_transl(1:nspins, 1)
949 matrix_s(1:1) => matrix_s_kp(1:1, 1)
950 mos(1:nspins) => kpoints%kp_env(1)%kpoint_env%mos(1:nspins, 1)
951
952 END SUBROUTINE get_gamma
953
954END MODULE mp2
955
Types and set/get functions for auxiliary density matrix methods.
Definition admm_types.F:15
Contains methods used in the context of density fitting.
Definition admm_utils.F:15
subroutine, public admm_uncorrect_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:127
subroutine, public admm_correct_for_eigenvalues(ispin, admm_env, ks_matrix)
...
Definition admm_utils.F:53
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind_set(atomic_kind_set, atom_of_kind, kind_of, natom_of_kind, maxatom, natom, nshell, fist_potential_present, shell_present, shell_adiabatic, shell_check_distance, damping_present)
Get attributes of an atomic kind set.
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public stein2024
integer, save, public stein2022
integer, save, public rybkin2016
integer, save, public delben2012
integer, save, public delben2015b
integer, save, public bussy2023
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
subroutine, public dbcsr_copy(matrix_b, matrix_a, name, keep_sparsity, keep_imaginary)
...
subroutine, public dbcsr_get_info(matrix, nblkrows_total, nblkcols_total, nfullrows_total, nfullcols_total, nblkrows_local, nblkcols_local, nfullrows_local, nfullcols_local, my_prow, my_pcol, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, row_blk_offset, col_blk_offset, distribution, name, matrix_type, group)
...
DBCSR operations in CP2K.
subroutine, public copy_dbcsr_to_fm(matrix, fm)
Copy a DBCSR matrix to a BLACS matrix.
Basic linear algebra operations for full matrices.
subroutine, public cp_fm_column_scale(matrixa, scaling)
scales column i of matrix a with scaling(i)
subroutine, public cp_fm_syrk(uplo, trans, k, alpha, matrix_a, ia, ja, beta, matrix_c)
performs a rank-k update of a symmetric matrix_c matrix_c = beta * matrix_c + alpha * matrix_a * tran...
subroutine, public cp_fm_uplo_to_full(matrix, work, uplo)
given a triangular matrix according to uplo, computes the corresponding full matrix
subroutine, public cp_fm_triangular_invert(matrix_a, uplo_tr)
inverts a triangular matrix
various cholesky decomposition related routines
subroutine, public cp_fm_cholesky_decompose(matrix, n, info_out)
used to replace a symmetric positive def. matrix M with its cholesky decomposition U: M = U^T * U,...
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
subroutine, public choose_eigv_solver(matrix, eigenvectors, eigenvalues, info)
Choose the Eigensolver depending on which library is available ELPA seems to be unstable for small sy...
Definition cp_fm_diag.F:239
represent the structure of a full matrix
subroutine, public cp_fm_struct_create(fmstruct, para_env, context, nrow_global, ncol_global, nrow_block, ncol_block, descriptor, first_p_pos, local_leading_dimension, template_fmstruct, square_blocks, force_block)
allocates and initializes a full matrix structure
subroutine, public cp_fm_struct_release(fmstruct)
releases a full matrix structure
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
subroutine, public cp_fm_create(matrix, matrix_struct, name, use_sp, nrow, ncol, set_zero)
creates a new full matrix with the given structure
subroutine, public cp_fm_set_all(matrix, alpha, beta)
set all elements of a matrix to the same value, and optionally the diagonal to a different one
subroutine, public cp_fm_get_submatrix(fm, target_m, start_row, start_col, n_rows, n_cols, transpose)
gets a submatrix of a full matrix op(target_m)(1:n_rows,1:n_cols) =fm(start_row:start_row+n_rows,...
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
Types for excited states potential energies.
Routines to calculate EXX in RPA and energy correction methods.
Definition hfx_exx.F:16
subroutine, public calculate_exx(qs_env, unit_nr, hfx_sections, x_data, do_gw, do_admm, calc_forces, reuse_hfx, do_im_time, e_ex_from_gw, e_admm_from_gw, t3)
...
Definition hfx_exx.F:106
Types and set/get functions for HFX.
Definition hfx_types.F:15
subroutine, public hfx_init_container(container, memory_usage, do_disk_storage)
This routine deletes all list entries in a container in order to deallocate the memory.
Definition hfx_types.F:2552
subroutine, public hfx_release_basis_types(basis_parameter)
...
Definition hfx_types.F:1806
subroutine, public alloc_containers(data, bin_size)
...
Definition hfx_types.F:2935
subroutine, public dealloc_containers(data, memory_usage)
...
Definition hfx_types.F:2903
subroutine, public hfx_create_basis_types(basis_parameter, basis_info, qs_kind_set, basis_type)
This routine allocates and initializes the basis_info and basis_parameter types
Definition hfx_types.F:1680
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public rpa_exchange_none
integer, parameter, public sigma_none
integer, parameter, public mp2_method_direct
integer, parameter, public rpa_exchange_sosex
integer, parameter, public mp2_ri_optimize_basis
integer, parameter, public do_eri_mme
integer, parameter, public cholesky_off
integer, parameter, public cholesky_inverse
integer, parameter, public ri_rpa_method_gpw
integer, parameter, public ri_mp2_method_gpw
integer, parameter, public mp2_method_gpw
integer, parameter, public ri_mp2_laplace
integer, parameter, public rpa_exchange_axk
integer, parameter, public do_eri_gpw
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public int_8
Definition kinds.F:54
integer, parameter, public dp
Definition kinds.F:34
Types and basic routines needed for a kpoint calculation.
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_memory(mem)
Returns the total amount of memory [bytes] in use, if known, zero otherwise.
Definition machine.F:452
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:136
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:153
Interface to the message passing library MPI.
Routines to calculate MP2 energy.
subroutine, public mp2_direct_energy(dimen, occ_i, occ_j, mp2_biel, mp2_env, c_i, auto_i, emp2, emp2_cou, emp2_ex, qs_env, para_env, unit_nr, c_j, auto_j)
...
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
Routines to optimize the RI-MP2 basis. Only exponents of non-contracted auxiliary basis basis are opt...
subroutine, public optimize_ri_basis_main(emp2, emp2_cou, emp2_ex, emp2_s, emp2_t, dimen, natom, homo, mp2_biel, mp2_env, c, auto, kind_of, qs_env, para_env, unit_nr, homo_beta, c_beta, auto_beta)
optimize RI-MP2 basis set
Types needed for MP2 calculations.
Definition mp2_types.F:14
Routines to calculate MP2 energy.
Definition mp2.F:14
subroutine, public mp2_main(qs_env, calc_forces)
the main entry point for MP2 calculations
Definition mp2.F:128
basic linear algebra operations for full matrixes
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, mimic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, sab_cneo, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, xcint_weights, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, rhoz_cneo_set, ecoul_1c, rho0_s_rs, rho0_s_gs, rhoz_cneo_s_rs, rhoz_cneo_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs, do_rixs, tb_tblite)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
Definition and initialisation of the mo data type.
Definition qs_mo_types.F:22
subroutine, public allocate_mo_set(mo_set, nao, nmo, nelectron, n_el_f, maxocc, flexible_electron_count)
Allocates a mo set and partially initializes it (nao,nmo,nelectron, and flexible_electron_count are v...
subroutine, public deallocate_mo_set(mo_set)
Deallocate a wavefunction data structure.
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.
subroutine, public init_mo_set(mo_set, fm_pool, fm_ref, fm_struct, name)
initializes an allocated mo_set. eigenvalues, mo_coeff, occupation_numbers are valid only after this ...
groups fairly general SCF methods, so that modules other than qs_scf can use them too split off from ...
subroutine, public eigensolver(matrix_ks_fm, mo_set, ortho, work, cholesky_method, do_level_shift, level_shift, matrix_u_fm, use_jacobi)
Diagonalise the Kohn-Sham matrix to get a new set of MO eigen- vectors and MO eigenvalues....
subroutine, public eigensolver_symm(matrix_ks_fm, mo_set, ortho, work, do_level_shift, level_shift, matrix_u_fm, use_jacobi, jacobi_threshold, ortho_red, work_red, matrix_ks_fm_red, matrix_u_fm_red)
...
module that contains the definitions of the scf types
Routines to calculate EXX within GW.
subroutine, public compute_vec_sigma_x_minus_vxc_gw(qs_env, mp2_env, mos_mp2, energy_ex, energy_xc_admm, t3, unit_nr)
...
parameters that control an scf iteration
stores some data used in wavefunction fitting
Definition admm_types.F:120
Provides all information about an atomic kind.
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
keeps the information about the structure of a full matrix
represent a full matrix
type of a logger, at the moment it contains just a print level starting at which level it should be l...
Contains information on the excited states energy.
stores some data used in construction of Kohn-Sham matrix
Definition hfx_types.F:511
Contains information about kpoints.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.