Line data Source code
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 function that build the dft section of the input
10 : !> \par History
11 : !> 10.2005 moved out of input_cp2k [fawzi]
12 : !> \author fawzi
13 : ! **************************************************************************************************
14 : MODULE input_cp2k_dft
15 : USE basis_set_types, ONLY: basis_sort_default,&
16 : basis_sort_zet
17 : USE bibliography, ONLY: &
18 : Andermatt2016, Andreussi2012, Avezac2005, Bengtsson1999, Blochl1995, Brelaz1979, &
19 : Chai2024a, Fattebert2002, Guidon2010, Iannuzzi2006, Kunert2003, Marek2025, Merlot2014, &
20 : Perdew1981, VandeVondele2005b, Yin2017
21 : USE cp_output_handling, ONLY: add_last_numeric,&
22 : cp_print_key_section_create,&
23 : high_print_level,&
24 : low_print_level,&
25 : medium_print_level,&
26 : silent_print_level
27 : USE cp_spline_utils, ONLY: pw_interp,&
28 : spline3_nopbc_interp,&
29 : spline3_pbc_interp
30 : USE cp_units, ONLY: cp_unit_to_cp2k
31 : USE input_constants, ONLY: &
32 : admm1_type, admm2_type, admmp_type, admmq_type, admms_type, do_admm_aux_exch_func_bee, &
33 : do_admm_aux_exch_func_bee_libxc, do_admm_aux_exch_func_default, &
34 : do_admm_aux_exch_func_default_libxc, do_admm_aux_exch_func_none, &
35 : do_admm_aux_exch_func_opt, do_admm_aux_exch_func_opt_libxc, do_admm_aux_exch_func_pbex, &
36 : do_admm_aux_exch_func_pbex_libxc, do_admm_aux_exch_func_sx_libxc, &
37 : do_admm_basis_projection, do_admm_blocked_projection, do_admm_blocking_purify_full, &
38 : do_admm_charge_constrained_projection, do_admm_exch_scaling_merlot, &
39 : do_admm_exch_scaling_none, do_admm_purify_cauchy, do_admm_purify_cauchy_subspace, &
40 : do_admm_purify_mcweeny, do_admm_purify_mo_diag, do_admm_purify_mo_no_diag, &
41 : do_admm_purify_none, do_admm_purify_none_dm, do_arnoldi, do_bch, do_cn, do_em, do_etrs, &
42 : do_exact, do_pade, do_taylor, ehrenfest, gaussian, kg_color_dsatur, kg_color_greedy, &
43 : kg_tnadd_atomic, kg_tnadd_embed, kg_tnadd_embed_ri, kg_tnadd_none, no_admm_type, &
44 : numerical, plus_u_lowdin, plus_u_mulliken, plus_u_mulliken_charges, real_time_propagation, &
45 : rel_dkh, rel_none, rel_pot_erfc, rel_pot_full, rel_sczora_mp, rel_trans_atom, &
46 : rel_trans_full, rel_trans_molecule, rel_zora, rel_zora_full, rel_zora_mp, &
47 : rtp_bse_ham_g0w0, rtp_bse_ham_ks, rtp_method_bse, rtp_method_tddft, sccs_andreussi, &
48 : sccs_derivative_cd3, sccs_derivative_cd5, sccs_derivative_cd7, sccs_derivative_fft, &
49 : sccs_fattebert_gygi, sic_ad, sic_eo, sic_list_all, sic_list_unpaired, sic_mauri_spz, &
50 : sic_mauri_us, sic_none, slater, use_mom_ref_coac, use_mom_ref_com, use_mom_ref_user, &
51 : use_mom_ref_zero, use_restart_wfn, use_rt_restart, use_scf_wfn, weight_type_mass, &
52 : weight_type_unit
53 : USE input_cp2k_almo, ONLY: create_almo_scf_section
54 : USE input_cp2k_as, ONLY: create_active_space_section
55 : USE input_cp2k_ec, ONLY: create_ec_section
56 : USE input_cp2k_exstate, ONLY: create_exstate_section
57 : USE input_cp2k_external, ONLY: create_ext_den_section,&
58 : create_ext_pot_section,&
59 : create_ext_vxc_section
60 : USE input_cp2k_field, ONLY: create_efield_section,&
61 : create_per_efield_section
62 : USE input_cp2k_harris, ONLY: create_harris_section
63 : USE input_cp2k_kpoints, ONLY: create_kpoint_set_section,&
64 : create_kpoints_section
65 : USE input_cp2k_loc, ONLY: create_localize_section
66 : USE input_cp2k_ls, ONLY: create_ls_scf_section
67 : USE input_cp2k_poisson, ONLY: create_poisson_section
68 : USE input_cp2k_print_dft, ONLY: create_print_dft_section
69 : USE input_cp2k_projection_rtp, ONLY: create_projection_rtp_section
70 : USE input_cp2k_qs, ONLY: create_lrigpw_section,&
71 : create_qs_section
72 : USE input_cp2k_rsgrid, ONLY: create_rsgrid_section
73 : USE input_cp2k_scf, ONLY: create_scf_section
74 : USE input_cp2k_smeagol, ONLY: create_dft_smeagol_section
75 : USE input_cp2k_transport, ONLY: create_transport_section
76 : USE input_cp2k_xas, ONLY: create_xas_section,&
77 : create_xas_tdp_section
78 : USE input_cp2k_xc, ONLY: create_xc_section
79 : USE input_keyword_types, ONLY: keyword_create,&
80 : keyword_release,&
81 : keyword_type
82 : USE input_section_types, ONLY: section_add_keyword,&
83 : section_add_subsection,&
84 : section_create,&
85 : section_release,&
86 : section_type
87 : USE input_val_types, ONLY: char_t,&
88 : integer_t,&
89 : lchar_t,&
90 : logical_t,&
91 : real_t
92 : USE kinds, ONLY: dp
93 : USE physcon, ONLY: evolt,&
94 : femtoseconds
95 : USE pw_spline_utils, ONLY: no_precond,&
96 : precond_spl3_1,&
97 : precond_spl3_2,&
98 : precond_spl3_3,&
99 : precond_spl3_aint,&
100 : precond_spl3_aint2
101 : USE string_utilities, ONLY: s2a
102 : #include "./base/base_uses.f90"
103 :
104 : IMPLICIT NONE
105 : PRIVATE
106 :
107 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_dft'
108 :
109 : PUBLIC :: create_dft_section
110 : PUBLIC :: create_bsse_section
111 : PUBLIC :: create_interp_section
112 : PUBLIC :: create_mgrid_section
113 :
114 : CONTAINS
115 :
116 : ! **************************************************************************************************
117 : !> \brief creates the dft section
118 : !> \param section the section to be created
119 : !> \author fawzi
120 : ! **************************************************************************************************
121 10278 : SUBROUTINE create_dft_section(section)
122 : TYPE(section_type), POINTER :: section
123 :
124 : TYPE(keyword_type), POINTER :: keyword
125 : TYPE(section_type), POINTER :: subsection
126 :
127 10278 : CPASSERT(.NOT. ASSOCIATED(section))
128 : CALL section_create(section, __LOCATION__, name="DFT", &
129 : description="Controls electronic-structure settings for Quickstep and related "// &
130 : "Gaussian-basis DFT methods.", &
131 10278 : n_keywords=3, n_subsections=4, repeats=.FALSE.)
132 :
133 10278 : NULLIFY (keyword)
134 : CALL keyword_create(keyword, __LOCATION__, name="BASIS_SET_FILE_NAME", &
135 : description="Name of a basis-set library file, optionally including a path. "// &
136 : "This keyword can be repeated to search several basis-set files.", &
137 : usage="BASIS_SET_FILE_NAME <FILENAME>", &
138 : type_of_var=lchar_t, repeats=.TRUE., &
139 10278 : default_lc_val="BASIS_SET", n_var=1)
140 10278 : CALL section_add_keyword(section, keyword)
141 10278 : CALL keyword_release(keyword)
142 :
143 : CALL keyword_create(keyword, __LOCATION__, name="POTENTIAL_FILE_NAME", &
144 : description="Name of the pseudopotential library file, optionally including a path. "// &
145 : "The potential selected for each kind is set with KIND%POTENTIAL.", &
146 : usage="POTENTIAL_FILE_NAME <FILENAME>", &
147 10278 : default_lc_val="POTENTIAL")
148 10278 : CALL section_add_keyword(section, keyword)
149 10278 : CALL keyword_release(keyword)
150 :
151 : CALL keyword_create(keyword, __LOCATION__, name="WFN_RESTART_FILE_NAME", &
152 : variants=["RESTART_FILE_NAME"], &
153 : description="Name of the wavefunction restart file, may include a path."// &
154 : " If no file is specified, the default is to open the file as generated by the wfn restart print key.", &
155 : usage="WFN_RESTART_FILE_NAME <FILENAME>", &
156 20556 : type_of_var=lchar_t)
157 10278 : CALL section_add_keyword(section, keyword)
158 10278 : CALL keyword_release(keyword)
159 :
160 : CALL keyword_create(keyword, __LOCATION__, &
161 : name="UKS", &
162 : variants=s2a("UNRESTRICTED_KOHN_SHAM", &
163 : "LSD", &
164 : "SPIN_POLARIZED"), &
165 : description="Requests a spin-polarized calculation using alpha "// &
166 : "and beta orbitals, i.e. no spin restriction is applied", &
167 : usage="LSD", &
168 : default_l_val=.FALSE., &
169 10278 : lone_keyword_l_val=.TRUE.)
170 10278 : CALL section_add_keyword(section, keyword)
171 10278 : CALL keyword_release(keyword)
172 : CALL keyword_create(keyword, __LOCATION__, &
173 : name="ROKS", &
174 : variants=["RESTRICTED_OPEN_KOHN_SHAM"], &
175 : description="Requests a restricted open Kohn-Sham calculation", &
176 : usage="ROKS", &
177 : default_l_val=.FALSE., &
178 20556 : lone_keyword_l_val=.TRUE.)
179 10278 : CALL section_add_keyword(section, keyword)
180 10278 : CALL keyword_release(keyword)
181 : CALL keyword_create(keyword, __LOCATION__, &
182 : name="MULTIPLICITY", &
183 : variants=["MULTIP"], &
184 : description="Two times the total spin plus one. "// &
185 : "Specify 3 for a triplet, 4 for a quartet, "// &
186 : "and so on. Default is 1 (singlet) for an "// &
187 : "even number and 2 (doublet) for an odd number "// &
188 : "of electrons.", &
189 : usage="MULTIPLICITY 3", &
190 20556 : default_i_val=0) ! this default value is just a flag to get the above
191 10278 : CALL section_add_keyword(section, keyword)
192 10278 : CALL keyword_release(keyword)
193 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE", &
194 : description="The total charge of the system", &
195 : usage="CHARGE -1", &
196 10278 : default_i_val=0)
197 10278 : CALL section_add_keyword(section, keyword)
198 10278 : CALL keyword_release(keyword)
199 :
200 : CALL keyword_create(keyword, __LOCATION__, &
201 : name="PLUS_U_METHOD", &
202 : description="Method employed for the calculation of the DFT+U contribution", &
203 : repeats=.FALSE., &
204 : enum_c_vals=s2a("LOWDIN", "MULLIKEN", "MULLIKEN_CHARGES"), &
205 : enum_i_vals=[plus_u_lowdin, plus_u_mulliken, plus_u_mulliken_charges], &
206 : enum_desc=s2a("Method based on Lowdin population analysis "// &
207 : "(computationally expensive, since the diagonalization of the "// &
208 : "overlap matrix is required, but possibly more robust than Mulliken)", &
209 : "Method based on Mulliken population analysis using the net AO and "// &
210 : "overlap populations (computationally cheap method)", &
211 : "Method based on Mulliken gross orbital populations (GOP)"), &
212 : n_var=1, &
213 : default_i_val=plus_u_mulliken, &
214 10278 : usage="PLUS_U_METHOD Lowdin")
215 10278 : CALL section_add_keyword(section, keyword)
216 10278 : CALL keyword_release(keyword)
217 :
218 : CALL keyword_create(keyword, __LOCATION__, &
219 : name="RELAX_MULTIPLICITY", &
220 : variants=["RELAX_MULTIP"], &
221 : description="Tolerance in Hartrees. Do not enforce the occupation "// &
222 : "of alpha and beta MOs due to the initially "// &
223 : "defined multiplicity, but rather follow the Aufbau principle. "// &
224 : "A value greater than zero activates this option. "// &
225 : "If alpha/beta MOs differ in energy less than this tolerance, "// &
226 : "then alpha-MO occupation is preferred even if it is higher "// &
227 : "in energy (within the tolerance). "// &
228 : "Such spin-symmetry broken (spin-polarized) occupation is used "// &
229 : "as SCF input, which (is assumed to) bias the SCF "// &
230 : "towards a spin-polarized solution. "// &
231 : "Thus, larger tolerance increases chances of ending up "// &
232 : "with spin-polarization. "// &
233 : "This option is only valid for unrestricted (i.e. spin polarised) "// &
234 : "Kohn-Sham (UKS) calculations. It also needs non-zero "// &
235 : "[ADDED_MOS](#CP2K_INPUT.FORCE_EVAL.DFT.SCF.ADDED_MOS) to actually affect the calculations, "// &
236 : "which is why it is not expected to work with [OT](#CP2K_INPUT.FORCE_EVAL.DFT.SCF.OT) "// &
237 : "and may raise errors when used with OT. "// &
238 : "For more details see [this discussion](https://github.com/cp2k/cp2k/issues/4389).", &
239 : usage="RELAX_MULTIPLICITY 0.00001", &
240 : repeats=.FALSE., &
241 20556 : default_r_val=0.0_dp)
242 10278 : CALL section_add_keyword(section, keyword)
243 10278 : CALL keyword_release(keyword)
244 :
245 : CALL keyword_create(keyword, __LOCATION__, name="SUBCELLS", &
246 : description="Read the grid size for subcell generation in the construction of "// &
247 : "neighbor lists.", usage="SUBCELLS 1.5", &
248 10278 : n_var=1, default_r_val=2.0_dp)
249 10278 : CALL section_add_keyword(section, keyword)
250 10278 : CALL keyword_release(keyword)
251 :
252 : CALL keyword_create(keyword, __LOCATION__, name="AUTO_BASIS", &
253 : description="Specify size of automatically generated auxiliary (RI) basis sets: "// &
254 : "Options={small,medium,large,huge}", &
255 : usage="AUTO_BASIS {basis_type} {basis_size}", &
256 30834 : type_of_var=char_t, repeats=.TRUE., n_var=-1, default_c_vals=["X", "X"])
257 10278 : CALL section_add_keyword(section, keyword)
258 10278 : CALL keyword_release(keyword)
259 :
260 : CALL keyword_create(keyword, __LOCATION__, &
261 : name="SURFACE_DIPOLE_CORRECTION", &
262 : variants=s2a("SURFACE_DIPOLE", &
263 : "SURF_DIP"), &
264 : description="For slab calculations with asymmetric geometries, activate the correction of "// &
265 : "the electrostatic potential with "// &
266 : "by compensating for the surface dipole. Implemented only for slabs with normal "// &
267 : "parallel to one Cartesian axis. The normal direction is given by the keyword SURF_DIP_DIR", &
268 : usage="SURF_DIP", &
269 : default_l_val=.FALSE., &
270 : lone_keyword_l_val=.TRUE., &
271 20556 : citations=[Bengtsson1999])
272 10278 : CALL section_add_keyword(section, keyword)
273 10278 : CALL keyword_release(keyword)
274 :
275 : CALL keyword_create(keyword, __LOCATION__, &
276 : name="SURF_DIP_DIR", &
277 : description="Cartesian axis parallel to surface normal.", &
278 : enum_c_vals=s2a("X", "Y", "Z"), &
279 : enum_i_vals=[1, 2, 3], &
280 : enum_desc=s2a("Along x", "Along y", "Along z"), &
281 : n_var=1, &
282 : default_i_val=3, &
283 10278 : usage="SURF_DIP_DIR Z")
284 10278 : CALL section_add_keyword(section, keyword)
285 10278 : CALL keyword_release(keyword)
286 :
287 : CALL keyword_create(keyword, __LOCATION__, &
288 : name="SURF_DIP_POS", &
289 : description="This keyword assigns an user defined position in Angstroms "// &
290 : "in the direction normal to the surface (given by SURF_DIP_DIR). "// &
291 : "The default value is -1.0_dp which appplies the correction at a position "// &
292 : "that has minimum electron density on the grid.", &
293 : usage="SURF_DIP_POS -1.0_dp", &
294 10278 : default_r_val=-1.0_dp)
295 10278 : CALL section_add_keyword(section, keyword)
296 10278 : CALL keyword_release(keyword)
297 :
298 : CALL keyword_create(keyword, __LOCATION__, &
299 : name="SURF_DIP_SWITCH", &
300 : description="WARNING: Experimental feature under development that will help the "// &
301 : "user to switch parameters to facilitate SCF convergence. In its current form the "// &
302 : "surface dipole correction is switched off if the calculation does not converge in "// &
303 : "(0.5*MAX_SCF + 1) outer_scf steps. "// &
304 : "The default value is .FALSE.", &
305 : usage="SURF_DIP_SWITCH .TRUE.", &
306 : default_l_val=.FALSE., &
307 10278 : lone_keyword_l_val=.TRUE.)
308 10278 : CALL section_add_keyword(section, keyword)
309 10278 : CALL keyword_release(keyword)
310 :
311 : CALL keyword_create(keyword, __LOCATION__, &
312 : name="CORE_CORR_DIP", &
313 : description="If the total CORE_CORRECTION is non-zero and surface dipole "// &
314 : "correction is switched on, presence of this keyword will adjust electron "// &
315 : "density via MO occupation to reflect the total CORE_CORRECTION. "// &
316 : "The default value is .FALSE.", &
317 : usage="CORE_CORR_DIP .TRUE.", &
318 : default_l_val=.FALSE., &
319 10278 : lone_keyword_l_val=.TRUE.)
320 10278 : CALL section_add_keyword(section, keyword)
321 10278 : CALL keyword_release(keyword)
322 :
323 : CALL keyword_create(keyword, __LOCATION__, &
324 : name="SORT_BASIS", &
325 : description="Sorts basis functions according to a selected criterion. "// &
326 : "Sorting by exponent can improve data locality for selected exact-exchange and RI workflows.", &
327 : enum_c_vals=s2a("DEFAULT", "EXP"), &
328 : enum_i_vals=[basis_sort_default, basis_sort_zet], &
329 : enum_desc=s2a("don't sort", "sort w.r.t. exponent"), &
330 : default_i_val=basis_sort_default, &
331 10278 : usage="SORT_BASIS EXP")
332 10278 : CALL section_add_keyword(section, keyword)
333 10278 : CALL keyword_release(keyword)
334 :
335 10278 : NULLIFY (subsection)
336 10278 : CALL create_scf_section(subsection)
337 10278 : CALL section_add_subsection(section, subsection)
338 10278 : CALL section_release(subsection)
339 :
340 10278 : CALL create_ls_scf_section(subsection)
341 10278 : CALL section_add_subsection(section, subsection)
342 10278 : CALL section_release(subsection)
343 :
344 10278 : CALL create_almo_scf_section(subsection)
345 10278 : CALL section_add_subsection(section, subsection)
346 10278 : CALL section_release(subsection)
347 :
348 10278 : CALL create_kg_section(subsection)
349 10278 : CALL section_add_subsection(section, subsection)
350 10278 : CALL section_release(subsection)
351 :
352 10278 : CALL create_harris_section(subsection)
353 10278 : CALL section_add_subsection(section, subsection)
354 10278 : CALL section_release(subsection)
355 :
356 10278 : CALL create_ec_section(subsection)
357 10278 : CALL section_add_subsection(section, subsection)
358 10278 : CALL section_release(subsection)
359 :
360 10278 : CALL create_exstate_section(subsection)
361 10278 : CALL section_add_subsection(section, subsection)
362 10278 : CALL section_release(subsection)
363 :
364 10278 : CALL create_admm_section(subsection)
365 10278 : CALL section_add_subsection(section, subsection)
366 10278 : CALL section_release(subsection)
367 :
368 10278 : CALL create_qs_section(subsection)
369 10278 : CALL section_add_subsection(section, subsection)
370 10278 : CALL section_release(subsection)
371 :
372 10278 : CALL create_mgrid_section(subsection, create_subsections=.TRUE.)
373 10278 : CALL section_add_subsection(section, subsection)
374 10278 : CALL section_release(subsection)
375 :
376 10278 : CALL create_xc_section(subsection)
377 10278 : CALL section_add_subsection(section, subsection)
378 10278 : CALL section_release(subsection)
379 :
380 10278 : CALL create_relativistic_section(subsection)
381 10278 : CALL section_add_subsection(section, subsection)
382 10278 : CALL section_release(subsection)
383 :
384 10278 : CALL create_sic_section(subsection)
385 10278 : CALL section_add_subsection(section, subsection)
386 10278 : CALL section_release(subsection)
387 :
388 10278 : CALL create_low_spin_roks_section(subsection)
389 10278 : CALL section_add_subsection(section, subsection)
390 10278 : CALL section_release(subsection)
391 :
392 10278 : CALL create_efield_section(subsection)
393 10278 : CALL section_add_subsection(section, subsection)
394 10278 : CALL section_release(subsection)
395 :
396 10278 : CALL create_per_efield_section(subsection)
397 10278 : CALL section_add_subsection(section, subsection)
398 10278 : CALL section_release(subsection)
399 :
400 10278 : CALL create_ext_pot_section(subsection)
401 10278 : CALL section_add_subsection(section, subsection)
402 10278 : CALL section_release(subsection)
403 :
404 10278 : CALL create_transport_section(subsection)
405 10278 : CALL section_add_subsection(section, subsection)
406 10278 : CALL section_release(subsection)
407 :
408 : ! ZMP sections to include the external density or v_xc potential
409 10278 : CALL create_ext_den_section(subsection)
410 10278 : CALL section_add_subsection(section, subsection)
411 10278 : CALL section_release(subsection)
412 :
413 10278 : CALL create_ext_vxc_section(subsection)
414 10278 : CALL section_add_subsection(section, subsection)
415 10278 : CALL section_release(subsection)
416 :
417 10278 : CALL create_poisson_section(subsection)
418 10278 : CALL section_add_subsection(section, subsection)
419 10278 : CALL section_release(subsection)
420 :
421 10278 : CALL create_kpoints_section(subsection)
422 10278 : CALL section_add_subsection(section, subsection)
423 10278 : CALL section_release(subsection)
424 :
425 10278 : CALL create_kpoint_set_section(subsection)
426 10278 : CALL section_add_subsection(section, subsection)
427 10278 : CALL section_release(subsection)
428 :
429 10278 : CALL create_implicit_solv_section(subsection)
430 10278 : CALL section_add_subsection(section, subsection)
431 10278 : CALL section_release(subsection)
432 :
433 10278 : CALL create_density_fitting_section(subsection)
434 10278 : CALL section_add_subsection(section, subsection)
435 10278 : CALL section_release(subsection)
436 :
437 10278 : CALL create_xas_section(subsection)
438 10278 : CALL section_add_subsection(section, subsection)
439 10278 : CALL section_release(subsection)
440 :
441 10278 : CALL create_xas_tdp_section(subsection)
442 10278 : CALL section_add_subsection(section, subsection)
443 10278 : CALL section_release(subsection)
444 :
445 10278 : CALL create_localize_section(subsection)
446 10278 : CALL section_add_subsection(section, subsection)
447 10278 : CALL section_release(subsection)
448 :
449 10278 : CALL create_rtp_section(subsection)
450 10278 : CALL section_add_subsection(section, subsection)
451 10278 : CALL section_release(subsection)
452 :
453 10278 : CALL create_print_dft_section(subsection)
454 10278 : CALL section_add_subsection(section, subsection)
455 10278 : CALL section_release(subsection)
456 :
457 10278 : CALL create_sccs_section(subsection)
458 10278 : CALL section_add_subsection(section, subsection)
459 10278 : CALL section_release(subsection)
460 :
461 10278 : CALL create_active_space_section(subsection)
462 10278 : CALL section_add_subsection(section, subsection)
463 10278 : CALL section_release(subsection)
464 :
465 10278 : CALL create_dft_smeagol_section(subsection)
466 10278 : CALL section_add_subsection(section, subsection)
467 10278 : CALL section_release(subsection)
468 :
469 10278 : CALL create_hairy_probes_section(subsection)
470 10278 : CALL section_add_subsection(section, subsection)
471 10278 : CALL section_release(subsection)
472 :
473 10278 : CALL create_pcc_section(subsection)
474 10278 : CALL section_add_subsection(section, subsection)
475 10278 : CALL section_release(subsection)
476 :
477 10278 : CALL create_paep_section(subsection)
478 10278 : CALL section_add_subsection(section, subsection)
479 10278 : CALL section_release(subsection)
480 :
481 10278 : END SUBROUTINE create_dft_section
482 :
483 : ! **************************************************************************************************
484 : !> \brief Hairy Probe DFT Model
485 : !> \param section ...
486 : !> \author Margherita Buraschi
487 : ! **************************************************************************************************
488 :
489 10278 : SUBROUTINE create_hairy_probes_section(section)
490 : TYPE(section_type), POINTER :: section
491 :
492 : TYPE(keyword_type), POINTER :: keyword
493 :
494 10278 : NULLIFY (keyword)
495 10278 : CPASSERT(.NOT. ASSOCIATED(section))
496 : CALL section_create(section, __LOCATION__, &
497 : name="HAIRY_PROBES", &
498 : description="Sets up a Hairy Probe calculation. ", &
499 10278 : n_keywords=0, n_subsections=0, repeats=.TRUE.)
500 :
501 : CALL keyword_create(keyword, __LOCATION__, &
502 : name="_SECTION_PARAMETERS_", &
503 : description="Controls the activation of hairy probe", &
504 : usage="&HAIRY_PROBES ON", &
505 : default_l_val=.FALSE., &
506 10278 : lone_keyword_l_val=.TRUE.)
507 10278 : CALL section_add_keyword(section, keyword)
508 10278 : CALL keyword_release(keyword)
509 :
510 : CALL keyword_create(keyword, __LOCATION__, name="ATOM_IDS", &
511 : description="Indexes of the atoms to which the probes are attached.", &
512 : usage="ATOM_IDS <INTEGER> .. <INTEGER>", &
513 10278 : type_of_var=integer_t, n_var=-1)
514 10278 : CALL section_add_keyword(section, keyword)
515 10278 : CALL keyword_release(keyword)
516 :
517 : CALL keyword_create(keyword, __LOCATION__, name="T", &
518 : description="Electronic temperature [K]", &
519 : usage="T <REAL>", &
520 : default_r_val=cp_unit_to_cp2k(value=300.0_dp, unit_str="K"), &
521 10278 : unit_str="K")
522 10278 : CALL section_add_keyword(section, keyword)
523 10278 : CALL keyword_release(keyword)
524 :
525 : CALL keyword_create(keyword, __LOCATION__, name="MU", &
526 : description="Chemical potential of the electrons in the probes [eV] ", &
527 : usage="MU <REAL>", &
528 : default_r_val=cp_unit_to_cp2k(value=0.0_dp, unit_str="eV"), &
529 10278 : unit_str="eV")
530 10278 : CALL section_add_keyword(section, keyword)
531 10278 : CALL keyword_release(keyword)
532 :
533 : CALL keyword_create(keyword, __LOCATION__, name="ALPHA", &
534 : description="Parameter for solution probes ", &
535 : usage="ALPHA <REAL>", &
536 10278 : default_r_val=1.0_dp)
537 10278 : CALL section_add_keyword(section, keyword)
538 10278 : CALL keyword_release(keyword)
539 :
540 : CALL keyword_create(keyword, __LOCATION__, name="EPS_HP", &
541 : description=" Tolerance for accuracy checks on occupation numbers "// &
542 : "calculated using hair-probes. ", &
543 : usage="EPS_HP <REAL>", &
544 10278 : default_r_val=1.0E-5_dp)
545 10278 : CALL section_add_keyword(section, keyword)
546 10278 : CALL keyword_release(keyword)
547 10278 : END SUBROUTINE create_hairy_probes_section
548 : !####################################################################################
549 :
550 : ! **************************************************************************************************
551 : !> \brief Implicit Solvation Model
552 : !> \param section ...
553 : !> \author tlaino
554 : ! **************************************************************************************************
555 10278 : SUBROUTINE create_implicit_solv_section(section)
556 : TYPE(section_type), POINTER :: section
557 :
558 : TYPE(keyword_type), POINTER :: keyword
559 : TYPE(section_type), POINTER :: print_key, subsection
560 :
561 10278 : NULLIFY (keyword, subsection, print_key)
562 10278 : CPASSERT(.NOT. ASSOCIATED(section))
563 : CALL section_create(section, __LOCATION__, name="SCRF", &
564 : description="Adds an implicit solvation model to the DFT calculation."// &
565 : " Know also as Self Consistent Reaction Field.", &
566 10278 : n_keywords=0, n_subsections=0, repeats=.FALSE.)
567 :
568 : CALL keyword_create(keyword, __LOCATION__, name="EPS_OUT", &
569 : description="Value of the dielectric constant outside the sphere", &
570 : usage="EPS_OUT <REAL>", &
571 10278 : default_r_val=1.0_dp)
572 10278 : CALL section_add_keyword(section, keyword)
573 10278 : CALL keyword_release(keyword)
574 :
575 : CALL keyword_create(keyword, __LOCATION__, name="LMAX", &
576 : description="Maximum value of L used in the multipole expansion", &
577 : usage="LMAX <INTEGER>", &
578 10278 : default_i_val=3)
579 10278 : CALL section_add_keyword(section, keyword)
580 10278 : CALL keyword_release(keyword)
581 :
582 10278 : CALL create_sphere_section(subsection)
583 10278 : CALL section_add_subsection(section, subsection)
584 10278 : CALL section_release(subsection)
585 :
586 : CALL cp_print_key_section_create(print_key, __LOCATION__, "program_run_info", &
587 : description="Controls the printing basic info about the method", &
588 10278 : print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
589 10278 : CALL section_add_subsection(section, print_key)
590 10278 : CALL section_release(print_key)
591 :
592 10278 : END SUBROUTINE create_implicit_solv_section
593 :
594 : ! **************************************************************************************************
595 : !> \brief Create Sphere cavity
596 : !> \param section ...
597 : !> \author tlaino
598 : ! **************************************************************************************************
599 10278 : SUBROUTINE create_sphere_section(section)
600 : TYPE(section_type), POINTER :: section
601 :
602 : TYPE(keyword_type), POINTER :: keyword
603 : TYPE(section_type), POINTER :: subsection
604 :
605 10278 : NULLIFY (keyword, subsection)
606 10278 : CPASSERT(.NOT. ASSOCIATED(section))
607 : CALL section_create(section, __LOCATION__, name="SPHERE", &
608 : description="Treats the implicit solvent environment like a sphere", &
609 10278 : n_keywords=0, n_subsections=0, repeats=.FALSE.)
610 :
611 : CALL keyword_create(keyword, __LOCATION__, name="RADIUS", &
612 : description="Value of the spherical cavity in the dielectric medium", &
613 : usage="RADIUS <REAL>", &
614 : unit_str="angstrom", &
615 10278 : type_of_var=real_t)
616 10278 : CALL section_add_keyword(section, keyword)
617 10278 : CALL keyword_release(keyword)
618 :
619 10278 : CALL create_center_section(subsection)
620 10278 : CALL section_add_subsection(section, subsection)
621 10278 : CALL section_release(subsection)
622 :
623 10278 : END SUBROUTINE create_sphere_section
624 :
625 : ! **************************************************************************************************
626 : !> \brief ...
627 : !> \param section ...
628 : !> \author tlaino
629 : ! **************************************************************************************************
630 10278 : SUBROUTINE create_center_section(section)
631 : TYPE(section_type), POINTER :: section
632 :
633 : TYPE(keyword_type), POINTER :: keyword
634 :
635 10278 : NULLIFY (keyword)
636 10278 : CPASSERT(.NOT. ASSOCIATED(section))
637 : CALL section_create(section, __LOCATION__, name="CENTER", &
638 : description="Defines the center of the sphere.", &
639 10278 : n_keywords=0, n_subsections=0, repeats=.FALSE.)
640 : CALL keyword_create(keyword, __LOCATION__, name="XYZ", &
641 : description="Coordinates of the center of the sphere", &
642 : usage="XYZ <REAL> <REAL> <REAL>", &
643 : unit_str="angstrom", &
644 10278 : type_of_var=real_t, n_var=3)
645 10278 : CALL section_add_keyword(section, keyword)
646 10278 : CALL keyword_release(keyword)
647 :
648 : CALL keyword_create(keyword, __LOCATION__, name="ATOM_LIST", &
649 : description="Defines a list of atoms to define the center of the sphere", &
650 : usage="ATOM_LIST <INTEGER> .. <INTEGER>", &
651 10278 : type_of_var=integer_t, n_var=-1)
652 10278 : CALL section_add_keyword(section, keyword)
653 10278 : CALL keyword_release(keyword)
654 :
655 : CALL keyword_create(keyword, __LOCATION__, name="WEIGHT_TYPE", &
656 : description="Defines the weight used to define the center of the sphere"// &
657 : " (if ATOM_LIST is provided)", &
658 : usage="WEIGHT_TYPE (UNIT|MASS)", &
659 : enum_c_vals=["UNIT", "MASS"], &
660 : enum_i_vals=[weight_type_unit, weight_type_mass], &
661 30834 : default_i_val=weight_type_unit)
662 10278 : CALL section_add_keyword(section, keyword)
663 10278 : CALL keyword_release(keyword)
664 :
665 : CALL keyword_create(keyword, __LOCATION__, name="FIXED", &
666 : description="Specify if the center of the sphere should be fixed or"// &
667 : " allowed to move", &
668 : usage="FIXED <LOGICAL>", &
669 10278 : default_l_val=.TRUE.)
670 10278 : CALL section_add_keyword(section, keyword)
671 10278 : CALL keyword_release(keyword)
672 :
673 10278 : END SUBROUTINE create_center_section
674 :
675 : ! **************************************************************************************************
676 : !> \brief ...
677 : !> \param section ...
678 : ! **************************************************************************************************
679 10278 : SUBROUTINE create_admm_section(section)
680 : TYPE(section_type), POINTER :: section
681 :
682 : TYPE(keyword_type), POINTER :: keyword
683 :
684 10278 : NULLIFY (keyword)
685 10278 : CPASSERT(.NOT. ASSOCIATED(section))
686 : CALL section_create(section, __LOCATION__, name="AUXILIARY_DENSITY_MATRIX_METHOD", &
687 : description="Controls the auxiliary density matrix method (ADMM), which evaluates "// &
688 : "Hartree-Fock exchange on a smaller auxiliary basis and adds an exchange correction.", &
689 : n_keywords=1, n_subsections=1, repeats=.FALSE., &
690 20556 : citations=[Guidon2010])
691 :
692 : CALL keyword_create( &
693 : keyword, __LOCATION__, &
694 : name="ADMM_TYPE", &
695 : description="Named ADMM variant from the literature. This shortcut sets METHOD, "// &
696 : "ADMM_PURIFICATION_METHOD, and EXCH_SCALING_MODEL consistently for the selected variant.", &
697 : enum_c_vals=s2a("NONE", "ADMM1", "ADMM2", "ADMMS", "ADMMP", "ADMMQ"), &
698 : enum_desc=s2a("No short name is used, use specific definitions (default)", &
699 : "ADMM1 method from Guidon2010", &
700 : "ADMM2 method from Guidon2010", &
701 : "ADMMS method from Merlot2014", &
702 : "ADMMP method from Merlot2014", &
703 : "ADMMQ method from Merlot2014"), &
704 : enum_i_vals=[no_admm_type, admm1_type, admm2_type, admms_type, admmp_type, admmq_type], &
705 : default_i_val=no_admm_type, &
706 30834 : citations=[Guidon2010, Merlot2014])
707 10278 : CALL section_add_keyword(section, keyword)
708 10278 : CALL keyword_release(keyword)
709 :
710 : CALL keyword_create( &
711 : keyword, __LOCATION__, &
712 : name="ADMM_PURIFICATION_METHOD", &
713 : description="Method that shall be used for wavefunction fitting. Use MO_DIAG for MD.", &
714 : enum_c_vals=s2a("NONE", "CAUCHY", "CAUCHY_SUBSPACE", "MO_DIAG", "MO_NO_DIAG", "MCWEENY", "NONE_DM"), &
715 : enum_i_vals=[do_admm_purify_none, do_admm_purify_cauchy, do_admm_purify_cauchy_subspace, &
716 : do_admm_purify_mo_diag, do_admm_purify_mo_no_diag, &
717 : do_admm_purify_mcweeny, do_admm_purify_none_dm], &
718 : enum_desc=s2a("Do not apply any purification", &
719 : "Perform purification via general Cauchy representation", &
720 : "Perform purification via Cauchy representation in occupied subspace", &
721 : "Calculate MO derivatives via Cauchy representation by diagonalization", &
722 : "Calculate MO derivatives via Cauchy representation by inversion", &
723 : "Perform original McWeeny purification via matrix multiplications", &
724 : "Do not apply any purification, works directly with density matrix"), &
725 10278 : default_i_val=do_admm_purify_mo_diag)
726 10278 : CALL section_add_keyword(section, keyword)
727 10278 : CALL keyword_release(keyword)
728 :
729 : CALL keyword_create( &
730 : keyword, __LOCATION__, &
731 : name="METHOD", &
732 : description="Method that shall be used for wavefunction fitting. Use BASIS_PROJECTION for MD.", &
733 : enum_c_vals=s2a("BASIS_PROJECTION", "BLOCKED_PROJECTION_PURIFY_FULL", "BLOCKED_PROJECTION", &
734 : "CHARGE_CONSTRAINED_PROJECTION"), &
735 : enum_i_vals=[do_admm_basis_projection, do_admm_blocking_purify_full, do_admm_blocked_projection, &
736 : do_admm_charge_constrained_projection], &
737 : enum_desc=s2a("Construct auxiliary density matrix from auxiliary basis.", &
738 : "Construct auxiliary density from a blocked Fock matrix,"// &
739 : " but use the original matrix for purification.", &
740 : "Construct auxiliary density from a blocked Fock matrix.", &
741 : "Construct auxiliary density from auxiliary basis enforcing charge constrain."), &
742 10278 : default_i_val=do_admm_basis_projection)
743 10278 : CALL section_add_keyword(section, keyword)
744 10278 : CALL keyword_release(keyword)
745 :
746 : CALL keyword_create( &
747 : keyword, __LOCATION__, &
748 : name="EXCH_SCALING_MODEL", &
749 : description="Scaling of the exchange correction calculated by the auxiliary density matrix.", &
750 : enum_c_vals=s2a("NONE", "MERLOT"), &
751 : enum_i_vals=[do_admm_exch_scaling_none, do_admm_exch_scaling_merlot], &
752 : enum_desc=s2a("No scaling is enabled, refers to methods ADMM1, ADMM2 or ADMMQ.", &
753 : "Exchange scaling according to Merlot (2014)"), &
754 10278 : default_i_val=do_admm_exch_scaling_none)
755 10278 : CALL section_add_keyword(section, keyword)
756 10278 : CALL keyword_release(keyword)
757 :
758 : CALL keyword_create( &
759 : keyword, __LOCATION__, &
760 : name="EXCH_CORRECTION_FUNC", &
761 : description="Exchange functional used for the ADMM correction. It should be chosen consistently "// &
762 : "with the exchange functional in the main XC setup. LibXC implementations require linking with LibXC.", &
763 : enum_c_vals=s2a("DEFAULT", "PBEX", "NONE", "OPTX", "BECKE88X", &
764 : "PBEX_LIBXC", "BECKE88X_LIBXC", "OPTX_LIBXC", "DEFAULT_LIBXC", "LDA_X_LIBXC"), &
765 : enum_i_vals=[do_admm_aux_exch_func_default, do_admm_aux_exch_func_pbex, &
766 : do_admm_aux_exch_func_none, do_admm_aux_exch_func_opt, do_admm_aux_exch_func_bee, &
767 : do_admm_aux_exch_func_pbex_libxc, do_admm_aux_exch_func_bee_libxc, &
768 : do_admm_aux_exch_func_opt_libxc, do_admm_aux_exch_func_default_libxc, &
769 : do_admm_aux_exch_func_sx_libxc], &
770 : enum_desc=s2a("Use PBE-based corrections according to the chosen interaction operator.", &
771 : "Use PBEX functional for exchange correction.", &
772 : "No correction: X(D)-x(d)-> 0.", &
773 : "Use OPTX functional for exchange correction.", &
774 : "Use Becke88X functional for exchange correction.", &
775 : "Use PBEX functional (LibXC implementation) for exchange correction.", &
776 : "Use Becke88X functional (LibXC implementation) for exchange correction.", &
777 : "Use OPTX functional (LibXC implementation) for exchange correction.", &
778 : "Use PBE-based corrections (LibXC where possible) to the chosen interaction operator.", &
779 : "Use Slater X functional (LibXC where possible) for exchange correction."), &
780 10278 : default_i_val=do_admm_aux_exch_func_default)
781 10278 : CALL section_add_keyword(section, keyword)
782 10278 : CALL keyword_release(keyword)
783 :
784 : CALL keyword_create(keyword, __LOCATION__, name="optx_a1", &
785 : description="OPTX a1 coefficient", &
786 10278 : default_r_val=1.05151_dp)
787 10278 : CALL section_add_keyword(section, keyword)
788 10278 : CALL keyword_release(keyword)
789 : CALL keyword_create(keyword, __LOCATION__, name="optx_a2", &
790 : description="OPTX a2 coefficient", &
791 10278 : default_r_val=1.43169_dp)
792 10278 : CALL section_add_keyword(section, keyword)
793 10278 : CALL keyword_release(keyword)
794 : CALL keyword_create(keyword, __LOCATION__, name="optx_gamma", &
795 : description="OPTX gamma coefficient", &
796 10278 : default_r_val=0.006_dp)
797 10278 : CALL section_add_keyword(section, keyword)
798 10278 : CALL keyword_release(keyword)
799 :
800 : CALL keyword_create(keyword, __LOCATION__, name="BLOCK_LIST", &
801 : description="Specifies a list of atoms.", &
802 : usage="BLOCK_LIST {integer} {integer} .. {integer}", &
803 10278 : n_var=-1, type_of_var=integer_t, repeats=.TRUE.)
804 10278 : CALL section_add_keyword(section, keyword)
805 10278 : CALL keyword_release(keyword)
806 :
807 : CALL keyword_create(keyword, __LOCATION__, name="EPS_FILTER", &
808 : description="Define accuracy of DBCSR operations", &
809 10278 : usage="EPS_FILTER", default_r_val=0.0_dp)
810 10278 : CALL section_add_keyword(section, keyword)
811 10278 : CALL keyword_release(keyword)
812 :
813 10278 : END SUBROUTINE create_admm_section
814 :
815 : ! **************************************************************************************************
816 : !> \brief ...
817 : !> \param section ...
818 : ! **************************************************************************************************
819 10278 : SUBROUTINE create_density_fitting_section(section)
820 : TYPE(section_type), POINTER :: section
821 :
822 : TYPE(keyword_type), POINTER :: keyword
823 : TYPE(section_type), POINTER :: print_key
824 :
825 10278 : NULLIFY (keyword, print_key)
826 10278 : CPASSERT(.NOT. ASSOCIATED(section))
827 : CALL section_create(section, __LOCATION__, name="DENSITY_FITTING", &
828 : description="Setup parameters for density fitting (Bloechl charges or density derived "// &
829 : "atomic point charges (DDAPC) charges)", &
830 : n_keywords=7, n_subsections=0, repeats=.FALSE., &
831 20556 : citations=[Blochl1995])
832 :
833 : CALL keyword_create(keyword, __LOCATION__, name="NUM_GAUSS", &
834 : description="Specifies the numbers of gaussian used to fit the QM density for each atomic site.", &
835 : usage="NUM_GAUSS {integer}", &
836 10278 : n_var=1, type_of_var=integer_t, default_i_val=3)
837 10278 : CALL section_add_keyword(section, keyword)
838 10278 : CALL keyword_release(keyword)
839 :
840 : CALL keyword_create(keyword, __LOCATION__, name="PFACTOR", &
841 : description="Specifies the progression factor for the gaussian exponent for each atomic site.", &
842 : usage="PFACTOR {real}", &
843 10278 : n_var=1, type_of_var=real_t, default_r_val=1.5_dp)
844 10278 : CALL section_add_keyword(section, keyword)
845 10278 : CALL keyword_release(keyword)
846 :
847 : CALL keyword_create(keyword, __LOCATION__, name="MIN_RADIUS", &
848 : description="Specifies the smallest radius of the gaussian used in the fit. All other radius are"// &
849 : " obtained with the progression factor.", &
850 : usage="MIN_RADIUS {real}", &
851 10278 : unit_str="angstrom", n_var=1, type_of_var=real_t, default_r_val=0.5_dp)
852 10278 : CALL section_add_keyword(section, keyword)
853 10278 : CALL keyword_release(keyword)
854 :
855 : CALL keyword_create(keyword, __LOCATION__, name="RADII", &
856 : description="Specifies all the radius of the gaussian used in the fit for each atomic site. The use"// &
857 : " of this keyword disables all other keywords of this section.", &
858 : usage="RADII {real} {real} .. {real}", &
859 10278 : unit_str="angstrom", n_var=-1, type_of_var=real_t)
860 10278 : CALL section_add_keyword(section, keyword)
861 10278 : CALL keyword_release(keyword)
862 :
863 : CALL keyword_create(keyword, __LOCATION__, name="GCUT", &
864 : description="Cutoff for charge fit in G-space.", &
865 : usage="GCUT {real}", &
866 10278 : n_var=1, type_of_var=real_t, default_r_val=SQRT(6.0_dp))
867 10278 : CALL section_add_keyword(section, keyword)
868 10278 : CALL keyword_release(keyword)
869 :
870 : CALL cp_print_key_section_create(print_key, __LOCATION__, "program_run_info", &
871 : description="Controls the printing of basic information during the run", &
872 10278 : print_level=high_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
873 :
874 : CALL keyword_create(keyword, __LOCATION__, name="CONDITION_NUMBER", &
875 : description="Prints information regarding the condition numbers of the A matrix (to be inverted)", &
876 : usage="CONDITION_NUMBER <LOGICAL>", &
877 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
878 10278 : CALL section_add_keyword(print_key, keyword)
879 10278 : CALL keyword_release(keyword)
880 :
881 10278 : CALL section_add_subsection(section, print_key)
882 10278 : CALL section_release(print_key)
883 :
884 10278 : END SUBROUTINE create_density_fitting_section
885 :
886 : ! **************************************************************************************************
887 : !> \brief creates the input section for the relativistic part
888 : !> \param section the section to create
889 : !> \author jens
890 : ! **************************************************************************************************
891 10278 : SUBROUTINE create_relativistic_section(section)
892 : TYPE(section_type), POINTER :: section
893 :
894 : TYPE(keyword_type), POINTER :: keyword
895 :
896 10278 : CPASSERT(.NOT. ASSOCIATED(section))
897 : CALL section_create(section, __LOCATION__, name="relativistic", &
898 : description="parameters needed and setup for relativistic calculations", &
899 10278 : n_keywords=5, n_subsections=0, repeats=.FALSE.)
900 :
901 10278 : NULLIFY (keyword)
902 :
903 : CALL keyword_create(keyword, __LOCATION__, name="method", &
904 : description="type of relativistic correction used", &
905 : usage="method (NONE|DKH|ZORA)", default_i_val=rel_none, &
906 : enum_c_vals=s2a("NONE", "DKH", "ZORA"), &
907 : enum_i_vals=[rel_none, rel_dkh, rel_zora], &
908 : enum_desc=s2a("Use no relativistic correction", &
909 : "Use Douglas-Kroll-Hess method", &
910 10278 : "Use ZORA method"))
911 10278 : CALL section_add_keyword(section, keyword)
912 10278 : CALL keyword_release(keyword)
913 :
914 : CALL keyword_create(keyword, __LOCATION__, name="DKH_order", &
915 : description="The order of the DKH transformation ", &
916 10278 : usage="DKH_order 2", default_i_val=2)
917 10278 : CALL section_add_keyword(section, keyword)
918 10278 : CALL keyword_release(keyword)
919 :
920 : CALL keyword_create(keyword, __LOCATION__, name="ZORA_type", &
921 : description="Type of ZORA method to be used", &
922 : usage="ZORA_type scMP", default_i_val=rel_zora_full, &
923 : enum_c_vals=s2a("FULL", "MP", "scMP"), &
924 : enum_desc=s2a("Full ZORA method (not implemented)", &
925 : "ZORA with atomic model potential", &
926 : "Scaled ZORA with atomic model potential"), &
927 10278 : enum_i_vals=[rel_zora_full, rel_zora_mp, rel_sczora_mp])
928 10278 : CALL section_add_keyword(section, keyword)
929 10278 : CALL keyword_release(keyword)
930 :
931 : CALL keyword_create(keyword, __LOCATION__, name="transformation", &
932 : description="Type of DKH transformation", &
933 : usage="transformation (FULL|MOLECULE|ATOM)", default_i_val=rel_trans_atom, &
934 : enum_c_vals=s2a("FULL", "MOLECULE", "ATOM"), &
935 : enum_i_vals=[rel_trans_full, rel_trans_molecule, rel_trans_atom], &
936 : enum_desc=s2a("Use full matrix transformation", &
937 : "Use transformation blocked by molecule", &
938 10278 : "Use atomic blocks"))
939 10278 : CALL section_add_keyword(section, keyword)
940 10278 : CALL keyword_release(keyword)
941 :
942 : CALL keyword_create(keyword, __LOCATION__, name="z_cutoff", &
943 : description="The minimal atomic number considered for atom transformation", &
944 10278 : usage="z_cutoff 50", default_i_val=1)
945 10278 : CALL section_add_keyword(section, keyword)
946 10278 : CALL keyword_release(keyword)
947 :
948 : CALL keyword_create(keyword, __LOCATION__, name="potential", &
949 : description="External potential used in DKH transformation, full 1/r or erfc(r)/r", &
950 : usage="POTENTIAL {FULL,ERFC}", default_i_val=rel_pot_erfc, &
951 : enum_c_vals=s2a("FULL", "ERFC"), &
952 10278 : enum_i_vals=[rel_pot_full, rel_pot_erfc])
953 10278 : CALL section_add_keyword(section, keyword)
954 10278 : CALL keyword_release(keyword)
955 :
956 10278 : END SUBROUTINE create_relativistic_section
957 :
958 : ! **************************************************************************************************
959 : !> \brief creates the KG section
960 : !> \param section ...
961 : !> \author Martin Haeufel [2012.07]
962 : ! **************************************************************************************************
963 10278 : SUBROUTINE create_kg_section(section)
964 : TYPE(section_type), POINTER :: section
965 :
966 : TYPE(keyword_type), POINTER :: keyword
967 : TYPE(section_type), POINTER :: print_key, subsection
968 :
969 10278 : CPASSERT(.NOT. ASSOCIATED(section))
970 : CALL section_create(section, __LOCATION__, name="KG_METHOD", &
971 : description="Specifies the parameters for a Kim-Gordon-like partitioning"// &
972 : " into molecular subunits", &
973 : n_keywords=0, n_subsections=1, repeats=.FALSE., &
974 41112 : citations=[Iannuzzi2006, Brelaz1979, Andermatt2016])
975 :
976 10278 : NULLIFY (keyword, subsection, print_key)
977 :
978 : ! add a XC section
979 10278 : CALL create_xc_section(subsection)
980 10278 : CALL section_add_subsection(section, subsection)
981 10278 : CALL section_release(subsection)
982 :
983 : ! add LRI section
984 10278 : CALL create_lrigpw_section(subsection)
985 10278 : CALL section_add_subsection(section, subsection)
986 10278 : CALL section_release(subsection)
987 :
988 : CALL keyword_create(keyword, __LOCATION__, name="COLORING_METHOD", &
989 : description="Which algorithm to use for coloring.", &
990 : usage="COLORING_METHOD GREEDY", &
991 : default_i_val=kg_color_dsatur, &
992 : enum_c_vals=s2a("DSATUR", "GREEDY"), &
993 : enum_desc=s2a("Maximum degree of saturation, relatively accurate", &
994 : "Greedy, fast coloring, less accurate"), &
995 10278 : enum_i_vals=[kg_color_dsatur, kg_color_greedy])
996 10278 : CALL section_add_keyword(section, keyword)
997 10278 : CALL keyword_release(keyword)
998 :
999 : CALL keyword_create(keyword, __LOCATION__, name="TNADD_METHOD", &
1000 : description="Algorithm to use for the calculation of the nonadditive kinetic energy.", &
1001 : usage="TNADD_METHOD ATOMIC", &
1002 : default_i_val=kg_tnadd_embed, &
1003 : enum_c_vals=s2a("EMBEDDING", "RI_EMBEDDING", "ATOMIC", "NONE"), &
1004 : enum_desc=s2a("Use full embedding potential (see Iannuzzi et al)", &
1005 : "Use full embedding potential with RI density fitting", &
1006 : "Use sum of atomic model potentials", &
1007 : "Do not use kinetic energy embedding"), &
1008 10278 : enum_i_vals=[kg_tnadd_embed, kg_tnadd_embed_ri, kg_tnadd_atomic, kg_tnadd_none])
1009 10278 : CALL section_add_keyword(section, keyword)
1010 10278 : CALL keyword_release(keyword)
1011 :
1012 : CALL keyword_create(keyword, __LOCATION__, name="INTEGRATION_GRID", &
1013 : description="Grid [small,medium,large,huge]to be used for the TNADD integration.", &
1014 : usage="INTEGRATION_GRID MEDIUM", &
1015 10278 : default_c_val="MEDIUM")
1016 10278 : CALL section_add_keyword(section, keyword)
1017 10278 : CALL keyword_release(keyword)
1018 :
1019 : CALL section_create(subsection, __LOCATION__, name="PRINT", &
1020 : description="Print section", &
1021 10278 : n_keywords=0, n_subsections=1, repeats=.FALSE.)
1022 :
1023 : CALL cp_print_key_section_create(print_key, __LOCATION__, "NEIGHBOR_LISTS", &
1024 : description="Controls the printing of the neighbor lists.", &
1025 10278 : print_level=low_print_level, filename="__STD_OUT__", unit_str="angstrom")
1026 :
1027 : CALL keyword_create(keyword, __LOCATION__, &
1028 : name="SAB_ORB_FULL", &
1029 : description="Activates the printing of the full orbital "// &
1030 : "orbital neighbor lists.", &
1031 : default_l_val=.FALSE., &
1032 10278 : lone_keyword_l_val=.TRUE.)
1033 10278 : CALL section_add_keyword(print_key, keyword)
1034 10278 : CALL keyword_release(keyword)
1035 :
1036 : CALL keyword_create(keyword, __LOCATION__, &
1037 : name="SAB_ORB_MOLECULAR", &
1038 : description="Activates the printing of the orbital "// &
1039 : "orbital neighbor lists for molecular subsets.", &
1040 : default_l_val=.FALSE., &
1041 10278 : lone_keyword_l_val=.TRUE.)
1042 10278 : CALL section_add_keyword(print_key, keyword)
1043 10278 : CALL keyword_release(keyword)
1044 :
1045 : CALL keyword_create(keyword, __LOCATION__, &
1046 : name="SAC_KIN", &
1047 : description="Activates the printing of the orbital "// &
1048 : "atomic potential neighbor list.", &
1049 : default_l_val=.FALSE., &
1050 10278 : lone_keyword_l_val=.TRUE.)
1051 10278 : CALL section_add_keyword(print_key, keyword)
1052 10278 : CALL keyword_release(keyword)
1053 :
1054 10278 : CALL section_add_subsection(subsection, print_key)
1055 10278 : CALL section_release(print_key)
1056 :
1057 10278 : CALL section_add_subsection(section, subsection)
1058 10278 : CALL section_release(subsection)
1059 :
1060 10278 : END SUBROUTINE create_kg_section
1061 :
1062 : ! **************************************************************************************************
1063 : !> \brief Create the BSSE section for counterpoise correction
1064 : !> \param section the section to create
1065 : !> \author teo
1066 : ! **************************************************************************************************
1067 10262 : SUBROUTINE create_bsse_section(section)
1068 : TYPE(section_type), POINTER :: section
1069 :
1070 : TYPE(keyword_type), POINTER :: keyword
1071 : TYPE(section_type), POINTER :: subsection
1072 :
1073 10262 : CPASSERT(.NOT. ASSOCIATED(section))
1074 : CALL section_create(section, __LOCATION__, name="BSSE", &
1075 : description="This section is used to set up the BSSE calculation. "// &
1076 : "It also requires that for each atomic kind X a kind X_ghost is present, "// &
1077 : "with the GHOST keyword specified, in addition to the other required fields.", &
1078 10262 : n_keywords=3, n_subsections=1, repeats=.FALSE.)
1079 :
1080 10262 : NULLIFY (keyword, subsection)
1081 : ! FRAGMENT SECTION
1082 : CALL section_create(subsection, __LOCATION__, name="FRAGMENT", &
1083 : description="Specify the atom number belonging to this fragment.", &
1084 10262 : n_keywords=2, n_subsections=0, repeats=.TRUE.)
1085 :
1086 : CALL keyword_create(keyword, __LOCATION__, name="LIST", &
1087 : description="Specifies a list of atoms.", &
1088 : usage="LIST {integer} {integer} .. {integer}", &
1089 10262 : repeats=.TRUE., n_var=-1, type_of_var=integer_t)
1090 10262 : CALL section_add_keyword(subsection, keyword)
1091 10262 : CALL keyword_release(keyword)
1092 :
1093 10262 : CALL section_add_subsection(section, subsection)
1094 10262 : CALL section_release(subsection)
1095 :
1096 : ! CONFIGURATION SECTION
1097 : CALL section_create(subsection, __LOCATION__, name="CONFIGURATION", &
1098 : description="Specify additional parameters for the combinatorial configurations. "// &
1099 : "Use this section to manually specify charge and multiplicity of the fragments "// &
1100 : "and their combinations.", &
1101 10262 : n_keywords=2, n_subsections=0, repeats=.TRUE.)
1102 :
1103 : CALL keyword_create(keyword, __LOCATION__, name="GLB_CONF", &
1104 : description="Specifies the global configuration using 1 or 0 for each fragment. "// &
1105 : "1 specifies the respective fragment as used, 0 as unused.", &
1106 : usage="GLB_CONF {integer} {integer} .. {integer}", &
1107 10262 : n_var=-1, type_of_var=integer_t)
1108 10262 : CALL section_add_keyword(subsection, keyword)
1109 10262 : CALL keyword_release(keyword)
1110 :
1111 : CALL keyword_create(keyword, __LOCATION__, name="SUB_CONF", &
1112 : description="Specifies the subconfiguration using 1 or 0 belonging to the global configuration. "// &
1113 : "1 specifies the respective fragment as real, 0 as ghost.", &
1114 : usage="SUB_CONF {integer} {integer} .. {integer}", &
1115 10262 : n_var=-1, type_of_var=integer_t)
1116 10262 : CALL section_add_keyword(subsection, keyword)
1117 10262 : CALL keyword_release(keyword)
1118 :
1119 : CALL keyword_create(keyword, __LOCATION__, &
1120 : name="MULTIPLICITY", &
1121 : variants=["MULTIP"], &
1122 : description="Specify for each fragment the multiplicity. Two times the total spin plus one. "// &
1123 : "Specify 3 for a triplet, 4 for a quartet,and so on. Default is 1 (singlet) for an "// &
1124 : "even number and 2 (doublet) for an odd number of electrons.", &
1125 : usage="MULTIPLICITY 3", &
1126 20524 : default_i_val=0) ! this default value is just a flag to get the above
1127 10262 : CALL section_add_keyword(subsection, keyword)
1128 10262 : CALL keyword_release(keyword)
1129 :
1130 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE", &
1131 : description="The total charge for each fragment.", &
1132 : usage="CHARGE -1", &
1133 10262 : default_i_val=0)
1134 10262 : CALL section_add_keyword(subsection, keyword)
1135 10262 : CALL keyword_release(keyword)
1136 10262 : CALL section_add_subsection(section, subsection)
1137 10262 : CALL section_release(subsection)
1138 :
1139 : CALL section_create(subsection, __LOCATION__, name="FRAGMENT_ENERGIES", &
1140 : description="This section contains the energies of the fragments already"// &
1141 : " computed. It is useful as a summary and specifically for restarting BSSE runs.", &
1142 10262 : n_keywords=2, n_subsections=0, repeats=.TRUE.)
1143 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1144 : description="The energy computed for each fragment", repeats=.TRUE., &
1145 10262 : usage="{REAL}", type_of_var=real_t)
1146 10262 : CALL section_add_keyword(subsection, keyword)
1147 10262 : CALL keyword_release(keyword)
1148 10262 : CALL section_add_subsection(section, subsection)
1149 10262 : CALL section_release(subsection)
1150 :
1151 10262 : CALL create_print_bsse_section(subsection)
1152 10262 : CALL section_add_subsection(section, subsection)
1153 10262 : CALL section_release(subsection)
1154 :
1155 10262 : END SUBROUTINE create_bsse_section
1156 :
1157 : ! **************************************************************************************************
1158 : !> \brief Create the print bsse section
1159 : !> \param section the section to create
1160 : !> \author teo
1161 : ! **************************************************************************************************
1162 10262 : SUBROUTINE create_print_bsse_section(section)
1163 : TYPE(section_type), POINTER :: section
1164 :
1165 : TYPE(section_type), POINTER :: print_key
1166 :
1167 10262 : CPASSERT(.NOT. ASSOCIATED(section))
1168 : CALL section_create(section, __LOCATION__, name="print", &
1169 : description="Section of possible print options in BSSE code.", &
1170 10262 : n_keywords=0, n_subsections=1, repeats=.FALSE.)
1171 :
1172 10262 : NULLIFY (print_key)
1173 : CALL cp_print_key_section_create(print_key, __LOCATION__, "PROGRAM_RUN_INFO", &
1174 : description="Controls the printing of information regarding the run.", &
1175 10262 : print_level=low_print_level, filename="__STD_OUT__")
1176 10262 : CALL section_add_subsection(section, print_key)
1177 10262 : CALL section_release(print_key)
1178 :
1179 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RESTART", &
1180 : description="Controls the dumping of the restart file during BSSE runs. "// &
1181 : "By default the restart is updated after each configuration calculation is "// &
1182 : "completed.", &
1183 : print_level=silent_print_level, common_iter_levels=0, &
1184 10262 : add_last=add_last_numeric, filename="")
1185 10262 : CALL section_add_subsection(section, print_key)
1186 10262 : CALL section_release(print_key)
1187 :
1188 10262 : END SUBROUTINE create_print_bsse_section
1189 :
1190 : ! **************************************************************************************************
1191 : !> \brief input section for optional parameters for RIGPW
1192 : !> \param section the section to create
1193 : !> \author JGH [06.2017]
1194 : ! **************************************************************************************************
1195 0 : SUBROUTINE create_rigpw_section(section)
1196 : TYPE(section_type), POINTER :: section
1197 :
1198 0 : CPASSERT(.NOT. ASSOCIATED(section))
1199 : CALL section_create(section, __LOCATION__, name="RIGPW", &
1200 : description="This section specifies optional parameters for RIGPW.", &
1201 0 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1202 :
1203 : ! CALL keyword_create(keyword, __LOCATION__, name="RI_OVERLAP_MATRIX", &
1204 : ! description="Specifies whether to calculate the inverse or the "// &
1205 : ! "pseudoinverse of the overlap matrix of the auxiliary "// &
1206 : ! "basis set. Calculating the pseudoinverse is necessary "// &
1207 : ! "for very large auxiliary basis sets, but more expensive. "// &
1208 : ! "Using the pseudoinverse, consistent forces are not "// &
1209 : ! "guaranteed yet.", &
1210 : ! usage="RI_OVERLAP_MATRIX INVERSE", &
1211 : ! enum_c_vals=s2a("INVERSE", "PSEUDO_INVERSE_SVD", "PSEUDO_INVERSE_DIAG", &
1212 : ! "AUTOSELECT"), &
1213 : ! enum_desc=s2a("Calculate inverse of the overlap matrix.", &
1214 : ! "Calculate the pseuodinverse of the overlap matrix "// &
1215 : ! "using singular value decomposition.", &
1216 : ! "Calculate the pseudoinverse of the overlap matrix "// &
1217 : ! "by prior diagonalization.", &
1218 : ! "Choose automatically for each pair whether to "// &
1219 : ! "calculate the inverse or pseudoinverse based on the "// &
1220 : ! "condition number of the overlap matrix for each pair. "// &
1221 : ! "Calculating the pseudoinverse is much more expensive."), &
1222 : ! enum_i_vals=(/do_lri_inv, do_lri_pseudoinv_svd, &
1223 : ! do_lri_pseudoinv_diag, do_lri_inv_auto/), &
1224 : ! default_i_val=do_lri_inv)
1225 : ! CALL section_add_keyword(section, keyword)
1226 : ! CALL keyword_release(keyword)
1227 :
1228 0 : END SUBROUTINE create_rigpw_section
1229 :
1230 : ! **************************************************************************************************
1231 : !> \brief creates the multigrid
1232 : !> \param section input section to create
1233 : !> \param create_subsections indicates whether or not subsections INTERPOLATOR and RS_GRID
1234 : !> should be created
1235 : !> \author fawzi
1236 : ! **************************************************************************************************
1237 30802 : SUBROUTINE create_mgrid_section(section, create_subsections)
1238 : TYPE(section_type), POINTER :: section
1239 : LOGICAL, INTENT(in) :: create_subsections
1240 :
1241 : TYPE(keyword_type), POINTER :: keyword
1242 : TYPE(section_type), POINTER :: subsection
1243 :
1244 30802 : CPASSERT(.NOT. ASSOCIATED(section))
1245 : CALL section_create(section, __LOCATION__, name="mgrid", &
1246 : description="Controls the multigrid used by GPW/GAPW to represent densities, "// &
1247 : "potentials, and Gaussian products on real-space grids.", &
1248 30802 : n_keywords=5, n_subsections=1, repeats=.FALSE.)
1249 30802 : NULLIFY (keyword)
1250 : CALL keyword_create(keyword, __LOCATION__, name="NGRIDS", &
1251 : description="Number of multigrid levels. Smooth Gaussian products can be mapped to "// &
1252 : "coarser levels, while sharper products require finer levels.", &
1253 30802 : usage="ngrids 1", default_i_val=4)
1254 30802 : CALL section_add_keyword(section, keyword)
1255 30802 : CALL keyword_release(keyword)
1256 :
1257 : CALL keyword_create(keyword, __LOCATION__, name="cutoff", &
1258 : description= &
1259 : "Plane-wave cutoff of the finest real-space grid level. "// &
1260 : "Increasing this value improves the grid representation, but it is "// &
1261 : "not a substitute for converging the Gaussian basis set. "// &
1262 : "Default value for SE or DFTB calculation is 1.0 [Ry].", &
1263 : usage="cutoff 300", &
1264 : default_r_val=cp_unit_to_cp2k(value=280.0_dp, unit_str="Ry"), &
1265 30802 : n_var=1, unit_str="Ry")
1266 30802 : CALL section_add_keyword(section, keyword)
1267 30802 : CALL keyword_release(keyword)
1268 :
1269 : CALL keyword_create(keyword, __LOCATION__, name="progression_factor", &
1270 : description="Factor used to derive the cutoff of coarser multigrid levels when "// &
1271 : "they are not given explicitly.", &
1272 30802 : usage="progression_factor <integer>", default_r_val=3._dp)
1273 30802 : CALL section_add_keyword(section, keyword)
1274 30802 : CALL keyword_release(keyword)
1275 :
1276 : CALL keyword_create(keyword, __LOCATION__, name="commensurate", &
1277 : description="If the grids should be commensurate. If true overrides "// &
1278 : "the progression factor and the cutoffs of the sub grids", &
1279 : usage="commensurate", default_l_val=.FALSE., &
1280 30802 : lone_keyword_l_val=.TRUE.)
1281 30802 : CALL section_add_keyword(section, keyword)
1282 30802 : CALL keyword_release(keyword)
1283 :
1284 : CALL keyword_create(keyword, __LOCATION__, name="realspace", &
1285 : description="If both rho and rho_gspace are needed ", &
1286 : usage="realspace", default_l_val=.FALSE., &
1287 30802 : lone_keyword_l_val=.TRUE.)
1288 30802 : CALL section_add_keyword(section, keyword)
1289 30802 : CALL keyword_release(keyword)
1290 :
1291 : CALL keyword_create(keyword, __LOCATION__, name="REL_CUTOFF", &
1292 : variants=["RELATIVE_CUTOFF"], &
1293 : description="Controls to which multigrid level a Gaussian product is mapped. "// &
1294 : "It is the reference cutoff for a Gaussian with exponent alpha=1. Larger values "// &
1295 : "keep more Gaussian products on finer grids and can be important for accurate "// &
1296 : "energies, forces, stress tensors, and variable-cell simulations.", &
1297 : usage="RELATIVE_CUTOFF real", default_r_val=20.0_dp, &
1298 61604 : unit_str="Ry")
1299 30802 : CALL section_add_keyword(section, keyword)
1300 30802 : CALL keyword_release(keyword)
1301 :
1302 : CALL keyword_create(keyword, __LOCATION__, name="MULTIGRID_SET", &
1303 : description="Activate a manual setting of the multigrids", &
1304 30802 : usage="MULTIGRID_SET", default_l_val=.FALSE.)
1305 30802 : CALL section_add_keyword(section, keyword)
1306 30802 : CALL keyword_release(keyword)
1307 :
1308 : CALL keyword_create(keyword, __LOCATION__, &
1309 : name="SKIP_LOAD_BALANCE_DISTRIBUTED", &
1310 : description="Skips load balancing on distributed multigrids. "// &
1311 : "Memory usage is O(p) so may be used "// &
1312 : "for all but the very largest runs.", &
1313 : usage="SKIP_LOAD_BALANCE_DISTRIBUTED", &
1314 : default_l_val=.FALSE., &
1315 30802 : lone_keyword_l_val=.TRUE.)
1316 : ! CALL keyword_create(keyword, __LOCATION__, name="SKIP_LOAD_BALANCE_DISTRIBUTED",&
1317 : ! description="Skip load balancing on distributed multigrids, which might be memory intensive."//&
1318 : ! "If not explicitly specified, runs using more than 1024 MPI tasks will default to .TRUE.",&
1319 : ! usage="SKIP_LOAD_BALANCE_DISTRIBUTED", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1320 :
1321 30802 : CALL section_add_keyword(section, keyword)
1322 30802 : CALL keyword_release(keyword)
1323 :
1324 : CALL keyword_create(keyword, __LOCATION__, name="MULTIGRID_CUTOFF", &
1325 : variants=["CUTOFF_LIST"], &
1326 : description="List of cutoff values to set up multigrids manually", &
1327 : usage="MULTIGRID_CUTOFF 200.0 100.0 ", &
1328 : n_var=-1, &
1329 : type_of_var=real_t, &
1330 61604 : unit_str="Ry")
1331 30802 : CALL section_add_keyword(section, keyword)
1332 30802 : CALL keyword_release(keyword)
1333 :
1334 30802 : IF (create_subsections) THEN
1335 10278 : NULLIFY (subsection)
1336 10278 : CALL create_rsgrid_section(subsection)
1337 10278 : CALL section_add_subsection(section, subsection)
1338 10278 : CALL section_release(subsection)
1339 :
1340 10278 : NULLIFY (subsection)
1341 10278 : CALL create_interp_section(subsection)
1342 10278 : CALL section_add_subsection(section, subsection)
1343 10278 : CALL section_release(subsection)
1344 : END IF
1345 30802 : END SUBROUTINE create_mgrid_section
1346 :
1347 : ! **************************************************************************************************
1348 : !> \brief creates the interpolation section
1349 : !> \param section ...
1350 : !> \author tlaino
1351 : ! **************************************************************************************************
1352 82112 : SUBROUTINE create_interp_section(section)
1353 : TYPE(section_type), POINTER :: section
1354 :
1355 : TYPE(keyword_type), POINTER :: keyword
1356 : TYPE(section_type), POINTER :: print_key
1357 :
1358 82112 : CPASSERT(.NOT. ASSOCIATED(section))
1359 : CALL section_create(section, __LOCATION__, name="interpolator", &
1360 : description="kind of interpolation used between the multigrids", &
1361 82112 : n_keywords=5, n_subsections=0, repeats=.FALSE.)
1362 :
1363 82112 : NULLIFY (keyword, print_key)
1364 :
1365 : CALL keyword_create(keyword, __LOCATION__, name="kind", &
1366 : description="the interpolator to use", &
1367 : usage="kind spline3", &
1368 : default_i_val=pw_interp, &
1369 : enum_c_vals=s2a("pw", "spline3_nopbc", "spline3"), &
1370 : enum_i_vals=[pw_interp, &
1371 82112 : spline3_nopbc_interp, spline3_pbc_interp])
1372 82112 : CALL section_add_keyword(section, keyword)
1373 82112 : CALL keyword_release(keyword)
1374 :
1375 : CALL keyword_create(keyword, __LOCATION__, name="safe_computation", &
1376 : description="if a non unrolled calculation is to be performed in parallel", &
1377 : usage="safe_computation OFF", &
1378 : default_l_val=.FALSE., &
1379 82112 : lone_keyword_l_val=.TRUE.)
1380 82112 : CALL section_add_keyword(section, keyword)
1381 82112 : CALL keyword_release(keyword)
1382 :
1383 : CALL keyword_create(keyword, __LOCATION__, name="aint_precond", &
1384 : description="the approximate inverse to use to get the starting point"// &
1385 : " for the linear solver of the spline3 methods", &
1386 : usage="aint_precond copy", &
1387 : default_i_val=precond_spl3_aint, &
1388 : enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_aint2", &
1389 : "spl3_nopbc_precond1", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1390 : enum_i_vals=[no_precond, precond_spl3_aint, precond_spl3_aint2, &
1391 82112 : precond_spl3_1, precond_spl3_2, precond_spl3_3])
1392 82112 : CALL section_add_keyword(section, keyword)
1393 82112 : CALL keyword_release(keyword)
1394 :
1395 : CALL keyword_create(keyword, __LOCATION__, name="precond", &
1396 : description="The preconditioner used"// &
1397 : " for the linear solver of the spline3 methods", &
1398 : usage="PRECOND copy", &
1399 : default_i_val=precond_spl3_3, &
1400 : enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_aint2", &
1401 : "spl3_nopbc_precond1", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1402 : enum_i_vals=[no_precond, precond_spl3_aint, precond_spl3_aint2, &
1403 82112 : precond_spl3_1, precond_spl3_2, precond_spl3_3])
1404 82112 : CALL section_add_keyword(section, keyword)
1405 82112 : CALL keyword_release(keyword)
1406 :
1407 : CALL keyword_create(keyword, __LOCATION__, name="eps_x", &
1408 : description="accuracy on the solution for spline3 the interpolators", &
1409 82112 : usage="eps_x 1.e-15", default_r_val=1.e-10_dp)
1410 82112 : CALL section_add_keyword(section, keyword)
1411 82112 : CALL keyword_release(keyword)
1412 :
1413 : CALL keyword_create(keyword, __LOCATION__, name="eps_r", &
1414 : description="accuracy on the residual for spline3 the interpolators", &
1415 82112 : usage="eps_r 1.e-15", default_r_val=1.e-10_dp)
1416 82112 : CALL section_add_keyword(section, keyword)
1417 82112 : CALL keyword_release(keyword)
1418 :
1419 : CALL keyword_create(keyword, __LOCATION__, name="max_iter", &
1420 : variants=['maxiter'], &
1421 : description="the maximum number of iterations", &
1422 164224 : usage="max_iter 200", default_i_val=100)
1423 82112 : CALL section_add_keyword(section, keyword)
1424 82112 : CALL keyword_release(keyword)
1425 :
1426 82112 : NULLIFY (print_key)
1427 : CALL cp_print_key_section_create(print_key, __LOCATION__, "conv_info", &
1428 : description="if convergence information about the linear solver"// &
1429 : " of the spline methods should be printed", &
1430 : print_level=medium_print_level, each_iter_names=s2a("SPLINE_FIND_COEFFS"), &
1431 : each_iter_values=[10], filename="__STD_OUT__", &
1432 82112 : add_last=add_last_numeric)
1433 82112 : CALL section_add_subsection(section, print_key)
1434 82112 : CALL section_release(print_key)
1435 :
1436 82112 : END SUBROUTINE create_interp_section
1437 :
1438 : ! **************************************************************************************************
1439 : !> \brief creates the sic (self interaction correction) section
1440 : !> \param section ...
1441 : !> \author fawzi
1442 : ! **************************************************************************************************
1443 10278 : SUBROUTINE create_sic_section(section)
1444 : TYPE(section_type), POINTER :: section
1445 :
1446 : TYPE(keyword_type), POINTER :: keyword
1447 :
1448 10278 : CPASSERT(.NOT. ASSOCIATED(section))
1449 : CALL section_create(section, __LOCATION__, name="sic", &
1450 : description="parameters for the self interaction correction", &
1451 : n_keywords=6, n_subsections=0, repeats=.FALSE., &
1452 41112 : citations=[VandeVondele2005b, Perdew1981, Avezac2005])
1453 :
1454 10278 : NULLIFY (keyword)
1455 :
1456 : CALL keyword_create(keyword, __LOCATION__, name="SIC_SCALING_A", &
1457 : description="Scaling of the coulomb term in sic [experimental]", &
1458 : usage="SIC_SCALING_A 0.5", &
1459 : citations=[VandeVondele2005b], &
1460 20556 : default_r_val=1.0_dp)
1461 10278 : CALL section_add_keyword(section, keyword)
1462 10278 : CALL keyword_release(keyword)
1463 :
1464 : CALL keyword_create(keyword, __LOCATION__, name="SIC_SCALING_B", &
1465 : description="Scaling of the xc term in sic [experimental]", &
1466 : usage="SIC_SCALING_B 0.5", &
1467 : citations=[VandeVondele2005b], &
1468 20556 : default_r_val=1.0_dp)
1469 10278 : CALL section_add_keyword(section, keyword)
1470 10278 : CALL keyword_release(keyword)
1471 :
1472 : CALL keyword_create(keyword, __LOCATION__, name="SIC_METHOD", &
1473 : description="Method used to remove the self interaction", &
1474 : usage="SIC_METHOD MAURI_US", &
1475 : default_i_val=sic_none, &
1476 : enum_c_vals=s2a("NONE", "MAURI_US", "MAURI_SPZ", "AD", "EXPLICIT_ORBITALS"), &
1477 : enum_i_vals=[sic_none, sic_mauri_us, sic_mauri_spz, sic_ad, sic_eo], &
1478 : enum_desc=s2a("Do not apply a sic correction", &
1479 : "Employ a (scaled) correction proposed by Mauri and co-workers"// &
1480 : " on the spin density / doublet unpaired orbital", &
1481 : "Employ a (scaled) Perdew-Zunger expression"// &
1482 : " on the spin density / doublet unpaired orbital", &
1483 : "The average density correction", &
1484 : "(scaled) Perdew-Zunger correction explicitly on a set of orbitals."), &
1485 41112 : citations=[VandeVondele2005b, Perdew1981, Avezac2005])
1486 10278 : CALL section_add_keyword(section, keyword)
1487 10278 : CALL keyword_release(keyword)
1488 :
1489 : CALL keyword_create(keyword, __LOCATION__, name="ORBITAL_SET", &
1490 : description="Type of orbitals treated with the SIC", &
1491 : usage="ORBITAL_SET ALL", &
1492 : default_i_val=sic_list_unpaired, &
1493 : enum_c_vals=s2a("UNPAIRED", "ALL"), &
1494 : enum_desc=s2a("correction for the unpaired orbitals only, requires a restricted open shell calculation", &
1495 : "correction for all orbitals, requires a LSD or ROKS calculation"), &
1496 10278 : enum_i_vals=[sic_list_unpaired, sic_list_all])
1497 10278 : CALL section_add_keyword(section, keyword)
1498 10278 : CALL keyword_release(keyword)
1499 :
1500 10278 : END SUBROUTINE create_sic_section
1501 :
1502 : ! **************************************************************************************************
1503 : !> \brief creates the low spin roks section
1504 : !> \param section ...
1505 : !> \author Joost VandeVondele
1506 : ! **************************************************************************************************
1507 10278 : SUBROUTINE create_low_spin_roks_section(section)
1508 : TYPE(section_type), POINTER :: section
1509 :
1510 : TYPE(keyword_type), POINTER :: keyword
1511 :
1512 10278 : CPASSERT(.NOT. ASSOCIATED(section))
1513 : CALL section_create(section, __LOCATION__, name="LOW_SPIN_ROKS", &
1514 : description="Specify the details of the low spin ROKS method. "// &
1515 : "In particular, one can specify various terms added to the energy of the high spin roks configuration"// &
1516 : " with a energy scaling factor, and a prescription of the spin state.", &
1517 10278 : n_keywords=6, n_subsections=0, repeats=.FALSE.)
1518 :
1519 10278 : NULLIFY (keyword)
1520 : CALL keyword_create(keyword, __LOCATION__, name="ENERGY_SCALING", &
1521 : description="The scaling factors for each term added to the total energy. "// &
1522 : "This list should contain one number for each term added to the total energy.", &
1523 : usage="ENERGY_SCALING 1.0 -1.0 ", &
1524 10278 : n_var=-1, type_of_var=real_t, repeats=.FALSE.)
1525 10278 : CALL section_add_keyword(section, keyword)
1526 10278 : CALL keyword_release(keyword)
1527 : CALL keyword_create( &
1528 : keyword, __LOCATION__, name="SPIN_CONFIGURATION", &
1529 : description="For each singly occupied orbital, specify if this should be an alpha (=1) or a beta (=2) orbital. "// &
1530 : "This keyword should be repeated, each repetition corresponding to an additional term.", &
1531 : usage="SPIN_CONFIGURATION 1 2", &
1532 10278 : n_var=-1, type_of_var=integer_t, repeats=.TRUE.)
1533 10278 : CALL section_add_keyword(section, keyword)
1534 10278 : CALL keyword_release(keyword)
1535 :
1536 10278 : END SUBROUTINE create_low_spin_roks_section
1537 :
1538 : ! **************************************************************************************************
1539 : !> \brief ...
1540 : !> \param section ...
1541 : ! **************************************************************************************************
1542 10278 : SUBROUTINE create_rtp_section(section)
1543 : TYPE(section_type), POINTER :: section
1544 :
1545 : TYPE(keyword_type), POINTER :: keyword
1546 : TYPE(section_type), POINTER :: print_key, print_section, subsection
1547 :
1548 10278 : NULLIFY (keyword)
1549 10278 : CPASSERT(.NOT. ASSOCIATED(section))
1550 : CALL section_create(section, __LOCATION__, name="REAL_TIME_PROPAGATION", &
1551 : description="Parameters needed to set up the real time propagation"// &
1552 : " for the electron dynamics. This currently works only in the NVE ensemble.", &
1553 : n_keywords=4, n_subsections=4, repeats=.FALSE., &
1554 30834 : citations=[Kunert2003, Andermatt2016])
1555 :
1556 : CALL keyword_create(keyword, __LOCATION__, name="MAX_ITER", &
1557 : description="Maximal number of iterations for the self consistent propagator loop.", &
1558 : usage="MAX_ITER 10", &
1559 10278 : default_i_val=10)
1560 10278 : CALL section_add_keyword(section, keyword)
1561 10278 : CALL keyword_release(keyword)
1562 :
1563 : CALL keyword_create(keyword, __LOCATION__, name="EPS_ITER", &
1564 : description="Convergence criterion for the self consistent propagator loop.", &
1565 : usage="EPS_ITER 1.0E-5", &
1566 10278 : default_r_val=1.0E-7_dp)
1567 10278 : CALL section_add_keyword(section, keyword)
1568 10278 : CALL keyword_release(keyword)
1569 :
1570 : CALL keyword_create(keyword, __LOCATION__, name="ASPC_ORDER", &
1571 : description="Speciefies how many steps will be used for extrapolation. "// &
1572 : "One will be always used which is means X(t+dt)=X(t)", &
1573 : usage="ASPC_ORDER 3", &
1574 10278 : default_i_val=3)
1575 10278 : CALL section_add_keyword(section, keyword)
1576 10278 : CALL keyword_release(keyword)
1577 :
1578 : CALL keyword_create(keyword, __LOCATION__, name="MAT_EXP", &
1579 : description="Which method should be used to calculate the exponential"// &
1580 : " in the propagator. It is recommended to use BCH when employing density_propagation "// &
1581 : "and ARNOLDI otherwise.", &
1582 : usage="MAT_EXP TAYLOR", default_i_val=do_arnoldi, &
1583 : enum_c_vals=s2a("TAYLOR", "PADE", "ARNOLDI", "BCH", "EXACT"), &
1584 : enum_i_vals=[do_taylor, do_pade, do_arnoldi, do_bch, do_exact], &
1585 : enum_desc=s2a("exponential is evaluated using scaling and squaring in combination"// &
1586 : " with a taylor expansion of the exponential.", &
1587 : "uses scaling and squaring together with the pade approximation", &
1588 : "uses arnoldi subspace algorithm to compute exp(H)*MO directly, can't be used in "// &
1589 : "combination with Crank Nicholson or density propagation", &
1590 : "Uses a Baker-Campbell-Hausdorff expansion to propagate the density matrix,"// &
1591 : " only works for density propagation", &
1592 : "Uses diagonalisation of the exponent matrices to determine the "// &
1593 10278 : "matrix exponential exactly. Only implemented for GWBSE."))
1594 10278 : CALL section_add_keyword(section, keyword)
1595 10278 : CALL keyword_release(keyword)
1596 :
1597 : CALL keyword_create(keyword, __LOCATION__, name="DENSITY_PROPAGATION", &
1598 : description="The density matrix is propagated instead of the molecular orbitals. "// &
1599 : "This can allow a linear scaling simulation. The density matrix is filtered with "// &
1600 : "the threshold based on the EPS_FILTER keyword from the LS_SCF section", &
1601 : usage="DENSITY_PROPAGATION .TRUE.", &
1602 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1603 10278 : CALL section_add_keyword(section, keyword)
1604 10278 : CALL keyword_release(keyword)
1605 :
1606 : CALL keyword_create(keyword, __LOCATION__, name="SC_CHECK_START", &
1607 : description="Speciefies how many iteration steps will be done without "// &
1608 : "a check for self consistency. Can save some time in big calculations.", &
1609 : usage="SC_CHECK_START 3", &
1610 10278 : default_i_val=0)
1611 10278 : CALL section_add_keyword(section, keyword)
1612 10278 : CALL keyword_release(keyword)
1613 :
1614 : CALL keyword_create(keyword, __LOCATION__, name="EXP_ACCURACY", &
1615 : description="Accuracy for the taylor and pade approximation. "// &
1616 : "This is only an upper bound bound since the norm used for the guess "// &
1617 : "is an upper bound for the needed one.", &
1618 : usage="EXP_ACCURACY 1.0E-6", &
1619 10278 : default_r_val=1.0E-9_dp)
1620 10278 : CALL section_add_keyword(section, keyword)
1621 10278 : CALL keyword_release(keyword)
1622 :
1623 : CALL keyword_create(keyword, __LOCATION__, name="PROPAGATOR", &
1624 : description="Which propagator should be used for the orbitals", &
1625 : usage="PROPAGATOR ETRS", default_i_val=do_etrs, &
1626 : enum_c_vals=s2a("ETRS", "CN", "EM"), &
1627 : enum_i_vals=[do_etrs, do_cn, do_em], &
1628 : enum_desc=s2a("enforced time reversible symmetry", &
1629 : "Crank Nicholson propagator", &
1630 10278 : "Exponential midpoint propagator"))
1631 10278 : CALL section_add_keyword(section, keyword)
1632 10278 : CALL keyword_release(keyword)
1633 :
1634 : CALL keyword_create(keyword, __LOCATION__, name="INITIAL_WFN", &
1635 : description="Controls the initial WFN used for propagation. "// &
1636 : "Note that some energy contributions may not be "// &
1637 : "initialized in the restart cases, for instance "// &
1638 : "electronic entropy energy in the case of smearing.", &
1639 : usage="INITIAL_WFN SCF_WFN", default_i_val=use_scf_wfn, &
1640 : enum_c_vals=s2a("SCF_WFN", "RESTART_WFN", "RT_RESTART"), &
1641 : enum_i_vals=[use_scf_wfn, use_restart_wfn, use_rt_restart], &
1642 : enum_desc=s2a("An SCF run is performed to get the initial state.", &
1643 : "A wavefunction from a previous SCF is propagated. Especially useful,"// &
1644 : " if electronic constraints or restraints are used in the previous calculation, "// &
1645 : "since these do not work in the rtp scheme.", &
1646 10278 : "use the wavefunction of a real time propagation/ehrenfest run"))
1647 10278 : CALL section_add_keyword(section, keyword)
1648 10278 : CALL keyword_release(keyword)
1649 :
1650 : CALL keyword_create(keyword, __LOCATION__, name="APPLY_WFN_MIX_INIT_RESTART", &
1651 : description="If set to True and in the case of INITIAL_WFN=RESTART_WFN, call the "// &
1652 : "DFT%PRINT%WFN_MIX section to mix the read initial wfn. The starting wave-function of the "// &
1653 : "RTP will be the mixed one. Setting this to True without a defined WFN_MIX section will "// &
1654 : "not do anything as defining a WFN_MIX section without this keyword for RTP run with "// &
1655 : "INITIAL_WFN=RESTART_WFN. Note that if INITIAL_WFN=SCF_WFN, this keyword is not needed to "// &
1656 : "apply the mixing defined in the WFN_MIX section. Default is False.", &
1657 : usage="APPLY_WFN_MIX_INIT_RESTART", &
1658 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1659 10278 : CALL section_add_keyword(section, keyword)
1660 10278 : CALL keyword_release(keyword)
1661 :
1662 : CALL keyword_create(keyword, __LOCATION__, name="APPLY_DELTA_PULSE", &
1663 : description="Applies a delta kick to the initial wfn (only RTP for now - the EMD"// &
1664 : " case is not yet implemented). Only work for INITIAL_WFN=SCF_WFN", &
1665 : usage="APPLY_DELTA_PULSE", &
1666 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1667 10278 : CALL section_add_keyword(section, keyword)
1668 10278 : CALL keyword_release(keyword)
1669 :
1670 : CALL keyword_create(keyword, __LOCATION__, name="APPLY_DELTA_PULSE_MAG", &
1671 : description="Applies a magnetic delta kick to the initial wfn (only RTP for now - the EMD"// &
1672 : " case is not yet implemented). Only work for INITIAL_WFN=SCF_WFN", &
1673 : usage="APPLY_DELTA_PULSE_MAG", &
1674 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1675 10278 : CALL section_add_keyword(section, keyword)
1676 10278 : CALL keyword_release(keyword)
1677 :
1678 : CALL keyword_create(keyword, __LOCATION__, name="VELOCITY_GAUGE", &
1679 : description="Perform propagation in the velocity gauge using the explicit vector potential"// &
1680 : " only a constant vector potential as of now (corresonding to a delta-pulse)."// &
1681 : " uses DELTA_PULSE_SCALE and DELTA_PULSE_DIRECTION to define the vector potential", &
1682 : usage="VELOCITY_GAUGE T", &
1683 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1684 10278 : CALL section_add_keyword(section, keyword)
1685 10278 : CALL keyword_release(keyword)
1686 :
1687 : CALL keyword_create(keyword, __LOCATION__, name="GAUGE_ORIG", &
1688 : description="Define gauge origin for magnetic perturbation", &
1689 : usage="GAUGE_ORIG COM", &
1690 : enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
1691 : enum_desc=s2a("Use Center of Mass", &
1692 : "Use Center of Atomic Charges", &
1693 : "Use User Defined Point (Keyword:REF_POINT)", &
1694 : "Use Origin of Coordinate System"), &
1695 : enum_i_vals=[use_mom_ref_com, &
1696 : use_mom_ref_coac, &
1697 : use_mom_ref_user, &
1698 : use_mom_ref_zero], &
1699 10278 : default_i_val=use_mom_ref_com)
1700 10278 : CALL section_add_keyword(section, keyword)
1701 10278 : CALL keyword_release(keyword)
1702 :
1703 : CALL keyword_create(keyword, __LOCATION__, name="GAUGE_ORIG_MANUAL", &
1704 : description="Manually defined gauge origin for magnetic perturbation [in Bohr!]", &
1705 : usage="GAUGE_ORIG_MANUAL x y z", &
1706 : repeats=.FALSE., &
1707 : n_var=3, default_r_vals=[0._dp, 0._dp, 0._dp], &
1708 : type_of_var=real_t, &
1709 10278 : unit_str='bohr')
1710 10278 : CALL section_add_keyword(section, keyword)
1711 10278 : CALL keyword_release(keyword)
1712 :
1713 : CALL keyword_create(keyword, __LOCATION__, name="VG_COM_NL", &
1714 : description="apply gauge transformed non-local potential term"// &
1715 : " only affects VELOCITY_GAUGE=.TRUE.", &
1716 : usage="VG_COM_NL T", &
1717 10278 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1718 10278 : CALL section_add_keyword(section, keyword)
1719 10278 : CALL keyword_release(keyword)
1720 :
1721 : CALL keyword_create(keyword, __LOCATION__, name="COM_NL", &
1722 : description="Include non-local commutator for periodic delta pulse."// &
1723 : " only affects PERIODIC=.TRUE.", &
1724 : usage="COM_NL", &
1725 10278 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1726 10278 : CALL section_add_keyword(section, keyword)
1727 10278 : CALL keyword_release(keyword)
1728 :
1729 : CALL keyword_create(keyword, __LOCATION__, name="LEN_REP", &
1730 : description="Use length representation delta pulse (in conjunction with PERIODIC T)."// &
1731 : " This corresponds to a 1st order perturbation in the length gauge."// &
1732 : " Note that this is NOT compatible with a periodic calculation!"// &
1733 : " Uses the reference point defined in DFT%PRINT%MOMENTS ", &
1734 : usage="LEN_REP T", &
1735 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1736 10278 : CALL section_add_keyword(section, keyword)
1737 10278 : CALL keyword_release(keyword)
1738 :
1739 : CALL keyword_create(keyword, __LOCATION__, name="PERIODIC", &
1740 : description="Apply a delta-kick that is compatible with periodic boundary conditions"// &
1741 : " for any value of DELTA_PULSE_SCALE. Uses perturbation theory for the preparation of"// &
1742 : " the initial wfn with the velocity operator as perturbation."// &
1743 : " If LEN_REP is .FALSE. this corresponds to a first order velocity gauge."// &
1744 : " Note that the pulse is only applied when INITIAL_WFN is set to SCF_WFN,"// &
1745 : " and not for restarts (RT_RESTART).", &
1746 : usage="PERIODIC", &
1747 10278 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1748 10278 : CALL section_add_keyword(section, keyword)
1749 10278 : CALL keyword_release(keyword)
1750 :
1751 : CALL keyword_create(keyword, __LOCATION__, name="LOCALIZE", &
1752 : description="Localise the Molecular orbitals each n steps "// &
1753 : "real-time propagated TDDFT, 0 means never localise", &
1754 10278 : usage="LOCALIZE", default_i_val=0)
1755 10278 : CALL section_add_keyword(section, keyword)
1756 10278 : CALL keyword_release(keyword)
1757 :
1758 : CALL keyword_create(keyword, __LOCATION__, name="DELTA_PULSE_DIRECTION", &
1759 : description="Direction of the applied electric field. The k vector is given as"// &
1760 : " 2*Pi*[i,j,k]*inv(h_mat), which for PERIODIC .FALSE. yields exp(ikr) periodic with"// &
1761 : " the unit cell, only if DELTA_PULSE_SCALE is set to unity. For an orthorhombic cell"// &
1762 : " [1,0,0] yields [2*Pi/L_x,0,0]. For small cells, this results in a very large kick.", &
1763 : usage="DELTA_PULSE_DIRECTION 1 1 1", n_var=3, default_i_vals=[1, 0, 0], &
1764 10278 : type_of_var=integer_t)
1765 10278 : CALL section_add_keyword(section, keyword)
1766 10278 : CALL keyword_release(keyword)
1767 :
1768 : CALL keyword_create(keyword, __LOCATION__, name="DELTA_PULSE_SCALE", &
1769 : description="Scale the k vector, which for PERIODIC .FALSE. results in exp(ikr) no"// &
1770 : " longer being periodic with the unit cell. The norm of k is the strength of the"// &
1771 : " applied electric field in atomic units.", &
1772 10278 : usage="DELTA_PULSE_SCALE 0.01 ", n_var=1, default_r_val=0.001_dp)
1773 10278 : CALL section_add_keyword(section, keyword)
1774 10278 : CALL keyword_release(keyword)
1775 :
1776 : CALL keyword_create(keyword, __LOCATION__, name="HFX_BALANCE_IN_CORE", &
1777 : description="If HFX is used, this keyword forces a redistribution/recalculation"// &
1778 : " of the integrals, balanced with respect to the in core steps.", &
1779 : usage="HFX_BALANCE_IN_CORE", &
1780 10278 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1781 10278 : CALL section_add_keyword(section, keyword)
1782 10278 : CALL keyword_release(keyword)
1783 :
1784 : CALL keyword_create(keyword, __LOCATION__, name="MCWEENY_MAX_ITER", &
1785 : description="Determines the maximum amount of McWeeny steps used after each converged"// &
1786 : " step in density propagation", &
1787 10278 : usage="MCWEENY_MAX_ITER 2", default_i_val=1)
1788 10278 : CALL section_add_keyword(section, keyword)
1789 10278 : CALL keyword_release(keyword)
1790 :
1791 : CALL keyword_create( &
1792 : keyword, __LOCATION__, name="ACCURACY_REFINEMENT", &
1793 : description="If using density propagation some parts should be calculated with a higher accuracy than the rest"// &
1794 : " to reduce numerical noise. This factor determines by how much the filtering threshold is"// &
1795 : " reduced for these calculations.", &
1796 10278 : usage="ACCURACY_REFINEMENT", default_i_val=100)
1797 10278 : CALL section_add_keyword(section, keyword)
1798 10278 : CALL keyword_release(keyword)
1799 :
1800 : CALL keyword_create(keyword, __LOCATION__, name="MCWEENY_EPS", &
1801 : description="Threshold after which McWeeny is terminated", &
1802 : usage="MCWEENY_EPS 0.00001", &
1803 10278 : default_r_val=0.0_dp)
1804 10278 : CALL section_add_keyword(section, keyword)
1805 10278 : CALL keyword_release(keyword)
1806 :
1807 10278 : NULLIFY (print_section)
1808 : CALL section_create(print_section, __LOCATION__, name="PRINT", &
1809 : description="Section of possible print options for an RTP runs", &
1810 10278 : repeats=.FALSE.)
1811 :
1812 10278 : NULLIFY (print_key)
1813 : CALL cp_print_key_section_create(print_key, __LOCATION__, "PROGRAM_RUN_INFO", &
1814 : description="Controls the printing within real time propagation and Eherenfest dynamics", &
1815 10278 : print_level=low_print_level, filename="__STD_OUT__")
1816 10278 : CALL section_add_subsection(print_section, print_key)
1817 10278 : CALL section_release(print_key)
1818 :
1819 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RESTART", &
1820 : description="Controls the dumping of the MO restart file during rtp. "// &
1821 : "By default keeps a short history of three restarts. "// &
1822 : "See also RESTART_HISTORY. In density propagation this controls the printing of "// &
1823 : "density matrix.", &
1824 : print_level=low_print_level, common_iter_levels=3, &
1825 : each_iter_names=s2a("MD"), each_iter_values=[20], &
1826 10278 : add_last=add_last_numeric, filename="RESTART")
1827 : CALL keyword_create(keyword, __LOCATION__, name="BACKUP_COPIES", &
1828 : description="Specifies the maximum number of backup copies.", &
1829 : usage="BACKUP_COPIES {int}", &
1830 10278 : default_i_val=1)
1831 10278 : CALL section_add_keyword(print_key, keyword)
1832 10278 : CALL keyword_release(keyword)
1833 10278 : CALL section_add_subsection(print_section, print_key)
1834 10278 : CALL section_release(print_key)
1835 :
1836 : CALL cp_print_key_section_create(print_key, __LOCATION__, "RESTART_HISTORY", &
1837 : description="Dumps unique MO restart files during the run keeping all of them. "// &
1838 : "In density propagation it dumps the density matrix instead", &
1839 : print_level=low_print_level, common_iter_levels=0, &
1840 : each_iter_names=s2a("MD"), &
1841 : each_iter_values=[500], &
1842 10278 : filename="RESTART")
1843 : CALL keyword_create(keyword, __LOCATION__, name="BACKUP_COPIES", &
1844 : description="Specifies the maximum number of backup copies.", &
1845 : usage="BACKUP_COPIES {int}", &
1846 10278 : default_i_val=1)
1847 10278 : CALL section_add_keyword(print_key, keyword)
1848 10278 : CALL keyword_release(keyword)
1849 10278 : CALL section_add_subsection(print_section, print_key)
1850 10278 : CALL section_release(print_key)
1851 :
1852 : CALL cp_print_key_section_create(print_key, __LOCATION__, "FIELD", &
1853 : description="Print the time-dependent field applied during an EMD simulation in "// &
1854 : "atomic unit.", &
1855 : print_level=high_print_level, common_iter_levels=1, &
1856 : each_iter_names=s2a("MD"), &
1857 : each_iter_values=[1], &
1858 10278 : filename="FIELD")
1859 10278 : CALL section_add_subsection(print_section, print_key)
1860 10278 : CALL section_release(print_key)
1861 :
1862 10278 : CALL create_projection_rtp_section(print_key)
1863 10278 : CALL section_add_subsection(print_section, print_key)
1864 10278 : CALL section_release(print_key)
1865 :
1866 : CALL cp_print_key_section_create(print_key, __LOCATION__, "CURRENT_INT", &
1867 : description="Print the integral of the current density (only if the"// &
1868 : " imaginary part of the density is NOT zero.", &
1869 : print_level=high_print_level, common_iter_levels=1, &
1870 : each_iter_names=s2a("MD"), &
1871 : each_iter_values=[1], &
1872 10278 : filename="rtp_j_int")
1873 10278 : CALL section_add_subsection(print_section, print_key)
1874 10278 : CALL section_release(print_key)
1875 :
1876 : CALL cp_print_key_section_create(print_key, __LOCATION__, "CURRENT", &
1877 : description="Print the current during an EMD simulation to cube files.", &
1878 : print_level=high_print_level, common_iter_levels=0, &
1879 : each_iter_names=s2a("MD"), &
1880 : each_iter_values=[20], &
1881 10278 : filename="current")
1882 : CALL keyword_create(keyword, __LOCATION__, name="BACKUP_COPIES", &
1883 : description="Specifies the maximum number of backup copies.", &
1884 : usage="BACKUP_COPIES {int}", &
1885 10278 : default_i_val=1)
1886 10278 : CALL section_add_keyword(print_key, keyword)
1887 10278 : CALL keyword_release(keyword)
1888 : CALL keyword_create(keyword, __LOCATION__, name="STRIDE", &
1889 : description="The stride (X,Y,Z) used to write the cube file "// &
1890 : "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1891 : " 1 number valid for all components.", &
1892 10278 : usage="STRIDE 2 2 2", n_var=-1, default_i_vals=[2, 2, 2], type_of_var=integer_t)
1893 10278 : CALL section_add_keyword(print_key, keyword)
1894 10278 : CALL keyword_release(keyword)
1895 :
1896 10278 : CALL section_add_subsection(print_section, print_key)
1897 10278 : CALL section_release(print_key)
1898 :
1899 : ! Marek : Add print option for ASCII density files - DEVELPMENT ONLY?
1900 : CALL cp_print_key_section_create(print_key, __LOCATION__, "DENSITY_MATRIX", &
1901 : description="Prints the density matrix at iterations in clear text to a file", &
1902 : print_level=high_print_level, common_iter_levels=0, &
1903 : each_iter_names=s2a("MD"), &
1904 : each_iter_values=[1], &
1905 10278 : filename="rho")
1906 10278 : CALL section_add_subsection(print_section, print_key)
1907 10278 : CALL section_release(print_key)
1908 : ! Marek : Moments ASCII print
1909 : CALL cp_print_key_section_create(print_key, __LOCATION__, "MOMENTS", &
1910 : description="Prints the time-dependent electronic moments at "// &
1911 : "iterations in clear text to a file.", &
1912 : print_level=high_print_level, common_iter_levels=1, &
1913 : each_iter_names=s2a("MD"), &
1914 : each_iter_values=[1], &
1915 10278 : filename="__STD_OUT__")
1916 : CALL keyword_create(keyword, __LOCATION__, name="REFERENCE", &
1917 : variants=s2a("REF"), &
1918 : description="Define the reference point for the calculation of the electrostatic moment.", &
1919 : usage="REFERENCE COM", &
1920 : enum_c_vals=s2a("COM", "COAC", "USER_DEFINED", "ZERO"), &
1921 : enum_desc=s2a("Use Center of Mass", &
1922 : "Use Center of Atomic Charges", &
1923 : "Use User Defined Point (Keyword:REFERENCE_POINT)", &
1924 : "Use Origin of Coordinate System"), &
1925 : enum_i_vals=[use_mom_ref_com, &
1926 : use_mom_ref_coac, &
1927 : use_mom_ref_user, &
1928 : use_mom_ref_zero], &
1929 10278 : default_i_val=use_mom_ref_coac)
1930 10278 : CALL section_add_keyword(print_key, keyword)
1931 10278 : CALL keyword_release(keyword)
1932 :
1933 : CALL keyword_create(keyword, __LOCATION__, name="REFERENCE_POINT", &
1934 : variants=s2a("REF_POINT"), &
1935 : description="Fixed reference point for the calculations of the electrostatic moment.", &
1936 : usage="REFERENCE_POINT x y z", &
1937 : repeats=.FALSE., &
1938 : n_var=3, default_r_vals=[0._dp, 0._dp, 0._dp], &
1939 : type_of_var=real_t, &
1940 10278 : unit_str='angstrom')
1941 10278 : CALL section_add_keyword(print_key, keyword)
1942 10278 : CALL keyword_release(keyword)
1943 10278 : CALL section_add_subsection(print_section, print_key)
1944 10278 : CALL section_release(print_key)
1945 : ! Marek : Fourier transform of MOMENTS ASCII print
1946 : CALL cp_print_key_section_create(print_key, __LOCATION__, "MOMENTS_FT", &
1947 : description="Prints the calculated Fourier transform of "// &
1948 : "time-dependent moments. For calculations with real time pulse (not delta kick) "// &
1949 : "can be supplied with starting time.", &
1950 : print_level=medium_print_level, common_iter_levels=0, &
1951 : each_iter_names=s2a("MD"), &
1952 : each_iter_values=[1], &
1953 10278 : filename="MOMENTS_FT")
1954 10278 : CALL section_add_subsection(print_section, print_key)
1955 10278 : CALL section_release(print_key)
1956 : ! Marek : Chosen element of (Fourier transformed) polarizability tensor (energy dependent) - text format
1957 : CALL cp_print_key_section_create(print_key, __LOCATION__, "POLARIZABILITY", &
1958 : description="Prints the chosen element of the energy dependent polarizability tensor "// &
1959 : "to a specified file. The tensor is calculated as ratio of "// &
1960 : "Fourier transform of the dipole "// &
1961 : "moment trace and Fourier transform of the applied field "// &
1962 : "(for delta kick, constant real field is applied.", &
1963 : print_level=medium_print_level, common_iter_levels=0, &
1964 : each_iter_names=s2a("MD"), &
1965 : each_iter_values=[1], &
1966 10278 : filename="POLARIZABILITY")
1967 : CALL keyword_create(keyword, __LOCATION__, "ELEMENT", &
1968 : description="Specifies the element of polarizability which is to be printed out "// &
1969 : "(indexing starts at 1). If not explicitly provided, RTBSE code tries to guess "// &
1970 : "the optimal values - for applied electric field (both delta pulse and RT field) "// &
1971 : "with only a single non-zero cartesian component, prints the 3 trivially available elements.", &
1972 10278 : type_of_var=integer_t, default_i_vals=[1, 1], n_var=2, usage="ELEMENT 1 1", repeats=.TRUE.)
1973 10278 : CALL section_add_keyword(print_key, keyword)
1974 10278 : CALL keyword_release(keyword)
1975 10278 : CALL section_add_subsection(print_section, print_key)
1976 10278 : CALL section_release(print_key)
1977 :
1978 : CALL cp_print_key_section_create(print_key, __LOCATION__, "E_CONSTITUENTS", &
1979 : description="Print the energy constituents (relevant to RTP) which make up "// &
1980 : "the Total Energy", &
1981 : print_level=high_print_level, common_iter_levels=1, &
1982 : each_iter_names=s2a("MD"), &
1983 : each_iter_values=[1], &
1984 10278 : filename="rtp")
1985 10278 : CALL section_add_subsection(print_section, print_key)
1986 10278 : CALL section_release(print_key)
1987 :
1988 10278 : CALL section_add_subsection(section, print_section)
1989 10278 : CALL section_release(print_section)
1990 :
1991 : ! RTBSE subsection
1992 10278 : NULLIFY (subsection)
1993 10278 : CALL create_rtbse_section(subsection)
1994 10278 : CALL section_add_subsection(section, subsection)
1995 10278 : CALL section_release(subsection)
1996 : ! FT subsection
1997 10278 : CALL create_ft_section(subsection)
1998 10278 : CALL section_add_subsection(section, subsection)
1999 10278 : CALL section_release(subsection)
2000 :
2001 10278 : END SUBROUTINE create_rtp_section
2002 : ! **************************************************************************************************
2003 : !> \brief Creates the subsection for specialized options of RTBSE code
2004 : !> \param section The created RTBSE section
2005 : !> \author Stepan Marek
2006 : ! **************************************************************************************************
2007 10278 : SUBROUTINE create_rtbse_section(section)
2008 : TYPE(section_type), POINTER :: section
2009 :
2010 : TYPE(keyword_type), POINTER :: keyword
2011 :
2012 10278 : NULLIFY (keyword)
2013 10278 : CPASSERT(.NOT. ASSOCIATED(section))
2014 :
2015 : CALL section_create(section, __LOCATION__, name="RTBSE", &
2016 : description="Controls options for the real-time Bethe-Salpeter (RTBSE) propagation. "// &
2017 : "Note that running RTBSE requires previous low-scaling "// &
2018 : "[GW](#CP2K_INPUT.FORCE_EVAL.PROPERTIES.BANDSTRUCTURE.GW) calculation. Also note that "// &
2019 : "designating this section as RTBSE run but choosing run type ENERGY leads to potential "// &
2020 : "deallocation errors. More details (including description of output files) is available in "// &
2021 : "the [methods](../../../../methods/properties/optical/rtbse) section of the documentation.", &
2022 20556 : repeats=.FALSE., citations=[Marek2025])
2023 :
2024 : ! Marek : Controlling flow to RTBSE
2025 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
2026 : description="Which method is used for the time propagation of electronic structure. "// &
2027 : "By default, use the TDDFT method. Can also choose RT-BSE method, which propagates the lesser "// &
2028 : "Green's function instead of density matrix/molecular orbitals.", &
2029 : usage="&RTBSE TDDFT", &
2030 : default_i_val=rtp_method_tddft, &
2031 : lone_keyword_i_val=rtp_method_bse, &
2032 : enum_c_vals=s2a("TDDFT", "RTBSE"), &
2033 : enum_i_vals=[rtp_method_tddft, rtp_method_bse], &
2034 : enum_desc=s2a("Use TDDFT for density matrix/MO propagation.", &
2035 10278 : "Use RT-BSE for Green's function propagation"))
2036 10278 : CALL section_add_keyword(section, keyword)
2037 10278 : CALL keyword_release(keyword)
2038 :
2039 : ! Marek : Development option - run GWBSE starting from the KS Hamiltonian
2040 : CALL keyword_create(keyword, __LOCATION__, name="RTBSE_HAMILTONIAN", &
2041 : description="Which Hamiltonian to use as the single-particle Hamiltonian"// &
2042 : " in the Green's propagator.", &
2043 : usage="RTBSE_HAMILTONIAN G0W0", &
2044 : default_i_val=rtp_bse_ham_g0w0, &
2045 : enum_c_vals=s2a("KS", "G0W0"), &
2046 : enum_i_vals=[rtp_bse_ham_ks, rtp_bse_ham_g0w0], &
2047 : enum_desc=s2a("Use Kohn-Sham Hamiltonian for Green's propagation.", &
2048 10278 : "Use G0W0 Hamiltonian for Green's function propagation"))
2049 10278 : CALL section_add_keyword(section, keyword)
2050 10278 : CALL keyword_release(keyword)
2051 :
2052 10278 : END SUBROUTINE create_rtbse_section
2053 : ! **************************************************************************************************
2054 : !> \brief Creates the subsection for Fourier transform options applicable to RTP output
2055 : !> \param ft_section The created FT section
2056 : !> \date 11.2025
2057 : !> \author Stepan Marek
2058 : ! **************************************************************************************************
2059 10278 : SUBROUTINE create_ft_section(ft_section)
2060 : TYPE(section_type), POINTER :: ft_section
2061 :
2062 : TYPE(keyword_type), POINTER :: keyword
2063 : TYPE(section_type), POINTER :: subsection
2064 :
2065 10278 : CPASSERT(.NOT. ASSOCIATED(ft_section))
2066 :
2067 : ! Create the section itself
2068 : CALL section_create(ft_section, __LOCATION__, &
2069 : name="FT", &
2070 : description="Define parameters for Fourier transforms used in RTP outputs.", &
2071 10278 : repeats=.FALSE.)
2072 :
2073 : ! Start time keyword
2074 10278 : NULLIFY (keyword)
2075 : CALL keyword_create(keyword, __LOCATION__, "START_TIME", &
2076 : description="The starting time from which damping is applied and from which on the trace is "// &
2077 : "considered for the Fourier transform (Fourier transform is used for the calculation of "// &
2078 : "MOMENTS_FT and POLARIZABILITY). Useful for real-time pulse - "// &
2079 : "one can specify the center of the pulse as the starting point.", &
2080 : type_of_var=real_t, &
2081 : unit_str="fs", &
2082 10278 : default_r_val=0.0_dp)
2083 10278 : CALL section_add_keyword(ft_section, keyword)
2084 10278 : CALL keyword_release(keyword)
2085 :
2086 : ! Damping keyword
2087 : CALL keyword_create(keyword, __LOCATION__, "DAMPING", &
2088 : description="Numerical Fourier transform (required for calculation of "// &
2089 : "MOMENTS_FT and POLARIZABILITY) can oscillate "// &
2090 : "when the final time trace values are far away from zero. "// &
2091 : "This keyword controls the exponential damping added to the Fourier transform "// &
2092 : "(Fourier transform is used for calculation of MOMENTS_FT and POLARIZABILITY). "// &
2093 : "For negative values (the default), calculates the damping at the run time so that the last point "// &
2094 : "in the time trace is reduced by factor e^(-4). When set manually, determines the time in which "// &
2095 : "the moments trace is reduced by factor of e^(-1), except when set to zero, in which case "// &
2096 : "the damping is not applied.", &
2097 : type_of_var=real_t, &
2098 : unit_str="fs", &
2099 10278 : default_r_val=-1.0_dp/femtoseconds)
2100 10278 : CALL section_add_keyword(ft_section, keyword)
2101 10278 : CALL keyword_release(keyword)
2102 :
2103 : ! Create the Padé subsection
2104 10278 : NULLIFY (subsection)
2105 : CALL section_create(subsection, __LOCATION__, name="PADE", &
2106 : description="Defines the parameters for the Padé interpolation of the "// &
2107 : "Fourier transforms used in the output of RTP. Only available with the GreenX library linked to CP2K.", &
2108 10278 : repeats=.FALSE.)
2109 :
2110 : ! Explicit presence of the section turns on the Padé interpolation
2111 : CALL keyword_create(keyword, __LOCATION__, "_SECTION_PARAMETERS_", &
2112 : description="Turns on the Padé interpolation", &
2113 : type_of_var=logical_t, &
2114 : default_l_val=.FALSE., &
2115 10278 : lone_keyword_l_val=.TRUE.)
2116 10278 : CALL section_add_keyword(subsection, keyword)
2117 10278 : CALL keyword_release(keyword)
2118 :
2119 : ! Minimum interpolated energy
2120 : CALL keyword_create(keyword, __LOCATION__, "E_MIN", &
2121 : description="The minimum energy of the Padé interpolation output.", &
2122 : type_of_var=real_t, &
2123 : unit_str="eV", &
2124 10278 : default_r_val=0.0_dp)
2125 10278 : CALL section_add_keyword(subsection, keyword)
2126 10278 : CALL keyword_release(keyword)
2127 :
2128 : ! Maximum interpolated energy
2129 : CALL keyword_create(keyword, __LOCATION__, "E_MAX", &
2130 : description="The maximum energy of the Padé interpolation output.", &
2131 : type_of_var=real_t, &
2132 : unit_str="eV", &
2133 10278 : default_r_val=100.0_dp)
2134 10278 : CALL section_add_keyword(subsection, keyword)
2135 10278 : CALL keyword_release(keyword)
2136 :
2137 : ! Energy resolution
2138 : CALL keyword_create(keyword, __LOCATION__, "E_STEP", &
2139 : description="The energy resolution of the Padé interpolation output.", &
2140 : type_of_var=real_t, &
2141 : unit_str="eV", &
2142 10278 : default_r_val=0.02_dp/evolt)
2143 10278 : CALL section_add_keyword(subsection, keyword)
2144 10278 : CALL keyword_release(keyword)
2145 :
2146 : ! Minimum fitting energy
2147 : CALL keyword_create(keyword, __LOCATION__, "FIT_E_MIN", &
2148 : description="The lower boundary in energy for the points "// &
2149 : "used in the fitting of Padé parameters. If negative, uses "// &
2150 : "value of E_MIN (default).", &
2151 : type_of_var=real_t, &
2152 : unit_str="eV", &
2153 10278 : default_r_val=-1.0_dp)
2154 10278 : CALL section_add_keyword(subsection, keyword)
2155 10278 : CALL keyword_release(keyword)
2156 :
2157 : ! Maximum fitting energy
2158 : CALL keyword_create(keyword, __LOCATION__, "FIT_E_MAX", &
2159 : description="The upper boundary in energy for the points "// &
2160 : "used in the fitting of Padé parameters. If negative, uses "// &
2161 : "the value of E_MAX (default).", &
2162 : type_of_var=real_t, &
2163 : unit_str="eV", &
2164 10278 : default_r_val=-1.0_dp)
2165 10278 : CALL section_add_keyword(subsection, keyword)
2166 10278 : CALL keyword_release(keyword)
2167 :
2168 : ! Add the Padé subsection
2169 10278 : CALL section_add_subsection(ft_section, subsection)
2170 10278 : CALL section_release(subsection)
2171 10278 : END SUBROUTINE create_ft_section
2172 :
2173 : ! **************************************************************************************************
2174 : !> \brief Create CP2K input section for the SCCS model
2175 : !> \param section ...
2176 : !> \par History:
2177 : !> - Creation (10.10.2013,MK)
2178 : !> \author Matthias Krack (MK)
2179 : !> \version 1.0
2180 : ! **************************************************************************************************
2181 10278 : SUBROUTINE create_sccs_section(section)
2182 :
2183 : TYPE(section_type), POINTER :: section
2184 :
2185 : TYPE(keyword_type), POINTER :: keyword
2186 : TYPE(section_type), POINTER :: subsection
2187 :
2188 10278 : CPASSERT(.NOT. ASSOCIATED(section))
2189 :
2190 : CALL section_create(section, __LOCATION__, &
2191 : name="SCCS", &
2192 : description="Define the parameters for self-consistent continuum solvation (SCCS) model", &
2193 : citations=[Fattebert2002, Andreussi2012, Yin2017], &
2194 : n_keywords=8, &
2195 : n_subsections=2, &
2196 41112 : repeats=.FALSE.)
2197 :
2198 10278 : NULLIFY (keyword)
2199 :
2200 : CALL keyword_create(keyword, __LOCATION__, &
2201 : name="_SECTION_PARAMETERS_", &
2202 : description="Controls the activation of the SCCS section", &
2203 : usage="&SCCS ON", &
2204 : default_l_val=.FALSE., &
2205 10278 : lone_keyword_l_val=.TRUE.)
2206 10278 : CALL section_add_keyword(section, keyword)
2207 10278 : CALL keyword_release(keyword)
2208 :
2209 : CALL keyword_create(keyword, __LOCATION__, &
2210 : name="ALPHA", &
2211 : description="Solvent specific tunable parameter for the calculation of "// &
2212 : "the repulsion term $G^\text{rep} = \alpha S$ "// &
2213 : "where $S$ is the (quantum) surface of the cavity", &
2214 : repeats=.FALSE., &
2215 : n_var=1, &
2216 : type_of_var=real_t, &
2217 : default_r_val=0.0_dp, &
2218 10278 : unit_str="mN/m")
2219 10278 : CALL section_add_keyword(section, keyword)
2220 10278 : CALL keyword_release(keyword)
2221 :
2222 : CALL keyword_create(keyword, __LOCATION__, &
2223 : name="BETA", &
2224 : description="Solvent specific tunable parameter for the calculation of "// &
2225 : "the dispersion term $G^\text{dis} = \beta V$ "// &
2226 : "where $V$ is the (quantum) volume of the cavity", &
2227 : repeats=.FALSE., &
2228 : n_var=1, &
2229 : type_of_var=real_t, &
2230 : default_r_val=0.0_dp, &
2231 10278 : unit_str="GPa")
2232 10278 : CALL section_add_keyword(section, keyword)
2233 10278 : CALL keyword_release(keyword)
2234 :
2235 : CALL keyword_create(keyword, __LOCATION__, &
2236 : name="DELTA_RHO", &
2237 : description="Numerical increment for the calculation of the (quantum) "// &
2238 : "surface of the solute cavity", &
2239 : repeats=.FALSE., &
2240 : n_var=1, &
2241 : type_of_var=real_t, &
2242 10278 : default_r_val=2.0E-5_dp)
2243 10278 : CALL section_add_keyword(section, keyword)
2244 10278 : CALL keyword_release(keyword)
2245 :
2246 : CALL keyword_create(keyword, __LOCATION__, &
2247 : name="DERIVATIVE_METHOD", &
2248 : description="Method for the calculation of the numerical derivatives on the real-space grids", &
2249 : usage="DERIVATIVE_METHOD cd5", &
2250 : repeats=.FALSE., &
2251 : n_var=1, &
2252 : default_i_val=sccs_derivative_fft, &
2253 : enum_c_vals=s2a("FFT", "CD3", "CD5", "CD7"), &
2254 : enum_i_vals=[sccs_derivative_fft, &
2255 : sccs_derivative_cd3, &
2256 : sccs_derivative_cd5, &
2257 : sccs_derivative_cd7], &
2258 : enum_desc=s2a("Fast Fourier transformation", &
2259 : "3-point stencil central differences", &
2260 : "5-point stencil central differences", &
2261 10278 : "7-point stencil central differences"))
2262 10278 : CALL section_add_keyword(section, keyword)
2263 10278 : CALL keyword_release(keyword)
2264 :
2265 : CALL keyword_create(keyword, __LOCATION__, &
2266 : name="RELATIVE_PERMITTIVITY", &
2267 : variants=s2a("DIELECTRIC_CONSTANT", "EPSILON_RELATIVE", "EPSILON_SOLVENT"), &
2268 : description="Relative permittivity (dielectric constant) of the solvent (medium)", &
2269 : repeats=.FALSE., &
2270 : n_var=1, &
2271 : type_of_var=real_t, &
2272 : default_r_val=80.0_dp, &
2273 10278 : usage="RELATIVE_PERMITTIVITY 78.36")
2274 10278 : CALL section_add_keyword(section, keyword)
2275 10278 : CALL keyword_release(keyword)
2276 :
2277 : CALL keyword_create(keyword, __LOCATION__, &
2278 : name="EPS_SCCS", &
2279 : variants=s2a("EPS_ITER", "TAU_POL"), &
2280 : description="Tolerance for the convergence of the polarisation density, "// &
2281 : "i.e. requested accuracy for the SCCS iteration cycle", &
2282 : repeats=.FALSE., &
2283 : n_var=1, &
2284 : type_of_var=real_t, &
2285 : default_r_val=1.0E-6_dp, &
2286 10278 : usage="EPS_ITER 1.0E-7")
2287 10278 : CALL section_add_keyword(section, keyword)
2288 10278 : CALL keyword_release(keyword)
2289 :
2290 : CALL keyword_create(keyword, __LOCATION__, &
2291 : name="EPS_SCF", &
2292 : description="The SCCS iteration cycle is activated only if the SCF iteration cycle "// &
2293 : "is converged to this threshold value", &
2294 : repeats=.FALSE., &
2295 : n_var=1, &
2296 : type_of_var=real_t, &
2297 : default_r_val=0.5_dp, &
2298 10278 : usage="EPS_SCF 1.0E-2")
2299 10278 : CALL section_add_keyword(section, keyword)
2300 10278 : CALL keyword_release(keyword)
2301 :
2302 : CALL keyword_create(keyword, __LOCATION__, &
2303 : name="GAMMA", &
2304 : variants=s2a("SURFACE_TENSION"), &
2305 : description="Surface tension of the solvent used for the calculation of "// &
2306 : "the cavitation term $G^\text{cav} = \gamma S$ "// &
2307 : "where $S$ is the (quantum) surface of the cavity", &
2308 : repeats=.FALSE., &
2309 : n_var=1, &
2310 : type_of_var=real_t, &
2311 : default_r_val=0.0_dp, &
2312 10278 : unit_str="mN/m")
2313 10278 : CALL section_add_keyword(section, keyword)
2314 10278 : CALL keyword_release(keyword)
2315 :
2316 : CALL keyword_create(keyword, __LOCATION__, &
2317 : name="MAX_ITER", &
2318 : description="Maximum number of SCCS iteration steps performed to converge "// &
2319 : "within the given tolerance", &
2320 : repeats=.FALSE., &
2321 : n_var=1, &
2322 : type_of_var=integer_t, &
2323 : default_i_val=100, &
2324 10278 : usage="MAX_ITER 50")
2325 10278 : CALL section_add_keyword(section, keyword)
2326 10278 : CALL keyword_release(keyword)
2327 :
2328 : CALL keyword_create(keyword, __LOCATION__, &
2329 : name="METHOD", &
2330 : description="Method used for the smoothing of the dielectric function", &
2331 : usage="METHOD Fattebert-Gygi", &
2332 : default_i_val=sccs_andreussi, &
2333 : enum_c_vals=s2a("ANDREUSSI", "FATTEBERT-GYGI"), &
2334 : enum_i_vals=[sccs_andreussi, sccs_fattebert_gygi], &
2335 : enum_desc=s2a("Smoothing function proposed by Andreussi et al.", &
2336 10278 : "Smoothing function proposed by Fattebert and Gygi"))
2337 10278 : CALL section_add_keyword(section, keyword)
2338 10278 : CALL keyword_release(keyword)
2339 :
2340 : CALL keyword_create(keyword, __LOCATION__, &
2341 : name="MIXING", &
2342 : variants=["ETA"], &
2343 : description="Mixing parameter (Hartree damping) employed during the iteration procedure", &
2344 : repeats=.FALSE., &
2345 : n_var=1, &
2346 : type_of_var=real_t, &
2347 : default_r_val=0.6_dp, &
2348 20556 : usage="MIXING 0.2")
2349 10278 : CALL section_add_keyword(section, keyword)
2350 10278 : CALL keyword_release(keyword)
2351 :
2352 10278 : NULLIFY (subsection)
2353 :
2354 : CALL section_create(subsection, __LOCATION__, &
2355 : name="ANDREUSSI", &
2356 : description="Define the parameters of the dielectric smoothing function proposed by "// &
2357 : "Andreussi et al.", &
2358 : citations=[Andreussi2012], &
2359 : n_keywords=2, &
2360 : n_subsections=0, &
2361 20556 : repeats=.FALSE.)
2362 :
2363 : CALL keyword_create(keyword, __LOCATION__, &
2364 : name="RHO_MAX", &
2365 : description="Maximum density value used for the smoothing of the dielectric function", &
2366 : repeats=.FALSE., &
2367 : n_var=1, &
2368 : type_of_var=real_t, &
2369 : default_r_val=0.0035_dp, &
2370 10278 : usage="RHO_MAX 0.01")
2371 10278 : CALL section_add_keyword(subsection, keyword)
2372 10278 : CALL keyword_release(keyword)
2373 :
2374 : CALL keyword_create(keyword, __LOCATION__, &
2375 : name="RHO_MIN", &
2376 : description="Minimum density value used for the smoothing of the dielectric function", &
2377 : repeats=.FALSE., &
2378 : n_var=1, &
2379 : type_of_var=real_t, &
2380 : default_r_val=0.0001_dp, &
2381 10278 : usage="RHO_MIN 0.0003")
2382 10278 : CALL section_add_keyword(subsection, keyword)
2383 10278 : CALL keyword_release(keyword)
2384 :
2385 10278 : CALL section_add_subsection(section, subsection)
2386 10278 : CALL section_release(subsection)
2387 :
2388 : CALL section_create(subsection, __LOCATION__, &
2389 : name="FATTEBERT-GYGI", &
2390 : description="Define the parameters of the dielectric smoothing function proposed by "// &
2391 : "Fattebert and Gygi", &
2392 : citations=[Fattebert2002], &
2393 : n_keywords=2, &
2394 : n_subsections=0, &
2395 20556 : repeats=.FALSE.)
2396 :
2397 : CALL keyword_create(keyword, __LOCATION__, &
2398 : name="BETA", &
2399 : description="Parameter β changes the width of the interface solute-solvent", &
2400 : repeats=.FALSE., &
2401 : n_var=1, &
2402 : type_of_var=real_t, &
2403 : default_r_val=1.7_dp, &
2404 10278 : usage="BETA 1.3")
2405 10278 : CALL section_add_keyword(subsection, keyword)
2406 10278 : CALL keyword_release(keyword)
2407 :
2408 : CALL keyword_create(keyword, __LOCATION__, &
2409 : name="RHO_ZERO", &
2410 : variants=["RHO0"], &
2411 : description="Parameter $\rho_0$ defines the critical density in the middle "// &
2412 : "of the interface solute-solvent", &
2413 : repeats=.FALSE., &
2414 : n_var=1, &
2415 : type_of_var=real_t, &
2416 : default_r_val=0.0006_dp, &
2417 20556 : usage="RHO_ZERO 0.0004")
2418 10278 : CALL section_add_keyword(subsection, keyword)
2419 10278 : CALL keyword_release(keyword)
2420 :
2421 10278 : CALL section_add_subsection(section, subsection)
2422 10278 : CALL section_release(subsection)
2423 :
2424 10278 : END SUBROUTINE create_sccs_section
2425 :
2426 : ! **************************************************************************************************
2427 : !> \brief Create CP2K input section for the planar counter charge density
2428 : !> \param section ...
2429 : !> \author Ziwei Chai
2430 : !> \version 1.0
2431 : ! **************************************************************************************************
2432 10278 : SUBROUTINE create_pcc_section(section)
2433 :
2434 : TYPE(section_type), POINTER :: section
2435 :
2436 : TYPE(keyword_type), POINTER :: keyword
2437 :
2438 10278 : CPASSERT(.NOT. ASSOCIATED(section))
2439 :
2440 : CALL section_create(section, __LOCATION__, &
2441 : name="PLANAR_COUNTER_CHARGE", &
2442 : description="Define the parameters for the planar counter charge density", &
2443 : citations=[Chai2024a], &
2444 : n_keywords=4, &
2445 : n_subsections=0, &
2446 20556 : repeats=.FALSE.)
2447 :
2448 10278 : NULLIFY (keyword)
2449 :
2450 : CALL keyword_create(keyword, __LOCATION__, &
2451 : name="_SECTION_PARAMETERS_", &
2452 : description="Controls the activation of the planar counter charge section", &
2453 : usage="&PLANAR_COUNTER_CHARGE ON", &
2454 : default_l_val=.FALSE., &
2455 10278 : lone_keyword_l_val=.TRUE.)
2456 10278 : CALL section_add_keyword(section, keyword)
2457 10278 : CALL keyword_release(keyword)
2458 :
2459 : CALL keyword_create(keyword, __LOCATION__, name="PARALLEL_PLANE", &
2460 : enum_c_vals=s2a('XY', 'YZ', 'XZ'), &
2461 : enum_i_vals=[3, 1, 2], &
2462 : description="The coordinate plane that the surface is parallel to.", &
2463 : enum_desc=s2a("Parallel to XY", "Parallel to YZ", "Parallel to XZ"), &
2464 : n_var=1, &
2465 : default_i_val=3, &
2466 10278 : usage="PARALLEL_PLANE XY")
2467 10278 : CALL section_add_keyword(section, keyword)
2468 10278 : CALL keyword_release(keyword)
2469 :
2470 : CALL keyword_create(keyword, __LOCATION__, &
2471 : name="DIST_EDGE", &
2472 : description="Controls the distance between the center of the Gaussian "// &
2473 : "and the cell boundary", &
2474 : usage="DIST_EDGE 1.0", &
2475 : unit_str="angstrom", &
2476 10278 : type_of_var=real_t)
2477 10278 : CALL section_add_keyword(section, keyword)
2478 10278 : CALL keyword_release(keyword)
2479 :
2480 : CALL keyword_create(keyword, __LOCATION__, &
2481 : name="GAU_C", &
2482 : description="Controls the spread of the Gaussian distribution", &
2483 : usage="GAU_C 0.1", &
2484 : unit_str="angstrom", &
2485 10278 : type_of_var=real_t)
2486 10278 : CALL section_add_keyword(section, keyword)
2487 10278 : CALL keyword_release(keyword)
2488 :
2489 10278 : END SUBROUTINE create_pcc_section
2490 :
2491 : ! **************************************************************************************************
2492 : !> \brief Create CP2K input section for calculating and printing the planar averaged
2493 : !> electrostatic potential (Hartree potential) for symmetric slab systems
2494 : !> \param section ...
2495 : !> \author Ziwei Chai
2496 : !> \version 1.0
2497 : ! **************************************************************************************************
2498 10278 : SUBROUTINE create_paep_section(section)
2499 :
2500 : TYPE(section_type), POINTER :: section
2501 :
2502 : TYPE(keyword_type), POINTER :: keyword
2503 :
2504 10278 : CPASSERT(.NOT. ASSOCIATED(section))
2505 :
2506 : CALL section_create(section, __LOCATION__, &
2507 : name="PLANAR_AVERAGED_V_HARTREE", &
2508 : description="Define the parameters for calculating and printing the planar "// &
2509 : "averaged electrostatic potential (Hartree potential) "// &
2510 : "for symmetric slab systems", &
2511 : n_keywords=2, &
2512 : n_subsections=0, &
2513 10278 : repeats=.FALSE.)
2514 :
2515 10278 : NULLIFY (keyword)
2516 :
2517 : CALL keyword_create(keyword, __LOCATION__, &
2518 : name="_SECTION_PARAMETERS_", &
2519 : description="Controls the activation of the planar averaged electrostatic "// &
2520 : "potential (Hartree potential) section", &
2521 : usage="&PLANAR_AVERAGED_V_HARTREE ON", &
2522 : default_l_val=.FALSE., &
2523 10278 : lone_keyword_l_val=.TRUE.)
2524 10278 : CALL section_add_keyword(section, keyword)
2525 10278 : CALL keyword_release(keyword)
2526 :
2527 : CALL keyword_create(keyword, __LOCATION__, name="PARALLEL_PLANE", &
2528 : enum_c_vals=s2a('XY', 'YZ', 'XZ'), &
2529 : enum_i_vals=[3, 1, 2], &
2530 : description="The coordinate plane that the surface is parallel to.", &
2531 : enum_desc=s2a("Parallel to XY", "Parallel to YZ", "Parallel to XZ"), &
2532 : n_var=1, &
2533 : default_i_val=3, &
2534 10278 : usage="PARALLEL_PLANE XY")
2535 10278 : CALL section_add_keyword(section, keyword)
2536 10278 : CALL keyword_release(keyword)
2537 :
2538 10278 : END SUBROUTINE create_paep_section
2539 :
2540 : END MODULE input_cp2k_dft
|