(git:1155b05)
Loading...
Searching...
No Matches
ec_environment.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 Energy correction environment setup and handling
10!> \par History
11!> 2019.09 created
12!> \author JGH
13! **************************************************************************************************
22 USE bibliography, ONLY: niklasson2003,&
24 cite_reference
25 USE cell_types, ONLY: cell_type
31 USE input_constants, ONLY: &
44 USE kinds, ONLY: dp
47 USE kpoint_types, ONLY: kpoint_create,&
62 USE qs_kind_types, ONLY: get_qs_kind,&
65 USE qs_rho_types, ONLY: qs_rho_type
71#include "./base/base_uses.f90"
72
73 IMPLICIT NONE
74
75 PRIVATE
76
77 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ec_environment'
78
79 PUBLIC :: ec_env_create
80 PUBLIC :: ec_write_input
81
82CONTAINS
83
84! **************************************************************************************************
85!> \brief Allocates and intitializes ec_env
86!> \param qs_env The QS environment
87!> \param ec_env The energy correction environment (the object to create)
88!> \param dft_section The DFT section
89!> \param ec_section The energy correction input section
90!> \par History
91!> 2019.09 created
92!> \author JGH
93! **************************************************************************************************
94 SUBROUTINE ec_env_create(qs_env, ec_env, dft_section, ec_section)
95 TYPE(qs_environment_type), POINTER :: qs_env
96 TYPE(energy_correction_type), POINTER :: ec_env
97 TYPE(section_vals_type), POINTER :: dft_section
98 TYPE(section_vals_type), OPTIONAL, POINTER :: ec_section
99
100 cpassert(.NOT. ASSOCIATED(ec_env))
101 ALLOCATE (ec_env)
102 CALL init_ec_env(qs_env, ec_env, dft_section, ec_section)
103
104 END SUBROUTINE ec_env_create
105
106! **************************************************************************************************
107!> \brief Initializes energy correction environment
108!> \param qs_env The QS environment
109!> \param ec_env The energy correction environment
110!> \param dft_section The DFT section
111!> \param ec_section The energy correction input section
112!> \par History
113!> 2019.09 created
114!> \author JGH
115! **************************************************************************************************
116 SUBROUTINE init_ec_env(qs_env, ec_env, dft_section, ec_section)
117 TYPE(qs_environment_type), POINTER :: qs_env
118 TYPE(energy_correction_type), POINTER :: ec_env
119 TYPE(section_vals_type), POINTER :: dft_section
120 TYPE(section_vals_type), OPTIONAL, POINTER :: ec_section
121
122 CHARACTER(LEN=*), PARAMETER :: routinen = 'init_ec_env'
123
124 INTEGER :: handle, ikind, maxlgto, nkind, unit_nr
125 LOGICAL :: explicit, gpw, gs_kpoints, paw_atom
126 REAL(kind=dp) :: eps_pgf_orb, etemp, rc
127 TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
128 TYPE(cell_type), POINTER :: cell
129 TYPE(cp_blacs_env_type), POINTER :: blacs_env
130 TYPE(dft_control_type), POINTER :: dft_control
131 TYPE(gto_basis_set_type), POINTER :: basis_set, harris_basis, &
132 harris_soft_basis
133 TYPE(mp_para_env_type), POINTER :: para_env
134 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
135 TYPE(qs_dispersion_type), POINTER :: dispersion_env
136 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
137 TYPE(qs_kind_type), POINTER :: qs_kind
138 TYPE(qs_rho_type), POINTER :: rho
139 TYPE(section_vals_type), POINTER :: ec_hfx_section, kp_section, nl_section, &
140 pp_section, section1, section2, &
141 xc_fun_section, xc_section
142
143 CALL timeset(routinen, handle)
144
145 NULLIFY (atomic_kind_set, dispersion_env, ec_env%ls_env, para_env)
146 NULLIFY (ec_env%sab_orb, ec_env%sac_ae, ec_env%sac_ppl, ec_env%sap_ppnl)
147 NULLIFY (ec_env%matrix_ks, ec_env%matrix_h, ec_env%matrix_s)
148 NULLIFY (ec_env%matrix_t, ec_env%matrix_p, ec_env%matrix_w)
149 NULLIFY (ec_env%task_list)
150 NULLIFY (ec_env%mao_coef)
151 NULLIFY (ec_env%force)
152 NULLIFY (ec_env%dispersion_env)
153 NULLIFY (ec_env%xc_section)
154 NULLIFY (ec_env%matrix_z)
155 NULLIFY (ec_env%matrix_hz)
156 NULLIFY (ec_env%matrix_wz)
157 NULLIFY (ec_env%z_admm)
158 NULLIFY (ec_env%p_env)
159 NULLIFY (ec_env%vxc_rspace)
160 NULLIFY (ec_env%vtau_rspace)
161 NULLIFY (ec_env%vadmm_rspace)
162 NULLIFY (ec_env%rhoout_r, ec_env%rhoz_r)
163 NULLIFY (ec_env%x_data)
164 ec_env%should_update = .true.
165 ec_env%mao = .false.
166 ec_env%do_ec_admm = .false.
167 ec_env%do_kpoints = .false.
168 ec_env%do_ec_hfx = .false.
169 ec_env%reuse_hfx = .false.
170
171 IF (qs_env%energy_correction) THEN
172
173 cpassert(PRESENT(ec_section))
174
175 ! get a useful output_unit
176 unit_nr = cp_logger_get_default_unit_nr(local=.false.)
177
178 CALL section_vals_val_get(ec_section, "ALGORITHM", &
179 i_val=ec_env%ks_solver)
180 CALL section_vals_val_get(ec_section, "ENERGY_FUNCTIONAL", &
181 i_val=ec_env%energy_functional)
182 CALL section_vals_val_get(ec_section, "FACTORIZATION", &
183 i_val=ec_env%factorization)
184 CALL section_vals_val_get(ec_section, "OT_INITIAL_GUESS", &
185 i_val=ec_env%ec_initial_guess)
186 CALL section_vals_val_get(ec_section, "EPS_DEFAULT", &
187 r_val=ec_env%eps_default)
188 CALL section_vals_val_get(ec_section, "HARRIS_BASIS", &
189 c_val=ec_env%basis)
190 CALL section_vals_val_get(ec_section, "ELECTRONIC_TEMPERATURE", &
191 r_val=etemp)
192 CALL section_vals_val_get(ec_section, "MAO", &
193 l_val=ec_env%mao)
194 CALL section_vals_val_get(ec_section, "MAO_MAX_ITER", &
195 i_val=ec_env%mao_max_iter)
196 CALL section_vals_val_get(ec_section, "MAO_EPS_GRAD", &
197 r_val=ec_env%mao_eps_grad)
198 CALL section_vals_val_get(ec_section, "MAO_EPS1", &
199 r_val=ec_env%mao_eps1)
200 CALL section_vals_val_get(ec_section, "MAO_IOLEVEL", &
201 i_val=ec_env%mao_iolevel)
202 ! Skip EC calculation if ground-state calculation did not converge
203 CALL section_vals_val_get(ec_section, "SKIP_EC", &
204 l_val=ec_env%skip_ec)
205 ! Debug output
206 CALL section_vals_val_get(ec_section, "DEBUG_FORCES", &
207 l_val=ec_env%debug_forces)
208 CALL section_vals_val_get(ec_section, "DEBUG_STRESS", &
209 l_val=ec_env%debug_stress)
210 CALL section_vals_val_get(ec_section, "DEBUG_EXTERNAL_METHOD", &
211 l_val=ec_env%debug_external)
212 ! ADMM
213 CALL section_vals_val_get(ec_section, "ADMM", l_val=ec_env%do_ec_admm)
214 ! EXTERNAL
215 CALL section_vals_val_get(ec_section, "EXTERNAL_RESPONSE_FILENAME", &
216 c_val=ec_env%exresp_fn)
217 CALL section_vals_val_get(ec_section, "EXTERNAL_RESPONSE_ERROR_FILENAME", &
218 c_val=ec_env%exresperr_fn)
219 CALL section_vals_val_get(ec_section, "EXTERNAL_RESULT_FILENAME", &
220 c_val=ec_env%exresult_fn)
221 CALL section_vals_val_get(ec_section, "ERROR_ESTIMATION", &
222 l_val=ec_env%do_error)
223 CALL section_vals_val_get(ec_section, "ERROR_ESTIMATION_METHOD", &
224 c_val=ec_env%error_method)
225 CALL uppercase(ec_env%error_method)
226 CALL section_vals_val_get(ec_section, "ERROR_CUTOFF", &
227 r_val=ec_env%error_cutoff)
228 CALL section_vals_val_get(ec_section, "ERROR_SUBSPACE_SIZE", &
229 i_val=ec_env%error_subspace)
230
231 ec_env%do_skip = .false.
232
233 ! smearing
234 IF (etemp > 0.0_dp) THEN
235 ec_env%smear%do_smear = .true.
236 ec_env%smear%method = smear_fermi_dirac
237 ec_env%smear%electronic_temperature = etemp
238 ec_env%smear%eps_fermi_dirac = 1.0e-5_dp
239 ec_env%smear%fixed_mag_mom = -100.0_dp
240 ELSE
241 ec_env%smear%do_smear = .false.
242 ec_env%smear%electronic_temperature = 0.0_dp
243 END IF
244
245 ! kpoints
246 ! Options: gs ec
247 ! gamma gamma
248 ! gamma KPec
249 ! KPgs KPgs
250 ! KPgs KPec
251 CALL get_qs_env(qs_env, do_kpoints=gs_kpoints)
252 kp_section => section_vals_get_subs_vals(ec_section, "KPOINTS")
253 CALL section_vals_get(kp_section, explicit=explicit)
254 ec_env%do_kpoints = gs_kpoints .OR. explicit
255 IF (ec_env%do_kpoints) THEN
256 IF (.NOT. explicit) THEN
257 kp_section => section_vals_get_subs_vals(dft_section, "KPOINTS")
258 END IF
259 CALL get_qs_env(qs_env, particle_set=particle_set, cell=cell)
260 CALL kpoint_create(ec_env%kpoints)
261 CALL read_kpoint_section(ec_env%kpoints, kp_section, cell%hmat)
262 CALL kpoint_initialize(ec_env%kpoints, particle_set, cell)
263 CALL get_qs_env(qs_env, para_env=para_env, blacs_env=blacs_env)
264 CALL kpoint_env_initialize(ec_env%kpoints, para_env, blacs_env)
265 ELSE
266 NULLIFY (ec_env%kpoints)
267 END IF
268
269 ! set basis
270 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, nkind=nkind)
271 CALL uppercase(ec_env%basis)
272 SELECT CASE (ec_env%basis)
273 CASE ("ORBITAL")
274 DO ikind = 1, nkind
275 qs_kind => qs_kind_set(ikind)
276 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
277 IF (ASSOCIATED(basis_set)) THEN
278 NULLIFY (harris_basis)
279 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
280 IF (ASSOCIATED(harris_basis)) THEN
281 CALL remove_basis_from_container(qs_kind%basis_sets, basis_type="HARRIS")
282 END IF
283 NULLIFY (harris_basis)
284 CALL copy_gto_basis_set(basis_set, harris_basis)
285 CALL add_basis_set_to_container(qs_kind%basis_sets, harris_basis, "HARRIS")
286 END IF
287 END DO
288 CASE ("PRIMITIVE")
289 DO ikind = 1, nkind
290 qs_kind => qs_kind_set(ikind)
291 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
292 IF (ASSOCIATED(basis_set)) THEN
293 NULLIFY (harris_basis)
294 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
295 IF (ASSOCIATED(harris_basis)) THEN
296 CALL remove_basis_from_container(qs_kind%basis_sets, basis_type="HARRIS")
297 END IF
298 NULLIFY (harris_basis)
299 CALL create_primitive_basis_set(basis_set, harris_basis)
300 CALL get_qs_env(qs_env, dft_control=dft_control)
301 eps_pgf_orb = dft_control%qs_control%eps_pgf_orb
302 CALL init_interaction_radii_orb_basis(harris_basis, eps_pgf_orb)
303 harris_basis%kind_radius = basis_set%kind_radius
304 CALL add_basis_set_to_container(qs_kind%basis_sets, harris_basis, "HARRIS")
305 END IF
306 END DO
307 CASE ("HARRIS")
308 DO ikind = 1, nkind
309 qs_kind => qs_kind_set(ikind)
310 NULLIFY (harris_basis)
311 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
312 IF (.NOT. ASSOCIATED(harris_basis)) THEN
313 cpwarn("Harris Basis not defined for all types of atoms.")
314 END IF
315 END DO
316 CASE DEFAULT
317 cpabort("Unknown basis set for energy correction (Harris functional)")
318 END SELECT
319 !
320 CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, basis_type="HARRIS")
321 CALL init_orbital_pointers(maxlgto + 1)
322 ! GAPW: Generate soft version of Harris basis
323 CALL get_qs_env(qs_env, dft_control=dft_control)
324 IF (dft_control%qs_control%gapw .OR. dft_control%qs_control%gapw_xc) THEN
325 eps_pgf_orb = dft_control%qs_control%eps_pgf_orb
326 DO ikind = 1, nkind
327 qs_kind => qs_kind_set(ikind)
328 NULLIFY (harris_basis)
329 CALL get_qs_kind(qs_kind, basis_set=harris_basis, basis_type="HARRIS")
330 CALL get_qs_kind(qs_kind, hard_radius=rc, gpw_type_forced=gpw)
331 NULLIFY (harris_soft_basis)
332 CALL allocate_gto_basis_set(harris_soft_basis)
333 CALL create_soft_basis(harris_basis, harris_soft_basis, &
334 dft_control%qs_control%gapw_control%eps_fit, &
335 rc, paw_atom, &
336 dft_control%qs_control%gapw_control%force_paw, gpw)
337 CALL init_interaction_radii_orb_basis(harris_soft_basis, eps_pgf_orb)
338 CALL add_basis_set_to_container(qs_kind%basis_sets, harris_soft_basis, "HARRIS_SOFT")
339 END DO
340 END IF
341 !
342 CALL uppercase(ec_env%basis)
343
344 ! Basis may only differ from ground-state if explicitly added
345 ec_env%basis_inconsistent = .false.
346 IF (ec_env%basis == "HARRIS") THEN
347 DO ikind = 1, nkind
348 qs_kind => qs_kind_set(ikind)
349 ! Basis sets of ground-state
350 CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
351 ! Basis sets of energy correction
352 CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
353
354 IF (basis_set%name /= harris_basis%name) THEN
355 ec_env%basis_inconsistent = .true.
356 END IF
357 END DO
358 END IF
359
360 !Density-corrected DFT must be performed with the same basis as ground-state
361 IF (ec_env%energy_functional == ec_functional_dc .AND. ec_env%basis_inconsistent) THEN
362 CALL cp_abort(__location__, &
363 "DC-DFT: Correction and ground state need to use the same basis. "// &
364 "Checked by comparing basis set names only.")
365 END IF
366 IF (ec_env%energy_functional == ec_functional_ext .AND. ec_env%basis_inconsistent) THEN
367 CALL cp_abort(__location__, &
368 "Exteranl Energy: Correction and ground state need to use the same basis. "// &
369 "Checked by comparing basis set names only.")
370 END IF
371 !
372 ! set functional
373 SELECT CASE (ec_env%energy_functional)
375 ec_env%ec_name = "Harris"
376 CASE (ec_functional_dc)
377 ec_env%ec_name = "DC-DFT"
378 CASE (ec_functional_ext)
379 ec_env%ec_name = "External Energy"
380 CASE DEFAULT
381 cpabort("unknown energy correction")
382 END SELECT
383 ! select the XC section
384 NULLIFY (xc_section)
385 xc_section => section_vals_get_subs_vals(dft_section, "XC")
386 section1 => section_vals_get_subs_vals(ec_section, "XC")
387 section2 => section_vals_get_subs_vals(ec_section, "XC%XC_FUNCTIONAL")
388 CALL section_vals_get(section2, explicit=explicit)
389 IF (explicit) THEN
390 CALL xc_functionals_expand(section2, section1)
391 ec_env%xc_section => section1
392 ELSE
393 ec_env%xc_section => xc_section
394 END IF
395 ! Check whether energy correction requires the kinetic energy density and rebuild rho if necessary
396 CALL get_qs_env(qs_env, dft_control=dft_control, rho=rho)
397 xc_fun_section => section_vals_get_subs_vals(ec_env%xc_section, "XC_FUNCTIONAL")
398 dft_control%use_kinetic_energy_density = dft_control%use_kinetic_energy_density .OR. &
399 xc_uses_kinetic_energy_density(xc_fun_section, dft_control%lsd)
400 ! Same for density gradient
401 dft_control%drho_by_collocation = dft_control%drho_by_collocation .OR. &
402 (xc_uses_norm_drho(xc_fun_section, dft_control%lsd) .AND. &
403 (section_get_ival(xc_section, "XC_GRID%XC_DERIV") == xc_deriv_collocate))
404 ! dispersion
405 ALLOCATE (dispersion_env)
406 NULLIFY (xc_section)
407 xc_section => ec_env%xc_section
408 CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, para_env=para_env)
409 CALL qs_dispersion_env_set(dispersion_env, xc_section)
410 IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
411 NULLIFY (pp_section)
412 pp_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%PAIR_POTENTIAL")
413 CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
414 ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
415 cpabort("nl-vdW functionals not available for EC calculations")
416 NULLIFY (nl_section)
417 nl_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%NON_LOCAL")
418 CALL qs_dispersion_nonloc_init(dispersion_env, para_env)
419 END IF
420 ec_env%dispersion_env => dispersion_env
421
422 ! Check if hybrid functional are used
423 ec_hfx_section => section_vals_get_subs_vals(ec_section, "XC%HF")
424 CALL section_vals_get(ec_hfx_section, explicit=ec_env%do_ec_hfx)
425
426 ! Initialize Harris LS solver environment
427 ec_env%use_ls_solver = .false.
428 ec_env%use_ls_solver = (ec_env%ks_solver == ec_matrix_sign) &
429 .OR. (ec_env%ks_solver == ec_matrix_trs4) &
430 .OR. (ec_env%ks_solver == ec_matrix_tc2)
431
432 IF (ec_env%use_ls_solver) THEN
433 CALL ec_ls_create(qs_env, ec_env)
434 END IF
435
436 ! check that Harris functional with electronic temperature uses diagonalization
437 IF (ec_env%energy_functional == ec_functional_harris) THEN
438 IF (ec_env%smear%do_smear .AND. ec_env%ks_solver /= ec_diagonalization) THEN
439 cpabort("Harris functional with Fermi-Dirac smearing needs diagonalization solver.")
440 END IF
441 END IF
442
443 END IF
444
445 CALL timestop(handle)
446
447 END SUBROUTINE init_ec_env
448
449! **************************************************************************************************
450!> \brief Initializes linear scaling environment for LS based solver of
451!> Harris energy functional and parses input section
452!> \param qs_env ...
453!> \param ec_env ...
454!> \par History
455!> 2020.10 created [Fabian Belleflamme]
456!> \author Fabian Belleflamme
457! **************************************************************************************************
458 SUBROUTINE ec_ls_create(qs_env, ec_env)
459 TYPE(qs_environment_type), POINTER :: qs_env
460 TYPE(energy_correction_type), POINTER :: ec_env
461
462 CHARACTER(LEN=*), PARAMETER :: routinen = 'ec_ls_create'
463
464 INTEGER :: handle
465 REAL(kind=dp) :: mu
466 TYPE(dft_control_type), POINTER :: dft_control
467 TYPE(ls_scf_env_type), POINTER :: ls_env
468 TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
469 TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
470 TYPE(section_vals_type), POINTER :: ec_section, input
471
472 CALL timeset(routinen, handle)
473
474 ALLOCATE (ec_env%ls_env)
475 ls_env => ec_env%ls_env
476
477 NULLIFY (dft_control, input, ls_env%para_env)
478
479 CALL get_qs_env(qs_env, &
480 dft_control=dft_control, &
481 input=input, &
482 molecule_set=molecule_set, &
483 particle_set=particle_set, &
484 para_env=ls_env%para_env, &
485 nelectron_spin=ls_env%nelectron_spin)
486
487 ! copy some basic stuff
488 ls_env%nspins = dft_control%nspins
489 ls_env%natoms = SIZE(particle_set, 1)
490 CALL ls_env%para_env%retain()
491
492 ! initialize block to group to defined molecules
493 ALLOCATE (ls_env%ls_mstruct%atom_to_molecule(ls_env%natoms))
494 CALL molecule_of_atom(molecule_set, atom_to_mol=ls_env%ls_mstruct%atom_to_molecule)
495
496 ls_env%do_transport = .false.
497 ls_env%do_pao = .false.
498 ls_env%ls_mstruct%do_pao = ls_env%do_pao
499 ls_env%do_pexsi = .false.
500 ls_env%has_unit_metric = .false.
501
502 ec_section => section_vals_get_subs_vals(input, "DFT%ENERGY_CORRECTION")
503 CALL section_vals_val_get(ec_section, "EPS_FILTER", r_val=ls_env%eps_filter)
504 CALL section_vals_val_get(ec_section, "MU", r_val=mu)
505 CALL section_vals_val_get(ec_section, "FIXED_MU", l_val=ls_env%fixed_mu)
506 ls_env%mu_spin = mu
507 CALL section_vals_val_get(ec_section, "S_PRECONDITIONER", i_val=ls_env%s_preconditioner_type)
508 CALL section_vals_val_get(ec_section, "MATRIX_CLUSTER_TYPE", i_val=ls_env%ls_mstruct%cluster_type)
509 CALL section_vals_val_get(ec_section, "S_INVERSION", i_val=ls_env%s_inversion_type)
510 CALL section_vals_val_get(ec_section, "CHECK_S_INV", l_val=ls_env%check_s_inv)
511 CALL section_vals_val_get(ec_section, "REPORT_ALL_SPARSITIES", l_val=ls_env%report_all_sparsities)
512 CALL section_vals_val_get(ec_section, "SIGN_METHOD", i_val=ls_env%sign_method)
513 CALL section_vals_val_get(ec_section, "SIGN_ORDER", i_val=ls_env%sign_order)
514 CALL section_vals_val_get(ec_section, "DYNAMIC_THRESHOLD", l_val=ls_env%dynamic_threshold)
515 CALL section_vals_val_get(ec_section, "NON_MONOTONIC", l_val=ls_env%non_monotonic)
516 CALL section_vals_val_get(ec_section, "S_SQRT_METHOD", i_val=ls_env%s_sqrt_method)
517 CALL section_vals_val_get(ec_section, "S_SQRT_ORDER", i_val=ls_env%s_sqrt_order)
518 CALL section_vals_val_get(ec_section, "EPS_LANCZOS", r_val=ls_env%eps_lanczos)
519 CALL section_vals_val_get(ec_section, "MAX_ITER_LANCZOS", i_val=ls_env%max_iter_lanczos)
520
521 SELECT CASE (ec_env%ks_solver)
522 CASE (ec_matrix_sign)
523 ! S inverse required for Sign matrix algorithm,
524 ! calculated either by Hotelling or multiplying S matrix sqrt inv
525 SELECT CASE (ls_env%s_inversion_type)
527 ls_env%needs_s_inv = .true.
528 ls_env%use_s_sqrt = .true.
530 ls_env%needs_s_inv = .true.
531 ls_env%use_s_sqrt = .false.
533 ls_env%needs_s_inv = .false.
534 ls_env%use_s_sqrt = .false.
535 CASE DEFAULT
536 cpabort("")
537 END SELECT
539 ls_env%needs_s_inv = .false.
540 ls_env%use_s_sqrt = .true.
541 CASE DEFAULT
542 cpabort("")
543 END SELECT
544
545 SELECT CASE (ls_env%s_preconditioner_type)
547 ls_env%has_s_preconditioner = .false.
548 CASE DEFAULT
549 ls_env%has_s_preconditioner = .true.
550 END SELECT
551
552 ! buffer for the history of matrices, not needed here
553 ls_env%extrapolation_order = 0
554 ls_env%scf_history%nstore = 0
555 ls_env%scf_history%istore = 0
556 ALLOCATE (ls_env%scf_history%matrix(ls_env%nspins, ls_env%scf_history%nstore))
557
558 NULLIFY (ls_env%mixing_store)
559
560 CALL timestop(handle)
561
562 END SUBROUTINE ec_ls_create
563
564! **************************************************************************************************
565!> \brief Print out the energy correction input section
566!>
567!> \param ec_env ...
568!> \par History
569!> 2020.10 created [Fabian Belleflamme]
570!> \author Fabian Belleflamme
571! **************************************************************************************************
572 SUBROUTINE ec_write_input(ec_env)
573 TYPE(energy_correction_type), POINTER :: ec_env
574
575 CHARACTER(LEN=*), PARAMETER :: routinen = 'ec_write_input'
576
577 INTEGER :: handle, unit_nr
578 TYPE(ls_scf_env_type), POINTER :: ls_env
579
580 CALL timeset(routinen, handle)
581
582 unit_nr = cp_logger_get_default_unit_nr(local=.false.)
583
584 IF (unit_nr > 0) THEN
585
586 WRITE (unit_nr, '(T2,A)') &
587 "!"//repeat("-", 29)//" Energy Correction "//repeat("-", 29)//"!"
588
589 ! Type of energy correction
590 SELECT CASE (ec_env%energy_functional)
592 WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "HARRIS FUNCTIONAL"
593 CASE (ec_functional_dc)
594 WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "DC-DFT"
595 CASE (ec_functional_ext)
596 WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "External"
597 END SELECT
598 WRITE (unit_nr, '()')
599
600 ! Energy correction parameters
601 WRITE (unit_nr, '(T2,A,T61,E20.3)') "eps_default:", ec_env%eps_default
602
603 CALL uppercase(ec_env%basis)
604 SELECT CASE (ec_env%basis)
605 CASE ("ORBITAL")
606 WRITE (unit_nr, '(T2,A,T61,A20)') "EC basis: ", "ORBITAL"
607 CASE ("PRIMITIVE")
608 WRITE (unit_nr, '(T2,A,T61,A20)') "EC basis: ", "PRIMITIVE"
609 CASE ("HARRIS")
610 WRITE (unit_nr, '(T2,A,T61,A20)') "EC Basis: ", "HARRIS"
611 END SELECT
612
613 ! Info how HFX in energy correction is treated
614 IF (ec_env%do_ec_hfx) THEN
615
616 WRITE (unit_nr, '(T2,A,T61,L20)') "DC-DFT with HFX", ec_env%do_ec_hfx
617 WRITE (unit_nr, '(T2,A,T61,L20)') "Reuse HFX integrals", ec_env%reuse_hfx
618 WRITE (unit_nr, '(T2,A,T61,L20)') "DC-DFT HFX with ADMM", ec_env%do_ec_admm
619
620 END IF ! ec_env%do_ec_hfx
621
622 IF (ec_env%do_kpoints) THEN
623 CALL write_kpoint_info(ec_env%kpoints, iounit=unit_nr)
624 END IF
625
626 ! Parameters for Harris functional solver
627 IF (ec_env%energy_functional == ec_functional_harris) THEN
628
629 ! Algorithm
630 SELECT CASE (ec_env%ks_solver)
631 CASE (ec_diagonalization)
632 WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "DIAGONALIZATION"
633 CASE (ec_ot_diag)
634 WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "OT DIAGONALIZATION"
635 CASE (ec_matrix_sign)
636 WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "MATRIX_SIGN"
637 CASE (ec_matrix_trs4)
638 WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "TRS4"
639 CALL cite_reference(niklasson2003)
640 CASE (ec_matrix_tc2)
641 WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "TC2"
642 CALL cite_reference(niklasson2014)
643 END SELECT
644 WRITE (unit_nr, '()')
645
646 ! MAO
647 IF (ec_env%mao) THEN
648 WRITE (unit_nr, '(T2,A,T61,L20)') "MAO:", ec_env%mao
649 WRITE (unit_nr, '(T2,A,T61,L20)') "MAO_IOLEVEL:", ec_env%mao_iolevel
650 WRITE (unit_nr, '(T2,A,T61,I20)') "MAO_MAX_ITER:", ec_env%mao_max_iter
651 WRITE (unit_nr, '(T2,A,T61,E20.3)') "MAO_EPS_GRAD:", ec_env%mao_eps_grad
652 WRITE (unit_nr, '(T2,A,T61,E20.3)') "MAO_EPS1:", ec_env%mao_eps1
653 WRITE (unit_nr, '()')
654 END IF
655
656 ! Parameters for linear response solver
657 IF (.NOT. ec_env%use_ls_solver) THEN
658
659 WRITE (unit_nr, '(T2,A)') "MO Solver"
660 WRITE (unit_nr, '()')
661
662 SELECT CASE (ec_env%ks_solver)
663 CASE (ec_diagonalization)
664
665 SELECT CASE (ec_env%factorization)
666 CASE (kg_cholesky)
667 WRITE (unit_nr, '(T2,A,T61,A20)') "Factorization: ", "CHOLESKY"
668 END SELECT
669
670 CASE (ec_ot_diag)
671
672 ! OT Diagonalization
673 ! Initial guess : 1) block diagonal initial guess
674 ! 2) GS-density matrix (might require trafo if basis diff)
675
676 SELECT CASE (ec_env%ec_initial_guess)
677 CASE (ec_ot_atomic)
678 WRITE (unit_nr, '(T2,A,T61,A20)') "OT Diag initial guess: ", "ATOMIC"
679 CASE (ec_ot_gs)
680 WRITE (unit_nr, '(T2,A,T61,A20)') "OT Diag initial guess: ", "GROUND STATE DM"
681 END SELECT
682
683 CASE DEFAULT
684 cpabort("Unknown Diagonalization algorithm for Harris functional")
685 END SELECT
686
687 ELSE
688
689 WRITE (unit_nr, '(T2,A)') "AO Solver"
690 WRITE (unit_nr, '()')
691
692 ls_env => ec_env%ls_env
693 WRITE (unit_nr, '(T2,A,T61,E20.3)') "eps_filter:", ls_env%eps_filter
694 WRITE (unit_nr, '(T2,A,T61,L20)') "fixed chemical potential (mu)", ls_env%fixed_mu
695 WRITE (unit_nr, '(T2,A,T61,L20)') "Computing inv(S):", ls_env%needs_s_inv
696 WRITE (unit_nr, '(T2,A,T61,L20)') "Computing sqrt(S):", ls_env%use_s_sqrt
697 WRITE (unit_nr, '(T2,A,T61,L20)') "Computing S preconditioner ", ls_env%has_s_preconditioner
698
699 IF (ls_env%use_s_sqrt) THEN
700 SELECT CASE (ls_env%s_sqrt_method)
701 CASE (ls_s_sqrt_ns)
702 WRITE (unit_nr, '(T2,A,T61,A20)') "S sqrt method:", "NEWTONSCHULZ"
703 CASE (ls_s_sqrt_proot)
704 WRITE (unit_nr, '(T2,A,T61,A20)') "S sqrt method:", "PROOT"
705 CASE DEFAULT
706 cpabort("Unknown sqrt method.")
707 END SELECT
708 WRITE (unit_nr, '(T2,A,T61,I20)') "S sqrt order:", ls_env%s_sqrt_order
709 END IF
710
711 SELECT CASE (ls_env%s_preconditioner_type)
713 WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "NONE"
715 WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "ATOMIC"
717 WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "MOLECULAR"
718 END SELECT
719
720 SELECT CASE (ls_env%ls_mstruct%cluster_type)
721 CASE (ls_cluster_atomic)
722 WRITE (unit_nr, '(T2,A,T61,A20)') "Cluster type", adjustr("ATOMIC")
724 WRITE (unit_nr, '(T2,A,T61,A20)') "Cluster type", adjustr("MOLECULAR")
725 CASE DEFAULT
726 cpabort("Unknown cluster type")
727 END SELECT
728
729 END IF
730
731 END IF ! if ec_functional_harris
732
733 WRITE (unit_nr, '(T2,A)') repeat("-", 79)
734 WRITE (unit_nr, '()')
735
736 END IF ! unit_nr
737
738 CALL timestop(handle)
739
740 END SUBROUTINE ec_write_input
741
742END MODULE ec_environment
Define the atomic kind types and their sub types.
subroutine, public remove_basis_from_container(container, inum, basis_type)
...
subroutine, public add_basis_set_to_container(container, basis_set, basis_set_type)
...
subroutine, public allocate_gto_basis_set(gto_basis_set)
...
subroutine, public copy_gto_basis_set(basis_set_in, basis_set_out)
...
subroutine, public create_primitive_basis_set(basis_set, pbasis, lmax)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public niklasson2014
integer, save, public niklasson2003
Handles all functions related to the CELL.
Definition cell_types.F:15
methods related to the blacs parallel environment
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
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
Types needed for a linear scaling quickstep SCF run based on the density matrix.
Types needed for a for a Energy Correction.
Energy correction environment setup and handling.
subroutine, public ec_write_input(ec_env)
Print out the energy correction input section.
subroutine, public ec_env_create(qs_env, ec_env, dft_section, ec_section)
Allocates and intitializes ec_env.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public ec_functional_harris
integer, parameter, public smear_fermi_dirac
integer, parameter, public ec_functional_dc
integer, parameter, public ls_s_preconditioner_molecular
integer, parameter, public ls_s_inversion_hotelling
integer, parameter, public xc_vdw_fun_nonloc
integer, parameter, public ls_cluster_molecular
integer, parameter, public ec_ot_atomic
integer, parameter, public ls_s_preconditioner_atomic
integer, parameter, public ec_ot_diag
integer, parameter, public ec_ot_gs
integer, parameter, public ls_s_inversion_none
integer, parameter, public ls_s_sqrt_proot
integer, parameter, public ls_s_sqrt_ns
integer, parameter, public ls_s_preconditioner_none
integer, parameter, public ec_matrix_tc2
integer, parameter, public ec_diagonalization
integer, parameter, public ec_matrix_trs4
integer, parameter, public ec_functional_ext
integer, parameter, public ec_matrix_sign
integer, parameter, public ls_s_inversion_sign_sqrt
integer, parameter, public xc_vdw_fun_pairpot
integer, parameter, public kg_cholesky
integer, parameter, public ls_cluster_atomic
checks the input and perform some automatic "magic" on it
subroutine, public xc_functionals_expand(functionals, xc_section)
expand a shortcutted functional section
objects that represent the structure of input sections and the data contained in an input section
integer function, public section_get_ival(section_vals, keyword_name)
...
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
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
Routines needed for kpoint calculation.
subroutine, public kpoint_initialize(kpoint, particle_set, cell)
Generate the kpoints and initialize the kpoint environment.
subroutine, public kpoint_env_initialize(kpoint, para_env, blacs_env, with_aux_fit)
Initialize the kpoint environment.
Types and basic routines needed for a kpoint calculation.
subroutine, public read_kpoint_section(kpoint, kpoint_section, a_vec)
Read the kpoint input section.
subroutine, public write_kpoint_info(kpoint, iounit, dft_section)
Write information on the kpoints to output.
subroutine, public kpoint_create(kpoint)
Create a kpoint environment.
Interface to the message passing library MPI.
Define the data structure for the molecule information.
subroutine, public molecule_of_atom(molecule_set, atom_to_mol)
finds for each atom the molecule it belongs to
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
Define the data structure for the particle information.
Calculation of non local dispersion functionals Some routines adapted from: Copyright (C) 2001-2009 Q...
subroutine, public qs_dispersion_nonloc_init(dispersion_env, para_env)
...
Calculation of dispersion using pair potentials.
subroutine, public qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
...
Definition of disperson types for DFT calculations.
Set disperson types for DFT calculations.
subroutine, public qs_dispersion_env_set(dispersion_env, xc_section)
...
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.
Calculate the interaction radii for the operator matrix calculation.
subroutine, public init_interaction_radii_orb_basis(orb_basis_set, eps_pgf_orb, eps_pgf_short)
...
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, cneo_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, monovalent, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg, cneo_potential_present, nkind_q, natom_q)
Get attributes of an atomic kind set.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public create_soft_basis(orb_basis, soft_basis, eps_fit, rc, paw_atom, paw_type_forced, gpw_r3d_rs_type_forced)
create the soft basis from a GTO basis
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
input constants for xc
integer, parameter, public xc_deriv_collocate
Exchange and Correlation functional calculations.
Definition xc.F:17
logical function, public xc_uses_norm_drho(xc_fun_section, lsd)
...
Definition xc.F:111
logical function, public xc_uses_kinetic_energy_density(xc_fun_section, lsd)
...
Definition xc.F:91
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
Definition cell_types.F:60
represent a blacs multidimensional parallel environment (for the mpi corrispective see cp_paratypes/m...
Contains information on the energy correction functional for KG.
stores all the informations relevant to an mpi environment
Provides all information about a quickstep kind.
keeps the density in various representations, keeping track of which ones are valid.