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