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