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 Energy correction environment setup and handling
10 : !> \par History
11 : !> 2019.09 created
12 : !> \author JGH
13 : ! **************************************************************************************************
14 : MODULE ec_environment
15 : USE atomic_kind_types, ONLY: atomic_kind_type
16 : USE basis_set_container_types, ONLY: add_basis_set_to_container,&
17 : remove_basis_from_container
18 : USE basis_set_types, ONLY: copy_gto_basis_set,&
19 : create_primitive_basis_set,&
20 : gto_basis_set_type
21 : USE bibliography, ONLY: Niklasson2003,&
22 : Niklasson2014,&
23 : cite_reference
24 : USE cp_control_types, ONLY: dft_control_type
25 : USE cp_log_handling, ONLY: cp_get_default_logger,&
26 : cp_logger_get_default_unit_nr,&
27 : cp_logger_type
28 : USE dm_ls_scf_types, ONLY: ls_scf_env_type
29 : USE ec_env_types, ONLY: energy_correction_type
30 : USE input_constants, ONLY: &
31 : ec_diagonalization, ec_functional_dc, ec_functional_ext, ec_functional_harris, &
32 : ec_matrix_sign, ec_matrix_tc2, ec_matrix_trs4, ec_ot_atomic, ec_ot_diag, ec_ot_gs, &
33 : kg_cholesky, ls_cluster_atomic, ls_cluster_molecular, ls_s_inversion_hotelling, &
34 : ls_s_inversion_none, ls_s_inversion_sign_sqrt, ls_s_preconditioner_atomic, &
35 : ls_s_preconditioner_molecular, ls_s_preconditioner_none, ls_s_sqrt_ns, ls_s_sqrt_proot, &
36 : xc_vdw_fun_nonloc, xc_vdw_fun_pairpot
37 : USE input_cp2k_check, ONLY: xc_functionals_expand
38 : USE input_section_types, ONLY: section_get_ival,&
39 : section_vals_get,&
40 : section_vals_get_subs_vals,&
41 : section_vals_type,&
42 : section_vals_val_get
43 : USE kinds, ONLY: dp
44 : USE message_passing, ONLY: mp_para_env_type
45 : USE molecule_types, ONLY: molecule_of_atom,&
46 : molecule_type
47 : USE orbital_pointers, ONLY: init_orbital_pointers
48 : USE particle_types, ONLY: particle_type
49 : USE qs_dispersion_nonloc, ONLY: qs_dispersion_nonloc_init
50 : USE qs_dispersion_pairpot, ONLY: qs_dispersion_pairpot_init
51 : USE qs_dispersion_types, ONLY: qs_dispersion_type
52 : USE qs_dispersion_utils, ONLY: qs_dispersion_env_set
53 : USE qs_environment_types, ONLY: get_qs_env,&
54 : qs_environment_type
55 : USE qs_interactions, ONLY: init_interaction_radii_orb_basis
56 : USE qs_kind_types, ONLY: get_qs_kind,&
57 : get_qs_kind_set,&
58 : qs_kind_type
59 : USE qs_rho_types, ONLY: qs_rho_type
60 : USE string_utilities, ONLY: uppercase
61 : USE xc, ONLY: xc_uses_kinetic_energy_density,&
62 : xc_uses_norm_drho
63 : USE xc_input_constants, ONLY: xc_deriv_collocate
64 : #include "./base/base_uses.f90"
65 :
66 : IMPLICIT NONE
67 :
68 : PRIVATE
69 :
70 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ec_environment'
71 :
72 : PUBLIC :: ec_env_create
73 : PUBLIC :: ec_write_input
74 :
75 : CONTAINS
76 :
77 : ! **************************************************************************************************
78 : !> \brief Allocates and intitializes ec_env
79 : !> \param qs_env The QS environment
80 : !> \param ec_env The energy correction environment (the object to create)
81 : !> \param dft_section The DFT section
82 : !> \param ec_section The energy correction input section
83 : !> \par History
84 : !> 2019.09 created
85 : !> \author JGH
86 : ! **************************************************************************************************
87 6686 : SUBROUTINE ec_env_create(qs_env, ec_env, dft_section, ec_section)
88 : TYPE(qs_environment_type), POINTER :: qs_env
89 : TYPE(energy_correction_type), POINTER :: ec_env
90 : TYPE(section_vals_type), POINTER :: dft_section
91 : TYPE(section_vals_type), OPTIONAL, POINTER :: ec_section
92 :
93 6686 : CPASSERT(.NOT. ASSOCIATED(ec_env))
94 6686 : ALLOCATE (ec_env)
95 6686 : CALL init_ec_env(qs_env, ec_env, dft_section, ec_section)
96 :
97 6686 : END SUBROUTINE ec_env_create
98 :
99 : ! **************************************************************************************************
100 : !> \brief Initializes energy correction environment
101 : !> \param qs_env The QS environment
102 : !> \param ec_env The energy correction environment
103 : !> \param dft_section The DFT section
104 : !> \param ec_section The energy correction input section
105 : !> \par History
106 : !> 2019.09 created
107 : !> \author JGH
108 : ! **************************************************************************************************
109 6686 : SUBROUTINE init_ec_env(qs_env, ec_env, dft_section, ec_section)
110 : TYPE(qs_environment_type), POINTER :: qs_env
111 : TYPE(energy_correction_type), POINTER :: ec_env
112 : TYPE(section_vals_type), POINTER :: dft_section
113 : TYPE(section_vals_type), OPTIONAL, POINTER :: ec_section
114 :
115 : CHARACTER(LEN=*), PARAMETER :: routineN = 'init_ec_env'
116 :
117 : INTEGER :: handle, ikind, maxlgto, nkind, unit_nr
118 : LOGICAL :: explicit
119 : REAL(KIND=dp) :: eps_pgf_orb
120 6686 : TYPE(atomic_kind_type), DIMENSION(:), POINTER :: atomic_kind_set
121 : TYPE(cp_logger_type), POINTER :: logger
122 : TYPE(dft_control_type), POINTER :: dft_control
123 : TYPE(gto_basis_set_type), POINTER :: basis_set, harris_basis
124 : TYPE(mp_para_env_type), POINTER :: para_env
125 : TYPE(qs_dispersion_type), POINTER :: dispersion_env
126 6686 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
127 : TYPE(qs_kind_type), POINTER :: qs_kind
128 : TYPE(qs_rho_type), POINTER :: rho
129 : TYPE(section_vals_type), POINTER :: ec_hfx_section, nl_section, pp_section, &
130 : section1, section2, xc_fun_section, &
131 : xc_section
132 :
133 6686 : CALL timeset(routineN, handle)
134 :
135 6686 : NULLIFY (atomic_kind_set, dispersion_env, ec_env%ls_env, para_env)
136 6686 : NULLIFY (ec_env%sab_orb, ec_env%sac_ppl, ec_env%sap_ppnl)
137 6686 : NULLIFY (ec_env%matrix_ks, ec_env%matrix_h, ec_env%matrix_s)
138 6686 : NULLIFY (ec_env%matrix_t, ec_env%matrix_p, ec_env%matrix_w)
139 6686 : NULLIFY (ec_env%task_list)
140 6686 : NULLIFY (ec_env%mao_coef)
141 6686 : NULLIFY (ec_env%force)
142 6686 : NULLIFY (ec_env%dispersion_env)
143 6686 : NULLIFY (ec_env%xc_section)
144 6686 : NULLIFY (ec_env%matrix_z)
145 6686 : NULLIFY (ec_env%matrix_hz)
146 6686 : NULLIFY (ec_env%matrix_wz)
147 6686 : NULLIFY (ec_env%z_admm)
148 6686 : NULLIFY (ec_env%p_env)
149 6686 : NULLIFY (ec_env%vxc_rspace)
150 6686 : NULLIFY (ec_env%vtau_rspace)
151 6686 : NULLIFY (ec_env%vadmm_rspace)
152 6686 : NULLIFY (ec_env%rhoout_r, ec_env%rhoz_r)
153 6686 : NULLIFY (ec_env%x_data)
154 6686 : ec_env%should_update = .TRUE.
155 6686 : ec_env%mao = .FALSE.
156 6686 : ec_env%do_ec_admm = .FALSE.
157 6686 : ec_env%do_ec_hfx = .FALSE.
158 6686 : ec_env%reuse_hfx = .FALSE.
159 :
160 6686 : IF (qs_env%energy_correction) THEN
161 :
162 236 : CPASSERT(PRESENT(ec_section))
163 : ! get a useful output_unit
164 236 : logger => cp_get_default_logger()
165 236 : IF (logger%para_env%is_source()) THEN
166 118 : unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
167 : ELSE
168 : unit_nr = -1
169 : END IF
170 :
171 : CALL section_vals_val_get(ec_section, "ALGORITHM", &
172 236 : i_val=ec_env%ks_solver)
173 : CALL section_vals_val_get(ec_section, "ENERGY_FUNCTIONAL", &
174 236 : i_val=ec_env%energy_functional)
175 : CALL section_vals_val_get(ec_section, "FACTORIZATION", &
176 236 : i_val=ec_env%factorization)
177 : CALL section_vals_val_get(ec_section, "OT_INITIAL_GUESS", &
178 236 : i_val=ec_env%ec_initial_guess)
179 : CALL section_vals_val_get(ec_section, "EPS_DEFAULT", &
180 236 : r_val=ec_env%eps_default)
181 : CALL section_vals_val_get(ec_section, "HARRIS_BASIS", &
182 236 : c_val=ec_env%basis)
183 : CALL section_vals_val_get(ec_section, "MAO", &
184 236 : l_val=ec_env%mao)
185 : CALL section_vals_val_get(ec_section, "MAO_MAX_ITER", &
186 236 : i_val=ec_env%mao_max_iter)
187 : CALL section_vals_val_get(ec_section, "MAO_EPS_GRAD", &
188 236 : r_val=ec_env%mao_eps_grad)
189 : CALL section_vals_val_get(ec_section, "MAO_EPS1", &
190 236 : r_val=ec_env%mao_eps1)
191 : CALL section_vals_val_get(ec_section, "MAO_IOLEVEL", &
192 236 : i_val=ec_env%mao_iolevel)
193 : ! Skip EC calculation if ground-state calculation did not converge
194 : CALL section_vals_val_get(ec_section, "SKIP_EC", &
195 236 : l_val=ec_env%skip_ec)
196 : ! Debug output
197 : CALL section_vals_val_get(ec_section, "DEBUG_FORCES", &
198 236 : l_val=ec_env%debug_forces)
199 : CALL section_vals_val_get(ec_section, "DEBUG_STRESS", &
200 236 : l_val=ec_env%debug_stress)
201 : CALL section_vals_val_get(ec_section, "DEBUG_EXTERNAL_METHOD", &
202 236 : l_val=ec_env%debug_external)
203 : ! ADMM
204 236 : CALL section_vals_val_get(ec_section, "ADMM", l_val=ec_env%do_ec_admm)
205 :
206 236 : ec_env%do_skip = .FALSE.
207 :
208 : ! set basis
209 236 : CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, nkind=nkind)
210 236 : CALL uppercase(ec_env%basis)
211 384 : SELECT CASE (ec_env%basis)
212 : CASE ("ORBITAL")
213 316 : DO ikind = 1, nkind
214 168 : qs_kind => qs_kind_set(ikind)
215 168 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
216 316 : IF (ASSOCIATED(basis_set)) THEN
217 168 : NULLIFY (harris_basis)
218 168 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
219 168 : IF (ASSOCIATED(harris_basis)) THEN
220 6 : CALL remove_basis_from_container(qs_kind%basis_sets, basis_type="HARRIS")
221 : END IF
222 168 : NULLIFY (harris_basis)
223 168 : CALL copy_gto_basis_set(basis_set, harris_basis)
224 168 : CALL add_basis_set_to_container(qs_kind%basis_sets, harris_basis, "HARRIS")
225 : END IF
226 : END DO
227 : CASE ("PRIMITIVE")
228 6 : DO ikind = 1, nkind
229 4 : qs_kind => qs_kind_set(ikind)
230 4 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
231 6 : IF (ASSOCIATED(basis_set)) THEN
232 4 : NULLIFY (harris_basis)
233 4 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
234 4 : IF (ASSOCIATED(harris_basis)) THEN
235 0 : CALL remove_basis_from_container(qs_kind%basis_sets, basis_type="HARRIS")
236 : END IF
237 4 : NULLIFY (harris_basis)
238 4 : CALL create_primitive_basis_set(basis_set, harris_basis)
239 4 : CALL get_qs_env(qs_env, dft_control=dft_control)
240 4 : eps_pgf_orb = dft_control%qs_control%eps_pgf_orb
241 4 : CALL init_interaction_radii_orb_basis(harris_basis, eps_pgf_orb)
242 4 : harris_basis%kind_radius = basis_set%kind_radius
243 4 : CALL add_basis_set_to_container(qs_kind%basis_sets, harris_basis, "HARRIS")
244 : END IF
245 : END DO
246 : CASE ("HARRIS")
247 212 : DO ikind = 1, nkind
248 126 : qs_kind => qs_kind_set(ikind)
249 126 : NULLIFY (harris_basis)
250 126 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
251 212 : IF (.NOT. ASSOCIATED(harris_basis)) THEN
252 0 : CPWARN("Harris Basis not defined for all types of atoms.")
253 : END IF
254 : END DO
255 : CASE DEFAULT
256 236 : CPABORT("Unknown basis set for energy correction (Harris functional)")
257 : END SELECT
258 : !
259 236 : CALL get_qs_kind_set(qs_kind_set, maxlgto=maxlgto, basis_type="HARRIS")
260 236 : CALL init_orbital_pointers(maxlgto + 1)
261 : !
262 236 : CALL uppercase(ec_env%basis)
263 :
264 : ! Basis may only differ from ground-state if explicitly added
265 236 : ec_env%basis_inconsistent = .FALSE.
266 236 : IF (ec_env%basis == "HARRIS") THEN
267 212 : DO ikind = 1, nkind
268 126 : qs_kind => qs_kind_set(ikind)
269 : ! Basis sets of ground-state
270 126 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=basis_set, basis_type="ORB")
271 : ! Basis sets of energy correction
272 126 : CALL get_qs_kind(qs_kind=qs_kind, basis_set=harris_basis, basis_type="HARRIS")
273 :
274 212 : IF (basis_set%name .NE. harris_basis%name) THEN
275 64 : ec_env%basis_inconsistent = .TRUE.
276 : END IF
277 : END DO
278 : END IF
279 :
280 : !Density-corrected DFT must be performed with the same basis as ground-state
281 236 : IF (ec_env%energy_functional == ec_functional_dc .AND. ec_env%basis_inconsistent) THEN
282 : CALL cp_abort(__LOCATION__, &
283 : "DC-DFT: Correction and ground state need to use the same basis. "// &
284 0 : "Checked by comparing basis set names only.")
285 : END IF
286 236 : IF (ec_env%energy_functional == ec_functional_ext .AND. ec_env%basis_inconsistent) THEN
287 : CALL cp_abort(__LOCATION__, &
288 : "Exteranl Energy: Correction and ground state need to use the same basis. "// &
289 0 : "Checked by comparing basis set names only.")
290 : END IF
291 : !
292 : ! set functional
293 382 : SELECT CASE (ec_env%energy_functional)
294 : CASE (ec_functional_harris)
295 146 : ec_env%ec_name = "Harris"
296 : CASE (ec_functional_dc)
297 86 : ec_env%ec_name = "DC-DFT"
298 : CASE (ec_functional_ext)
299 4 : ec_env%ec_name = "External Energy"
300 : CASE DEFAULT
301 236 : CPABORT("unknown energy correction")
302 : END SELECT
303 : ! select the XC section
304 236 : NULLIFY (xc_section)
305 236 : xc_section => section_vals_get_subs_vals(dft_section, "XC")
306 236 : section1 => section_vals_get_subs_vals(ec_section, "XC")
307 236 : section2 => section_vals_get_subs_vals(ec_section, "XC%XC_FUNCTIONAL")
308 236 : CALL section_vals_get(section2, explicit=explicit)
309 236 : IF (explicit) THEN
310 232 : CALL xc_functionals_expand(section2, section1)
311 232 : ec_env%xc_section => section1
312 : ELSE
313 4 : ec_env%xc_section => xc_section
314 : END IF
315 : ! Check whether energy correction requires the kinetic energy density and rebuild rho if necessary
316 236 : CALL get_qs_env(qs_env, dft_control=dft_control, rho=rho)
317 236 : xc_fun_section => section_vals_get_subs_vals(ec_env%xc_section, "XC_FUNCTIONAL")
318 : dft_control%use_kinetic_energy_density = dft_control%use_kinetic_energy_density .OR. &
319 236 : xc_uses_kinetic_energy_density(xc_fun_section, dft_control%lsd)
320 : ! Same for density gradient
321 : dft_control%drho_by_collocation = dft_control%drho_by_collocation .OR. &
322 : (xc_uses_norm_drho(xc_fun_section, dft_control%lsd) .AND. &
323 236 : (section_get_ival(xc_section, "XC_GRID%XC_DERIV") == xc_deriv_collocate))
324 : ! dispersion
325 1180 : ALLOCATE (dispersion_env)
326 : NULLIFY (xc_section)
327 236 : xc_section => ec_env%xc_section
328 236 : CALL get_qs_env(qs_env, atomic_kind_set=atomic_kind_set, para_env=para_env)
329 236 : CALL qs_dispersion_env_set(dispersion_env, xc_section)
330 236 : IF (dispersion_env%type == xc_vdw_fun_pairpot) THEN
331 0 : NULLIFY (pp_section)
332 0 : pp_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%PAIR_POTENTIAL")
333 0 : CALL qs_dispersion_pairpot_init(atomic_kind_set, qs_kind_set, dispersion_env, pp_section, para_env)
334 236 : ELSE IF (dispersion_env%type == xc_vdw_fun_nonloc) THEN
335 0 : CPABORT("nl-vdW functionals not available for EC calculations")
336 0 : NULLIFY (nl_section)
337 0 : nl_section => section_vals_get_subs_vals(xc_section, "VDW_POTENTIAL%NON_LOCAL")
338 0 : CALL qs_dispersion_nonloc_init(dispersion_env, para_env)
339 : END IF
340 236 : ec_env%dispersion_env => dispersion_env
341 :
342 : ! Check if hybrid functional are used
343 236 : ec_hfx_section => section_vals_get_subs_vals(ec_section, "XC%HF")
344 236 : CALL section_vals_get(ec_hfx_section, explicit=ec_env%do_ec_hfx)
345 :
346 : ! Initialize Harris LS solver environment
347 236 : ec_env%use_ls_solver = .FALSE.
348 : ec_env%use_ls_solver = (ec_env%ks_solver .EQ. ec_matrix_sign) &
349 : .OR. (ec_env%ks_solver .EQ. ec_matrix_trs4) &
350 236 : .OR. (ec_env%ks_solver .EQ. ec_matrix_tc2)
351 :
352 236 : IF (ec_env%use_ls_solver) THEN
353 22 : CALL ec_ls_create(qs_env, ec_env)
354 : END IF
355 :
356 : END IF
357 :
358 6686 : CALL timestop(handle)
359 :
360 6686 : END SUBROUTINE init_ec_env
361 :
362 : ! **************************************************************************************************
363 : !> \brief Initializes linear scaling environment for LS based solver of
364 : !> Harris energy functional and parses input section
365 : !> \param qs_env ...
366 : !> \param ec_env ...
367 : !> \par History
368 : !> 2020.10 created [Fabian Belleflamme]
369 : !> \author Fabian Belleflamme
370 : ! **************************************************************************************************
371 22 : SUBROUTINE ec_ls_create(qs_env, ec_env)
372 : TYPE(qs_environment_type), POINTER :: qs_env
373 : TYPE(energy_correction_type), POINTER :: ec_env
374 :
375 : CHARACTER(LEN=*), PARAMETER :: routineN = 'ec_ls_create'
376 :
377 : INTEGER :: handle
378 : REAL(KIND=dp) :: mu
379 : TYPE(dft_control_type), POINTER :: dft_control
380 : TYPE(ls_scf_env_type), POINTER :: ls_env
381 22 : TYPE(molecule_type), DIMENSION(:), POINTER :: molecule_set
382 22 : TYPE(particle_type), DIMENSION(:), POINTER :: particle_set
383 : TYPE(section_vals_type), POINTER :: ec_section, input
384 :
385 22 : CALL timeset(routineN, handle)
386 :
387 858 : ALLOCATE (ec_env%ls_env)
388 22 : ls_env => ec_env%ls_env
389 :
390 22 : NULLIFY (dft_control, input, ls_env%para_env)
391 :
392 : CALL get_qs_env(qs_env, &
393 : dft_control=dft_control, &
394 : input=input, &
395 : molecule_set=molecule_set, &
396 : particle_set=particle_set, &
397 : para_env=ls_env%para_env, &
398 22 : nelectron_spin=ls_env%nelectron_spin)
399 :
400 : ! copy some basic stuff
401 22 : ls_env%nspins = dft_control%nspins
402 22 : ls_env%natoms = SIZE(particle_set, 1)
403 22 : CALL ls_env%para_env%retain()
404 :
405 : ! initialize block to group to defined molecules
406 66 : ALLOCATE (ls_env%ls_mstruct%atom_to_molecule(ls_env%natoms))
407 22 : CALL molecule_of_atom(molecule_set, atom_to_mol=ls_env%ls_mstruct%atom_to_molecule)
408 :
409 22 : ls_env%do_transport = .FALSE.
410 22 : ls_env%do_pao = .FALSE.
411 22 : ls_env%ls_mstruct%do_pao = ls_env%do_pao
412 22 : ls_env%do_pexsi = .FALSE.
413 22 : ls_env%has_unit_metric = .FALSE.
414 :
415 22 : ec_section => section_vals_get_subs_vals(input, "DFT%ENERGY_CORRECTION")
416 22 : CALL section_vals_val_get(ec_section, "EPS_FILTER", r_val=ls_env%eps_filter)
417 22 : CALL section_vals_val_get(ec_section, "MU", r_val=mu)
418 22 : CALL section_vals_val_get(ec_section, "FIXED_MU", l_val=ls_env%fixed_mu)
419 66 : ls_env%mu_spin = mu
420 22 : CALL section_vals_val_get(ec_section, "S_PRECONDITIONER", i_val=ls_env%s_preconditioner_type)
421 22 : CALL section_vals_val_get(ec_section, "MATRIX_CLUSTER_TYPE", i_val=ls_env%ls_mstruct%cluster_type)
422 22 : CALL section_vals_val_get(ec_section, "S_INVERSION", i_val=ls_env%s_inversion_type)
423 22 : CALL section_vals_val_get(ec_section, "CHECK_S_INV", l_val=ls_env%check_s_inv)
424 22 : CALL section_vals_val_get(ec_section, "REPORT_ALL_SPARSITIES", l_val=ls_env%report_all_sparsities)
425 22 : CALL section_vals_val_get(ec_section, "SIGN_METHOD", i_val=ls_env%sign_method)
426 22 : CALL section_vals_val_get(ec_section, "SIGN_ORDER", i_val=ls_env%sign_order)
427 22 : CALL section_vals_val_get(ec_section, "DYNAMIC_THRESHOLD", l_val=ls_env%dynamic_threshold)
428 22 : CALL section_vals_val_get(ec_section, "NON_MONOTONIC", l_val=ls_env%non_monotonic)
429 22 : CALL section_vals_val_get(ec_section, "S_SQRT_METHOD", i_val=ls_env%s_sqrt_method)
430 22 : CALL section_vals_val_get(ec_section, "S_SQRT_ORDER", i_val=ls_env%s_sqrt_order)
431 22 : CALL section_vals_val_get(ec_section, "EPS_LANCZOS", r_val=ls_env%eps_lanczos)
432 22 : CALL section_vals_val_get(ec_section, "MAX_ITER_LANCZOS", i_val=ls_env%max_iter_lanczos)
433 :
434 24 : SELECT CASE (ec_env%ks_solver)
435 : CASE (ec_matrix_sign)
436 : ! S inverse required for Sign matrix algorithm,
437 : ! calculated either by Hotelling or multiplying S matrix sqrt inv
438 24 : SELECT CASE (ls_env%s_inversion_type)
439 : CASE (ls_s_inversion_sign_sqrt)
440 2 : ls_env%needs_s_inv = .TRUE.
441 2 : ls_env%use_s_sqrt = .TRUE.
442 : CASE (ls_s_inversion_hotelling)
443 0 : ls_env%needs_s_inv = .TRUE.
444 0 : ls_env%use_s_sqrt = .FALSE.
445 : CASE (ls_s_inversion_none)
446 0 : ls_env%needs_s_inv = .FALSE.
447 0 : ls_env%use_s_sqrt = .FALSE.
448 : CASE DEFAULT
449 2 : CPABORT("")
450 : END SELECT
451 : CASE (ec_matrix_trs4, ec_matrix_tc2)
452 20 : ls_env%needs_s_inv = .FALSE.
453 20 : ls_env%use_s_sqrt = .TRUE.
454 : CASE DEFAULT
455 22 : CPABORT("")
456 : END SELECT
457 :
458 22 : SELECT CASE (ls_env%s_preconditioner_type)
459 : CASE (ls_s_preconditioner_none)
460 0 : ls_env%has_s_preconditioner = .FALSE.
461 : CASE DEFAULT
462 22 : ls_env%has_s_preconditioner = .TRUE.
463 : END SELECT
464 :
465 : ! buffer for the history of matrices, not needed here
466 22 : ls_env%extrapolation_order = 0
467 22 : ls_env%scf_history%nstore = 0
468 22 : ls_env%scf_history%istore = 0
469 44 : ALLOCATE (ls_env%scf_history%matrix(ls_env%nspins, ls_env%scf_history%nstore))
470 :
471 22 : NULLIFY (ls_env%mixing_store)
472 :
473 22 : CALL timestop(handle)
474 :
475 44 : END SUBROUTINE ec_ls_create
476 :
477 : ! **************************************************************************************************
478 : !> \brief Print out the energy correction input section
479 : !>
480 : !> \param ec_env ...
481 : !> \par History
482 : !> 2020.10 created [Fabian Belleflamme]
483 : !> \author Fabian Belleflamme
484 : ! **************************************************************************************************
485 236 : SUBROUTINE ec_write_input(ec_env)
486 : TYPE(energy_correction_type), POINTER :: ec_env
487 :
488 : CHARACTER(LEN=*), PARAMETER :: routineN = 'ec_write_input'
489 :
490 : INTEGER :: handle, unit_nr
491 : TYPE(cp_logger_type), POINTER :: logger
492 : TYPE(ls_scf_env_type), POINTER :: ls_env
493 :
494 236 : CALL timeset(routineN, handle)
495 :
496 236 : logger => cp_get_default_logger()
497 236 : IF (logger%para_env%is_source()) THEN
498 118 : unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
499 : ELSE
500 : unit_nr = -1
501 : END IF
502 :
503 118 : IF (unit_nr > 0) THEN
504 :
505 : WRITE (unit_nr, '(T2,A)') &
506 118 : "!"//REPEAT("-", 29)//" Energy Correction "//REPEAT("-", 29)//"!"
507 :
508 : ! Type of energy correction
509 191 : SELECT CASE (ec_env%energy_functional)
510 : CASE (ec_functional_harris)
511 73 : WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "HARRIS FUNCTIONAL"
512 : CASE (ec_functional_dc)
513 43 : WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "DC-DFT"
514 : CASE (ec_functional_ext)
515 118 : WRITE (unit_nr, '(T2,A,T61,A20)') "Energy Correction: ", "External"
516 : END SELECT
517 118 : WRITE (unit_nr, '()')
518 :
519 : ! Energy correction parameters
520 118 : WRITE (unit_nr, '(T2,A,T61,E20.3)') "eps_default:", ec_env%eps_default
521 :
522 118 : CALL uppercase(ec_env%basis)
523 192 : SELECT CASE (ec_env%basis)
524 : CASE ("ORBITAL")
525 74 : WRITE (unit_nr, '(T2,A,T61,A20)') "EC basis: ", "ORBITAL"
526 : CASE ("PRIMITIVE")
527 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "EC basis: ", "PRIMITIVE"
528 : CASE ("HARRIS")
529 118 : WRITE (unit_nr, '(T2,A,T61,A20)') "EC Basis: ", "HARRIS"
530 : END SELECT
531 :
532 : ! Info how HFX in energy correction is treated
533 118 : IF (ec_env%do_ec_hfx) THEN
534 :
535 8 : WRITE (unit_nr, '(T2,A,T61,L20)') "DC-DFT with HFX", ec_env%do_ec_hfx
536 8 : WRITE (unit_nr, '(T2,A,T61,L20)') "Reuse HFX integrals", ec_env%reuse_hfx
537 8 : WRITE (unit_nr, '(T2,A,T61,L20)') "DC-DFT HFX with ADMM", ec_env%do_ec_admm
538 :
539 : END IF ! ec_env%do_ec_hfx
540 :
541 : ! Parameters for Harris functional solver
542 118 : IF (ec_env%energy_functional == ec_functional_harris) THEN
543 :
544 : ! Algorithm
545 133 : SELECT CASE (ec_env%ks_solver)
546 : CASE (ec_diagonalization)
547 60 : WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "DIAGONALIZATION"
548 : CASE (ec_ot_diag)
549 2 : WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "OT DIAGONALIZATION"
550 : CASE (ec_matrix_sign)
551 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "MATRIX_SIGN"
552 : CASE (ec_matrix_trs4)
553 9 : WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "TRS4"
554 9 : CALL cite_reference(Niklasson2003)
555 : CASE (ec_matrix_tc2)
556 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "Algorithm: ", "TC2"
557 74 : CALL cite_reference(Niklasson2014)
558 : END SELECT
559 73 : WRITE (unit_nr, '()')
560 :
561 : ! MAO
562 73 : IF (ec_env%mao) THEN
563 2 : WRITE (unit_nr, '(T2,A,T61,L20)') "MAO:", ec_env%mao
564 2 : WRITE (unit_nr, '(T2,A,T61,L20)') "MAO_IOLEVEL:", ec_env%mao_iolevel
565 2 : WRITE (unit_nr, '(T2,A,T61,I20)') "MAO_MAX_ITER:", ec_env%mao_max_iter
566 2 : WRITE (unit_nr, '(T2,A,T61,E20.3)') "MAO_EPS_GRAD:", ec_env%mao_eps_grad
567 2 : WRITE (unit_nr, '(T2,A,T61,E20.3)') "MAO_EPS1:", ec_env%mao_eps1
568 2 : WRITE (unit_nr, '()')
569 : END IF
570 :
571 : ! Parameters for linear response solver
572 73 : IF (.NOT. ec_env%use_ls_solver) THEN
573 :
574 62 : WRITE (unit_nr, '(T2,A)') "MO Solver"
575 62 : WRITE (unit_nr, '()')
576 :
577 122 : SELECT CASE (ec_env%ks_solver)
578 : CASE (ec_diagonalization)
579 :
580 60 : SELECT CASE (ec_env%factorization)
581 : CASE (kg_cholesky)
582 60 : WRITE (unit_nr, '(T2,A,T61,A20)') "Factorization: ", "CHOLESKY"
583 : END SELECT
584 :
585 : CASE (ec_ot_diag)
586 :
587 : ! OT Diagonalization
588 : ! Initial guess : 1) block diagonal initial guess
589 : ! 2) GS-density matrix (might require trafo if basis diff)
590 :
591 2 : SELECT CASE (ec_env%ec_initial_guess)
592 : CASE (ec_ot_atomic)
593 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "OT Diag initial guess: ", "ATOMIC"
594 : CASE (ec_ot_gs)
595 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "OT Diag initial guess: ", "GROUND STATE DM"
596 : END SELECT
597 :
598 : CASE DEFAULT
599 62 : CPABORT("Unknown Diagonalization algorithm for Harris functional")
600 : END SELECT
601 :
602 : ELSE
603 :
604 11 : WRITE (unit_nr, '(T2,A)') "AO Solver"
605 11 : WRITE (unit_nr, '()')
606 :
607 11 : ls_env => ec_env%ls_env
608 11 : WRITE (unit_nr, '(T2,A,T61,E20.3)') "eps_filter:", ls_env%eps_filter
609 11 : WRITE (unit_nr, '(T2,A,T61,L20)') "fixed chemical potential (mu)", ls_env%fixed_mu
610 11 : WRITE (unit_nr, '(T2,A,T61,L20)') "Computing inv(S):", ls_env%needs_s_inv
611 11 : WRITE (unit_nr, '(T2,A,T61,L20)') "Computing sqrt(S):", ls_env%use_s_sqrt
612 11 : WRITE (unit_nr, '(T2,A,T61,L20)') "Computing S preconditioner ", ls_env%has_s_preconditioner
613 :
614 11 : IF (ls_env%use_s_sqrt) THEN
615 21 : SELECT CASE (ls_env%s_sqrt_method)
616 : CASE (ls_s_sqrt_ns)
617 10 : WRITE (unit_nr, '(T2,A,T61,A20)') "S sqrt method:", "NEWTONSCHULZ"
618 : CASE (ls_s_sqrt_proot)
619 1 : WRITE (unit_nr, '(T2,A,T61,A20)') "S sqrt method:", "PROOT"
620 : CASE DEFAULT
621 11 : CPABORT("Unknown sqrt method.")
622 : END SELECT
623 11 : WRITE (unit_nr, '(T2,A,T61,I20)') "S sqrt order:", ls_env%s_sqrt_order
624 : END IF
625 :
626 11 : SELECT CASE (ls_env%s_preconditioner_type)
627 : CASE (ls_s_preconditioner_none)
628 0 : WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "NONE"
629 : CASE (ls_s_preconditioner_atomic)
630 11 : WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "ATOMIC"
631 : CASE (ls_s_preconditioner_molecular)
632 11 : WRITE (unit_nr, '(T2,A,T61,A20)') "S preconditioner type ", "MOLECULAR"
633 : END SELECT
634 :
635 22 : SELECT CASE (ls_env%ls_mstruct%cluster_type)
636 : CASE (ls_cluster_atomic)
637 11 : WRITE (unit_nr, '(T2,A,T61,A20)') "Cluster type", ADJUSTR("ATOMIC")
638 : CASE (ls_cluster_molecular)
639 0 : WRITE (unit_nr, '(T2,A,T61,A20)') "Cluster type", ADJUSTR("MOLECULAR")
640 : CASE DEFAULT
641 11 : CPABORT("Unknown cluster type")
642 : END SELECT
643 :
644 : END IF
645 :
646 : END IF ! if ec_functional_harris
647 :
648 118 : WRITE (unit_nr, '(T2,A)') REPEAT("-", 79)
649 118 : WRITE (unit_nr, '()')
650 :
651 : END IF ! unit_nr
652 :
653 236 : CALL timestop(handle)
654 :
655 236 : END SUBROUTINE ec_write_input
656 :
657 : END MODULE ec_environment
|