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