Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief builds the subsystem section of the input
10 : !> \par History
11 : !> 10.2005 split input_cp2k [fawzi]
12 : !> \author teo & fawzi
13 : ! **************************************************************************************************
14 : MODULE input_cp2k_subsys
15 :
16 : USE bibliography, ONLY: Goedecker1996, &
17 : Guidon2010, &
18 : Hartwigsen1998, &
19 : Krack2005, &
20 : VandeVondele2005a, &
21 : VandeVondele2007
22 : USE cell_types, ONLY: &
23 : cell_sym_cubic, cell_sym_hexagonal_gamma_120, cell_sym_hexagonal_gamma_60, &
24 : cell_sym_monoclinic, cell_sym_monoclinic_gamma_ab, cell_sym_none, cell_sym_orthorhombic, &
25 : cell_sym_rhombohedral, cell_sym_tetragonal_ab, cell_sym_tetragonal_ac, &
26 : cell_sym_tetragonal_bc, cell_sym_triclinic, use_perd_none, use_perd_x, use_perd_xy, &
27 : use_perd_xyz, use_perd_xz, use_perd_y, use_perd_yz, use_perd_z
28 : USE cp_output_handling, ONLY: cp_print_key_section_create, debug_print_level, &
29 : high_print_level, medium_print_level
30 : USE cp_units, ONLY: cp_unit_to_cp2k
31 : USE input_constants, ONLY: do_add, do_bondparm_covalent, do_bondparm_vdw, &
32 : do_cell_cif, do_cell_cp2k, do_cell_xsc, &
33 : do_cell_extxyz, do_cell_pdb, &
34 : do_conn_amb7, do_conn_g87, do_conn_g96, &
35 : do_conn_generate, do_conn_mol_set, do_conn_off, &
36 : do_conn_psf, do_conn_psf_u, do_conn_user, &
37 : do_coord_cif, do_coord_cp2k, do_coord_crd, &
38 : do_coord_g96, do_coord_off, do_coord_pdb, &
39 : do_coord_xtl, do_coord_xyz, do_remove, &
40 : do_skip_11, do_skip_12, do_skip_13, do_skip_14, &
41 : dump_pdb, gaussian
42 : USE input_cp2k_colvar, ONLY: create_colvar_section
43 : USE input_cp2k_mm, ONLY: create_neighbor_lists_section
44 : USE input_keyword_types, ONLY: keyword_create, keyword_release, keyword_type
45 : USE input_section_types, ONLY: section_add_keyword, section_add_subsection, &
46 : section_create, section_release, section_type
47 : USE input_val_types, ONLY: char_t, integer_t, lchar_t, real_t
48 : USE kinds, ONLY: dp
49 : USE physcon, ONLY: bohr
50 : USE string_utilities, ONLY: newline, s2a
51 : #include "./base/base_uses.f90"
52 :
53 : IMPLICIT NONE
54 : PRIVATE
55 :
56 : LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .TRUE.
57 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_subsys'
58 :
59 : PUBLIC :: create_subsys_section, &
60 : create_cell_section, &
61 : create_structure_data_section, &
62 : create_rng_section, &
63 : create_basis_section
64 :
65 : CONTAINS
66 :
67 : ! **************************************************************************************************
68 : !> \brief creates the cell section
69 : !> \param section ...
70 : !> \param periodic ...
71 : !> \author Ole Schuett
72 : ! **************************************************************************************************
73 20614 : SUBROUTINE create_cell_section(section, periodic)
74 : TYPE(section_type), POINTER :: section
75 : INTEGER, INTENT(IN), OPTIONAL :: periodic
76 :
77 : TYPE(section_type), POINTER :: subsection
78 :
79 20614 : CPASSERT(.NOT. ASSOCIATED(section))
80 : CALL section_create(section, __LOCATION__, "CELL", &
81 : description="Input parameters needed to set up the simulation cell. "// &
82 : "Simple products and fractions combined with functions of a single "// &
83 : "number can be used like 2/3, 0.3*COS(60) or -SQRT(3)/2. The functions "// &
84 : "COS, EXP, LOG, LOG10, SIN, SQRT, and TAN are available."//newline//newline// &
85 : "Cell settings are parsed in the following precedence order:"//newline// &
86 : "1. The external file set by CELL_FILE_NAME with a CELL_FILE_FORMAT;"//newline// &
87 : "2. The lengths and angles of cell vectors set by ABC and ALPHA_BETA_GAMMA;"//newline// &
88 : "3. The vectors set by A, B, C together;"//newline// &
89 : "4. If none above exist, the external file set by TOPOLOGY/COORD_FILE_NAME with "// &
90 : "suitable TOPOLOGY/COORD_FILE_FORMAT may also be parsed for FORCE_EVAL/SUBSYS/CELL "// &
91 20614 : "but not for FORCE_EVAL/QMMM/CELL.")
92 20614 : CALL create_cell_section_low(section, periodic)
93 :
94 20614 : NULLIFY (subsection)
95 : CALL section_create(subsection, __LOCATION__, "CELL_REF", &
96 : description="Input parameters needed to set up the reference cell for "// &
97 : "FORCE_EVAL/SUBSYS/CELL. This option can be used to keep the FFT grid "// &
98 : "fixed while running a cell optimization or NpT molecular dynamics. "// &
99 20614 : "Check the &CELL section for further details.")
100 20614 : CALL create_cell_section_low(subsection, periodic)
101 20614 : CALL section_add_subsection(section, subsection)
102 20614 : CALL section_release(subsection)
103 :
104 20614 : END SUBROUTINE create_cell_section
105 :
106 : ! **************************************************************************************************
107 : !> \brief populates cell section with keywords
108 : !> \param section ...
109 : !> \param periodic ...
110 : !> \author teo
111 : ! **************************************************************************************************
112 41228 : SUBROUTINE create_cell_section_low(section, periodic)
113 : TYPE(section_type), POINTER :: section
114 : INTEGER, INTENT(IN), OPTIONAL :: periodic
115 :
116 : INTEGER :: my_periodic
117 : TYPE(keyword_type), POINTER :: keyword
118 :
119 41228 : my_periodic = use_perd_xyz
120 41228 : IF (PRESENT(periodic)) my_periodic = periodic
121 :
122 41228 : NULLIFY (keyword)
123 : CALL keyword_create(keyword, __LOCATION__, name="A", &
124 : description="Specify the Cartesian components for the cell vector A. "// &
125 : "This defines the first column of the h matrix. "// &
126 : "Ignored if the keywords ABC or CELL_FILE_NAME are used.", &
127 : usage="A 10.000 0.000 0.000", unit_str="angstrom", &
128 41228 : n_var=3, type_of_var=real_t, repeats=.FALSE.)
129 41228 : CALL section_add_keyword(section, keyword)
130 41228 : CALL keyword_release(keyword)
131 :
132 : CALL keyword_create(keyword, __LOCATION__, name="B", &
133 : description="Specify the Cartesian components for the cell vector B. "// &
134 : "This defines the second column of the h matrix. "// &
135 : "Ignored if the keywords ABC or CELL_FILE_NAME are used.", &
136 : usage="B 0.000 10.000 0.000", unit_str="angstrom", &
137 41228 : n_var=3, type_of_var=real_t, repeats=.FALSE.)
138 41228 : CALL section_add_keyword(section, keyword)
139 41228 : CALL keyword_release(keyword)
140 :
141 : CALL keyword_create(keyword, __LOCATION__, name="C", &
142 : description="Specify the Cartesian components for the cell vector C. "// &
143 : "This defines the third column of the h matrix. "// &
144 : "Ignored if the keywords ABC or CELL_FILE_NAME are used.", &
145 : usage="C 0.000 0.000 10.000", unit_str="angstrom", &
146 41228 : n_var=3, type_of_var=real_t, repeats=.FALSE.)
147 41228 : CALL section_add_keyword(section, keyword)
148 41228 : CALL keyword_release(keyword)
149 :
150 : CALL keyword_create(keyword, __LOCATION__, name="ABC", &
151 : description="Specify the lengths of the cell vectors A, B, and C, which"// &
152 : " defines the diagonal elements of h matrix for an orthorhombic cell."// &
153 : " For non-orthorhombic cells it is possible either to specify the angles "// &
154 : "ALPHA, BETA, GAMMA via ALPHA_BETA_GAMMA keyword or alternatively use the keywords "// &
155 : "A, B, and C. The convention is that A lies along the X-axis, B is in the XY plane. "// &
156 : "Ignored if CELL_FILE_NAME is used.", &
157 : usage="ABC 10.000 10.000 10.000", unit_str="angstrom", &
158 41228 : n_var=3, type_of_var=real_t, repeats=.FALSE.)
159 41228 : CALL section_add_keyword(section, keyword)
160 41228 : CALL keyword_release(keyword)
161 :
162 : CALL keyword_create(keyword, __LOCATION__, name="ALPHA_BETA_GAMMA", &
163 : variants=["ANGLES"], &
164 : description="Specify the angles between the vectors A, B and C when using the ABC keyword. "// &
165 : "The convention is that A lies along the X-axis, B is in the XY plane. "// &
166 : "ALPHA is the angle between B and C, BETA is the angle between A and C and "// &
167 : "GAMMA is the angle between A and B.", &
168 : usage="ALPHA_BETA_GAMMA [deg] 90.0 90.0 120.0", unit_str="deg", &
169 : n_var=3, default_r_vals=[cp_unit_to_cp2k(value=90.0_dp, unit_str="deg"), &
170 : cp_unit_to_cp2k(value=90.0_dp, unit_str="deg"), &
171 : cp_unit_to_cp2k(value=90.0_dp, unit_str="deg")], &
172 206140 : repeats=.FALSE.)
173 41228 : CALL section_add_keyword(section, keyword)
174 41228 : CALL keyword_release(keyword)
175 :
176 : CALL keyword_create(keyword, __LOCATION__, name="CELL_FILE_NAME", &
177 : description="The external file from which cell is parsed ", &
178 : repeats=.FALSE., usage="CELL_FILE_NAME <CHARACTER>", &
179 41228 : type_of_var=lchar_t)
180 41228 : CALL section_add_keyword(section, keyword)
181 41228 : CALL keyword_release(keyword)
182 :
183 : CALL keyword_create(keyword, __LOCATION__, name="CELL_FILE_FORMAT", &
184 : description="Format of the external file from which "// &
185 : "cell is parsed. If the format specifies a cell by "// &
186 : "lengths and angles of three vectors, then a cell "// &
187 : "matrix is constructed with the convention that A "// &
188 : "lies along the X-axis, B is in the XY plane. ALPHA "// &
189 : "is the angle between B and C, BETA is the angle "// &
190 : "between A and C, and GAMMA is the angle between A and B.", &
191 : usage="CELL_FILE_FORMAT (CP2K|CIF|XSC|EXTXYZ|XYZ|PDB)", &
192 : enum_c_vals=s2a("CP2K", "CIF", "XSC", "EXTXYZ", "XYZ", "PDB"), &
193 : enum_i_vals=[do_cell_cp2k, do_cell_cif, do_cell_xsc, do_cell_extxyz, do_cell_extxyz, do_cell_pdb], &
194 : enum_desc=s2a("Cell info in the CP2K native format", &
195 : "Cell info from CIF file (from fields `_cell_length_a` or `_cell.length_a`, etc)", &
196 : "Cell info in the XSC format (NAMD)", &
197 : "Cell info as `lattice=...` field in the comment line of Extended XYZ format", &
198 : "Alias for Extended XYZ", &
199 : "Cell info in the `CRYST1` record of PDB format"), &
200 41228 : default_i_val=do_cell_cp2k)
201 41228 : CALL section_add_keyword(section, keyword)
202 41228 : CALL keyword_release(keyword)
203 :
204 : CALL keyword_create(keyword, __LOCATION__, name="PERIODIC", &
205 : description="Specify the directions for which periodic boundary conditions (PBC) will be applied. "// &
206 : "Important notice: This applies to the generation of the pair lists as well as to the "// &
207 : "application of the PBCs to positions. "// &
208 : "See the POISSON section to specify the periodicity used for the electrostatics. "// &
209 : "Typically the settings should be the same.", &
210 : usage="PERIODIC (x|y|z|xy|xz|yz|xyz|none)", &
211 : enum_c_vals=s2a("x", "y", "z", "xy", "xz", "yz", "xyz", "none"), &
212 : enum_i_vals=[use_perd_x, use_perd_y, use_perd_z, &
213 : use_perd_xy, use_perd_xz, use_perd_yz, &
214 : use_perd_xyz, use_perd_none], &
215 41228 : default_i_val=my_periodic)
216 41228 : CALL section_add_keyword(section, keyword)
217 41228 : CALL keyword_release(keyword)
218 :
219 : CALL keyword_create(keyword, __LOCATION__, name="MULTIPLE_UNIT_CELL", &
220 : description="Specifies the numbers of repetition in space (X, Y, Z) of the defined cell, "// &
221 : "assuming it as a unit cell. This keyword affects only the CELL specification. The same keyword "// &
222 : "in SUBSYS%TOPOLOGY%MULTIPLE_UNIT_CELL should be modified in order to affect the coordinates "// &
223 : "specification.", usage="MULTIPLE_UNIT_CELL 1 1 1", &
224 41228 : n_var=3, default_i_vals=[1, 1, 1], repeats=.FALSE.)
225 41228 : CALL section_add_keyword(section, keyword)
226 41228 : CALL keyword_release(keyword)
227 :
228 : CALL keyword_create( &
229 : keyword, __LOCATION__, name="SYMMETRY", &
230 : description="Imposes an initial cell symmetry, according to the convention "// &
231 : "that A lies along the X-axis, B is in the XY plane. After the "// &
232 : "input cell information is parsed, the symmetry is enforced by "// &
233 : "reconstructing the cell matrix from lengths and angles of the "// &
234 : "cell vectors, taking averages if necessary. This process does "// &
235 : "not affect input atomic coordinates; in case a space group is "// &
236 : "to be detected and preserved for an optimization task, atomic "// &
237 : "coordinates should correspond to cell vectors already obeying "// &
238 : "the convention mentioned above.", &
239 : usage="SYMMETRY monoclinic", &
240 : enum_desc=s2a("No cell symmetry", &
241 : "Triclinic (a ≠ b ≠ c ≠ a, α ≠ β ≠ γ ≠ α ≠ 90°)", &
242 : "Monoclinic (a ≠ b ≠ c, α = γ = 90°, β ≠ 90°)", &
243 : "Monoclinic (a = b ≠ c, α = β = 90°, γ ≠ 90°)", &
244 : "Orthorhombic (a ≠ b ≠ c, α = β = γ = 90°)", &
245 : "Tetragonal (a = b ≠ c, α = β = γ = 90°)", &
246 : "Tetragonal (a = c ≠ b, α = β = γ = 90°)", &
247 : "Tetragonal (a ≠ b = c, α = β = γ = 90°)", &
248 : "Tetragonal (alias for TETRAGONAL_AB)", &
249 : "Rhombohedral (a = b = c, α = β = γ ≠ 90°)", &
250 : "Hexagonal (alias for HEXAGONAL_GAMMA_60)", &
251 : "Hexagonal (a = b ≠ c, α = β = 90°, γ = 60°)", &
252 : "Hexagonal (a = b ≠ c, α = β = 90°, γ = 120°)", &
253 : "Cubic (a = b = c, α = β = γ = 90°)"), &
254 : enum_c_vals=s2a("NONE", "TRICLINIC", "MONOCLINIC", "MONOCLINIC_GAMMA_AB", "ORTHORHOMBIC", &
255 : "TETRAGONAL_AB", "TETRAGONAL_AC", "TETRAGONAL_BC", "TETRAGONAL", "RHOMBOHEDRAL", &
256 : "HEXAGONAL", "HEXAGONAL_GAMMA_60", "HEXAGONAL_GAMMA_120", "CUBIC"), &
257 : enum_i_vals=[cell_sym_none, cell_sym_triclinic, cell_sym_monoclinic, cell_sym_monoclinic_gamma_ab, &
258 : cell_sym_orthorhombic, cell_sym_tetragonal_ab, cell_sym_tetragonal_ac, cell_sym_tetragonal_bc, &
259 : cell_sym_tetragonal_ab, cell_sym_rhombohedral, cell_sym_hexagonal_gamma_60, &
260 : cell_sym_hexagonal_gamma_60, cell_sym_hexagonal_gamma_120, cell_sym_cubic], &
261 41228 : default_i_val=cell_sym_none)
262 41228 : CALL section_add_keyword(section, keyword)
263 41228 : CALL keyword_release(keyword)
264 :
265 41228 : END SUBROUTINE create_cell_section_low
266 :
267 : ! **************************************************************************************************
268 : !> \brief Creates the random number restart section
269 : !> \param section the section to create
270 : !> \author teo
271 : ! **************************************************************************************************
272 227826 : SUBROUTINE create_rng_section(section)
273 : TYPE(section_type), POINTER :: section
274 :
275 : TYPE(keyword_type), POINTER :: keyword
276 :
277 227826 : CPASSERT(.NOT. ASSOCIATED(section))
278 : CALL section_create(section, __LOCATION__, name="RNG_INIT", &
279 : description="Information to initialize the parallel random number generator streams", &
280 227826 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
281 227826 : NULLIFY (keyword)
282 :
283 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
284 : description="Specify an initial RNG stream record", repeats=.TRUE., &
285 227826 : usage="{RNG record string}", type_of_var=lchar_t)
286 227826 : CALL section_add_keyword(section, keyword)
287 227826 : CALL keyword_release(keyword)
288 :
289 227826 : END SUBROUTINE create_rng_section
290 :
291 : ! **************************************************************************************************
292 : !> \brief creates the structure of a subsys, i.e. a full set of
293 : !> atoms+mol+bounds+cell
294 : !> \param section the section to create
295 : !> \author fawzi
296 : ! **************************************************************************************************
297 9582 : SUBROUTINE create_subsys_section(section)
298 : TYPE(section_type), POINTER :: section
299 :
300 : TYPE(keyword_type), POINTER :: keyword
301 : TYPE(section_type), POINTER :: subsection
302 :
303 9582 : CPASSERT(.NOT. ASSOCIATED(section))
304 : CALL section_create(section, __LOCATION__, name="subsys", &
305 : description="a subsystem: coordinates, topology, molecules and cell", &
306 9582 : n_keywords=1, n_subsections=9, repeats=.FALSE.)
307 :
308 9582 : NULLIFY (keyword)
309 : CALL keyword_create(keyword, __LOCATION__, name="SEED", &
310 : description="Initial seed for the (pseudo)random number generator for the "// &
311 : "Wiener process employed by the Langevin dynamics. Exactly 1 or 6 positive "// &
312 : "integer values are expected. A single value is replicated to fill up the "// &
313 : "full seed array with 6 numbers.", &
314 : n_var=-1, &
315 : type_of_var=integer_t, &
316 : usage="SEED {INTEGER} .. {INTEGER}", &
317 9582 : default_i_vals=[12345])
318 9582 : CALL section_add_keyword(section, keyword)
319 9582 : CALL keyword_release(keyword)
320 :
321 9582 : NULLIFY (subsection)
322 :
323 9582 : CALL create_rng_section(subsection)
324 9582 : CALL section_add_subsection(section, subsection)
325 9582 : CALL section_release(subsection)
326 :
327 9582 : CALL create_cell_section(subsection)
328 9582 : CALL section_add_subsection(section, subsection)
329 9582 : CALL section_release(subsection)
330 :
331 9582 : CALL create_coord_section(subsection)
332 9582 : CALL section_add_subsection(section, subsection)
333 9582 : CALL section_release(subsection)
334 :
335 9582 : CALL create_velocity_section(subsection)
336 9582 : CALL section_add_subsection(section, subsection)
337 9582 : CALL section_release(subsection)
338 :
339 9582 : CALL create_kind_section(subsection)
340 9582 : CALL section_add_subsection(section, subsection)
341 9582 : CALL section_release(subsection)
342 :
343 9582 : CALL create_topology_section(subsection)
344 9582 : CALL section_add_subsection(section, subsection)
345 9582 : CALL section_release(subsection)
346 :
347 9582 : CALL create_colvar_section(section=subsection)
348 9582 : CALL section_add_subsection(section, subsection)
349 9582 : CALL section_release(subsection)
350 :
351 9582 : CALL create_multipole_section(subsection)
352 9582 : CALL section_add_subsection(section, subsection)
353 9582 : CALL section_release(subsection)
354 :
355 9582 : CALL create_shell_coord_section(subsection)
356 9582 : CALL section_add_subsection(section, subsection)
357 9582 : CALL section_release(subsection)
358 :
359 9582 : CALL create_shell_vel_section(subsection)
360 9582 : CALL section_add_subsection(section, subsection)
361 9582 : CALL section_release(subsection)
362 9582 : CALL create_core_coord_section(subsection)
363 9582 : CALL section_add_subsection(section, subsection)
364 9582 : CALL section_release(subsection)
365 :
366 9582 : CALL create_core_vel_section(subsection)
367 9582 : CALL section_add_subsection(section, subsection)
368 9582 : CALL section_release(subsection)
369 :
370 9582 : CALL create_subsys_print_section(subsection)
371 9582 : CALL section_add_subsection(section, subsection)
372 9582 : CALL section_release(subsection)
373 :
374 9582 : END SUBROUTINE create_subsys_section
375 :
376 : ! **************************************************************************************************
377 : !> \brief Creates the subsys print section
378 : !> \param section the section to create
379 : !> \author teo
380 : ! **************************************************************************************************
381 9582 : SUBROUTINE create_subsys_print_section(section)
382 : TYPE(section_type), POINTER :: section
383 :
384 : TYPE(keyword_type), POINTER :: keyword
385 : TYPE(section_type), POINTER :: print_key
386 :
387 9582 : NULLIFY (print_key, keyword)
388 9582 : CPASSERT(.NOT. ASSOCIATED(section))
389 : CALL section_create(section, __LOCATION__, name="print", &
390 : description="Controls printings related to the subsys", &
391 9582 : n_keywords=0, n_subsections=9, repeats=.FALSE.)
392 :
393 : CALL cp_print_key_section_create(print_key, __LOCATION__, "atomic_coordinates", &
394 : description="controls the output of the atomic coordinates when setting up the"// &
395 : " force environment. For printing coordinates during MD or GEO refer to the keyword"// &
396 : " trajectory.", unit_str="angstrom", &
397 9582 : print_level=medium_print_level, filename="__STD_OUT__")
398 9582 : CALL section_add_subsection(section, print_key)
399 9582 : CALL section_release(print_key)
400 :
401 9582 : CALL create_structure_data_section(print_key)
402 9582 : CALL section_add_subsection(section, print_key)
403 9582 : CALL section_release(print_key)
404 :
405 : CALL cp_print_key_section_create(print_key, __LOCATION__, "INTERATOMIC_DISTANCES", &
406 : description="Controls the printout of the interatomic distances when setting up the "// &
407 : "force environment", unit_str="angstrom", &
408 9582 : print_level=debug_print_level, filename="__STD_OUT__")
409 : CALL keyword_create(keyword, __LOCATION__, name="CHECK_INTERATOMIC_DISTANCES", &
410 : description="Minimum allowed distance between two atoms. "// &
411 : "A warning is printed, if a smaller interatomic distance is encountered. "// &
412 : "The check is disabled for the threshold value 0 which is the default "// &
413 : "for systems with more than 2000 atoms (otherwise 0.5 A). "// &
414 : "The run is aborted, if an interatomic distance is smaller than the absolute "// &
415 : "value of a negative threshold value.", &
416 9582 : default_r_val=0.5_dp*bohr, unit_str="angstrom")
417 9582 : CALL section_add_keyword(print_key, keyword)
418 9582 : CALL keyword_release(keyword)
419 9582 : CALL section_add_subsection(section, print_key)
420 9582 : CALL section_release(print_key)
421 :
422 : CALL cp_print_key_section_create(print_key, __LOCATION__, "topology_info", description= &
423 : "controls the printing of information in the topology settings", &
424 9582 : print_level=high_print_level, filename="__STD_OUT__")
425 : CALL keyword_create(keyword, __LOCATION__, name="xtl_info", &
426 : description="Prints information when parsing XTL files.", &
427 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
428 9582 : CALL section_add_keyword(print_key, keyword)
429 9582 : CALL keyword_release(keyword)
430 : CALL keyword_create(keyword, __LOCATION__, name="cif_info", &
431 : description="Prints information when parsing CIF files.", &
432 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
433 9582 : CALL section_add_keyword(print_key, keyword)
434 9582 : CALL keyword_release(keyword)
435 : CALL keyword_create(keyword, __LOCATION__, name="pdb_info", &
436 : description="Prints information when parsing PDB files.", &
437 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
438 9582 : CALL section_add_keyword(print_key, keyword)
439 9582 : CALL keyword_release(keyword)
440 : CALL keyword_create(keyword, __LOCATION__, name="xyz_info", &
441 : description="Prints information when parsing XYZ files.", &
442 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
443 9582 : CALL section_add_keyword(print_key, keyword)
444 9582 : CALL keyword_release(keyword)
445 : CALL keyword_create(keyword, __LOCATION__, name="psf_info", &
446 : description="Prints information when parsing PSF files.", &
447 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
448 9582 : CALL section_add_keyword(print_key, keyword)
449 9582 : CALL keyword_release(keyword)
450 : CALL keyword_create(keyword, __LOCATION__, name="amber_info", &
451 : description="Prints information when parsing ABER topology files.", &
452 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
453 9582 : CALL section_add_keyword(print_key, keyword)
454 9582 : CALL keyword_release(keyword)
455 : CALL keyword_create(keyword, __LOCATION__, name="g96_info", &
456 : description="Prints information when parsing G96 files.", &
457 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
458 9582 : CALL section_add_keyword(print_key, keyword)
459 9582 : CALL keyword_release(keyword)
460 : CALL keyword_create(keyword, __LOCATION__, name="crd_info", &
461 : description="Prints information when parsing CRD files.", &
462 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
463 9582 : CALL section_add_keyword(print_key, keyword)
464 9582 : CALL keyword_release(keyword)
465 : CALL keyword_create(keyword, __LOCATION__, name="gtop_info", &
466 : description="Prints information when parsing GROMOS topology files.", &
467 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
468 9582 : CALL section_add_keyword(print_key, keyword)
469 9582 : CALL keyword_release(keyword)
470 : CALL keyword_create(keyword, __LOCATION__, name="util_info", &
471 : description="Prints information regarding topology utilities", &
472 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
473 9582 : CALL section_add_keyword(print_key, keyword)
474 9582 : CALL keyword_release(keyword)
475 : CALL keyword_create(keyword, __LOCATION__, name="generate_info", &
476 : description="Prints information regarding topology generation", &
477 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
478 9582 : CALL section_add_keyword(print_key, keyword)
479 9582 : CALL keyword_release(keyword)
480 9582 : CALL section_add_subsection(section, print_key)
481 9582 : CALL section_release(print_key)
482 :
483 : CALL cp_print_key_section_create(print_key, __LOCATION__, "cell", &
484 : description="controls the output of the cell parameters", &
485 : print_level=medium_print_level, filename="__STD_OUT__", &
486 9582 : unit_str="angstrom")
487 9582 : CALL section_add_subsection(section, print_key)
488 9582 : CALL section_release(print_key)
489 :
490 : CALL cp_print_key_section_create(print_key, __LOCATION__, "kinds", &
491 : description="controls the output of information on the kinds", &
492 9582 : print_level=medium_print_level, filename="__STD_OUT__")
493 : CALL keyword_create(keyword, __LOCATION__, name="potential", &
494 : description="If the printkey is activated controls the printing of the"// &
495 : " fist_potential, gth_potential, sgp_potential or all electron"// &
496 : " potential information", &
497 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
498 9582 : CALL section_add_keyword(print_key, keyword)
499 9582 : CALL keyword_release(keyword)
500 : CALL keyword_create(keyword, __LOCATION__, name="basis_set", &
501 : description="If the printkey is activated controls the printing of basis set information", &
502 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
503 9582 : CALL section_add_keyword(print_key, keyword)
504 9582 : CALL keyword_release(keyword)
505 : CALL keyword_create(keyword, __LOCATION__, name="se_parameters", &
506 : description="If the printkey is activated controls the printing of the semi-empirical parameters.", &
507 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
508 9582 : CALL section_add_keyword(print_key, keyword)
509 9582 : CALL keyword_release(keyword)
510 9582 : CALL section_add_subsection(section, print_key)
511 9582 : CALL section_release(print_key)
512 :
513 : CALL cp_print_key_section_create(print_key, __LOCATION__, "SYMMETRY", &
514 : description="controls the output of symmetry information", &
515 9582 : print_level=debug_print_level + 1, filename="__STD_OUT__")
516 : CALL keyword_create(keyword, __LOCATION__, name="MOLECULE", &
517 : description="Assume the system is an isolated molecule", &
518 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
519 9582 : CALL section_add_keyword(print_key, keyword)
520 9582 : CALL keyword_release(keyword)
521 : CALL keyword_create(keyword, __LOCATION__, name="EPS_GEO", &
522 : description="Accuracy required for symmetry detection", &
523 9582 : default_r_val=1.0E-4_dp)
524 9582 : CALL section_add_keyword(print_key, keyword)
525 9582 : CALL keyword_release(keyword)
526 : CALL keyword_create(keyword, __LOCATION__, name="STANDARD_ORIENTATION", &
527 : description="Print molecular coordinates in standard orientation", &
528 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
529 9582 : CALL section_add_keyword(print_key, keyword)
530 9582 : CALL keyword_release(keyword)
531 : CALL keyword_create(keyword, __LOCATION__, name="INERTIA", &
532 : description="Print molecular inertia tensor", &
533 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
534 9582 : CALL section_add_keyword(print_key, keyword)
535 9582 : CALL keyword_release(keyword)
536 : CALL keyword_create(keyword, __LOCATION__, name="SYMMETRY_ELEMENTS", &
537 : description="Print symmetry elements", &
538 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
539 9582 : CALL section_add_keyword(print_key, keyword)
540 9582 : CALL keyword_release(keyword)
541 : CALL keyword_create(keyword, __LOCATION__, name="ALL", &
542 : description="Print all symmetry information", &
543 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
544 9582 : CALL section_add_keyword(print_key, keyword)
545 9582 : CALL keyword_release(keyword)
546 : CALL keyword_create(keyword, __LOCATION__, name="ROTATION_MATRICES", &
547 : description="All the rotation matrices of the point group", &
548 9582 : default_l_val=.FALSE.)
549 9582 : CALL section_add_keyword(print_key, keyword)
550 9582 : CALL keyword_release(keyword)
551 : CALL keyword_create(keyword, __LOCATION__, name="CHECK_SYMMETRY", &
552 : description="Check if calculated symmetry has expected value."// &
553 : " Use either Schoenfliess or Hermann-Maugin symbols", &
554 9582 : default_c_val="NONE")
555 9582 : CALL section_add_keyword(print_key, keyword)
556 9582 : CALL keyword_release(keyword)
557 9582 : CALL section_add_subsection(section, print_key)
558 9582 : CALL section_release(print_key)
559 :
560 : CALL cp_print_key_section_create(print_key, __LOCATION__, "molecules", &
561 : description="controls the output of information on the molecules", &
562 9582 : print_level=medium_print_level, filename="__STD_OUT__")
563 9582 : CALL section_add_subsection(section, print_key)
564 9582 : CALL section_release(print_key)
565 :
566 : CALL cp_print_key_section_create(print_key, __LOCATION__, "radii", &
567 : description="controls the output of radii information", unit_str="angstrom", &
568 9582 : print_level=high_print_level, filename="__STD_OUT__")
569 :
570 : CALL keyword_create(keyword, __LOCATION__, name="core_charges_radii", &
571 : description="If the printkey is activated controls the printing of the radii of the core charges", &
572 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
573 9582 : CALL section_add_keyword(print_key, keyword)
574 9582 : CALL keyword_release(keyword)
575 :
576 : CALL keyword_create(keyword, __LOCATION__, name="pgf_radii", &
577 : description="If the printkey is activated controls the printing of the core gaussian radii", &
578 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
579 9582 : CALL section_add_keyword(print_key, keyword)
580 9582 : CALL keyword_release(keyword)
581 :
582 : CALL keyword_create(keyword, __LOCATION__, name="set_radii", &
583 : description="If the printkey is activated controls the printing of the set_radii", &
584 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
585 9582 : CALL section_add_keyword(print_key, keyword)
586 9582 : CALL keyword_release(keyword)
587 :
588 : CALL keyword_create(keyword, __LOCATION__, name="kind_radii", &
589 : description="If the printkey is activated controls the printing of the kind_radii", &
590 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
591 9582 : CALL section_add_keyword(print_key, keyword)
592 9582 : CALL keyword_release(keyword)
593 :
594 : CALL keyword_create(keyword, __LOCATION__, name="core_charge_radii", &
595 : description="If the printkey is activated controls the printing of the core_charge_radii", &
596 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
597 9582 : CALL section_add_keyword(print_key, keyword)
598 9582 : CALL keyword_release(keyword)
599 :
600 : CALL keyword_create(keyword, __LOCATION__, name="ppl_radii", &
601 : description="If the printkey is activated controls the printing of the "// &
602 : "pseudo potential local radii", &
603 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
604 9582 : CALL section_add_keyword(print_key, keyword)
605 9582 : CALL keyword_release(keyword)
606 :
607 : CALL keyword_create(keyword, __LOCATION__, name="ppnl_radii", &
608 : description="If the printkey is activated controls the printing of the "// &
609 : "pseudo potential non local radii", &
610 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
611 9582 : CALL section_add_keyword(print_key, keyword)
612 9582 : CALL keyword_release(keyword)
613 :
614 : CALL keyword_create(keyword, __LOCATION__, name="gapw_prj_radii", &
615 : description="If the printkey is activated controls the printing of the gapw projector radii", &
616 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
617 9582 : CALL section_add_keyword(print_key, keyword)
618 9582 : CALL keyword_release(keyword)
619 :
620 9582 : CALL section_add_subsection(section, print_key)
621 9582 : CALL section_release(print_key)
622 :
623 9582 : END SUBROUTINE create_subsys_print_section
624 :
625 : ! **************************************************************************************************
626 : !> \brief Creates the multipole section
627 : !> \param section the section to create
628 : !> \author teo
629 : ! **************************************************************************************************
630 9582 : SUBROUTINE create_multipole_section(section)
631 : TYPE(section_type), POINTER :: section
632 :
633 : TYPE(keyword_type), POINTER :: keyword
634 : TYPE(section_type), POINTER :: subsection
635 :
636 9582 : CPASSERT(.NOT. ASSOCIATED(section))
637 : CALL section_create(section, __LOCATION__, name="multipoles", &
638 : description="Specifies the dipoles and quadrupoles for particles.", &
639 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
640 :
641 9582 : NULLIFY (keyword, subsection)
642 : CALL section_create(subsection, __LOCATION__, name="dipoles", &
643 : description="Specifies the dipoles of the particles.", &
644 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
645 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
646 : description="The dipole components for each atom in the format: "// &
647 : "$D_x \ D_y \ D_z$", &
648 : repeats=.TRUE., usage="{Real} {Real} {Real}", &
649 9582 : type_of_var=real_t, n_var=3)
650 9582 : CALL section_add_keyword(subsection, keyword)
651 9582 : CALL keyword_release(keyword)
652 9582 : CALL section_add_subsection(section, subsection)
653 9582 : CALL section_release(subsection)
654 :
655 : CALL section_create(subsection, __LOCATION__, name="quadrupoles", &
656 : description="Specifies the quadrupoles of the particles.", &
657 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
658 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
659 : description="The quadrupole components for each atom in the format: "// &
660 : "$Q_{xx} \ Q_{xy} \ Q_{xz} \ Q_{yy} \ Q_{yz} \ Q_{zz}$", &
661 : repeats=.TRUE., usage="{Real} {Real} {Real} {Real} {Real} {Real}", &
662 9582 : type_of_var=real_t, n_var=6)
663 9582 : CALL section_add_keyword(subsection, keyword)
664 9582 : CALL keyword_release(keyword)
665 9582 : CALL section_add_subsection(section, subsection)
666 9582 : CALL section_release(subsection)
667 :
668 9582 : END SUBROUTINE create_multipole_section
669 :
670 : ! **************************************************************************************************
671 : !> \brief creates structure data section for output.. both subsys (for initialization)
672 : !> and motion section..
673 : !> \param print_key ...
674 : ! **************************************************************************************************
675 19164 : SUBROUTINE create_structure_data_section(print_key)
676 : TYPE(section_type), POINTER :: print_key
677 :
678 : TYPE(keyword_type), POINTER :: keyword
679 :
680 19164 : CPASSERT(.NOT. ASSOCIATED(print_key))
681 :
682 19164 : NULLIFY (keyword)
683 :
684 : CALL cp_print_key_section_create(print_key, __LOCATION__, name="STRUCTURE_DATA", &
685 : description="Request the printing of special structure data during a structure "// &
686 : "optimization (in MOTION%PRINT) or when setting up a subsys (in SUBSYS%PRINT).", &
687 19164 : print_level=high_print_level, filename="__STD_OUT__", unit_str="angstrom")
688 :
689 : CALL keyword_create(keyword, __LOCATION__, name="POSITION", variants=["POS"], &
690 : description="Print the position vectors in Cartesian coordinates of the atoms specified "// &
691 : "by a list of their indices", &
692 : usage="POSITION {integer} {integer} {integer}..{integer}", n_var=-1, repeats=.TRUE., &
693 38328 : type_of_var=integer_t)
694 19164 : CALL section_add_keyword(print_key, keyword)
695 19164 : CALL keyword_release(keyword)
696 :
697 : CALL keyword_create(keyword, __LOCATION__, name="POSITION_SCALED", variants=["POS_SCALED"], &
698 : description="Print the position vectors in scaled coordinates of the atoms specified "// &
699 : "by a list of their indices", &
700 : usage="POSITION_SCALED {integer} {integer} {integer}..{integer}", n_var=-1, repeats=.TRUE., &
701 38328 : type_of_var=integer_t)
702 19164 : CALL section_add_keyword(print_key, keyword)
703 19164 : CALL keyword_release(keyword)
704 :
705 : CALL keyword_create(keyword, __LOCATION__, name="DISTANCE", variants=["DIS"], &
706 : description="Print the distance between the atoms a and b specified by their indices", &
707 : usage="DISTANCE {integer} {integer}", n_var=2, repeats=.TRUE., &
708 38328 : type_of_var=integer_t)
709 19164 : CALL section_add_keyword(print_key, keyword)
710 19164 : CALL keyword_release(keyword)
711 :
712 : CALL keyword_create(keyword, __LOCATION__, name="ANGLE", variants=["ANG"], &
713 : description="Print the angle formed by the atoms specified by their indices", &
714 : usage="ANGLE {integer} {integer} {integer}", n_var=3, repeats=.TRUE., &
715 38328 : type_of_var=integer_t)
716 19164 : CALL section_add_keyword(print_key, keyword)
717 19164 : CALL keyword_release(keyword)
718 :
719 : CALL keyword_create(keyword, __LOCATION__, name="DIHEDRAL_ANGLE", variants=s2a("DIHEDRAL", "DIH"), &
720 : description="Print the dihedral angle between the planes defined by the atoms (a,b,c) and "// &
721 : "the atoms (b,c,d) specified by their indices", &
722 : usage="DIHEDRAL_ANGLE {integer} {integer} {integer} {integer}", n_var=4, &
723 19164 : repeats=.TRUE., type_of_var=integer_t)
724 19164 : CALL section_add_keyword(print_key, keyword)
725 19164 : CALL keyword_release(keyword)
726 :
727 19164 : END SUBROUTINE create_structure_data_section
728 :
729 : ! **************************************************************************************************
730 : !> \brief Creates the velocity section
731 : !> \param section the section to create
732 : !> \author teo
733 : ! **************************************************************************************************
734 9582 : SUBROUTINE create_velocity_section(section)
735 : TYPE(section_type), POINTER :: section
736 :
737 : TYPE(keyword_type), POINTER :: keyword
738 :
739 9582 : CPASSERT(.NOT. ASSOCIATED(section))
740 : CALL section_create(section, __LOCATION__, name="velocity", &
741 : description="The velocities for simple systems or "// &
742 : "the centroid mode in PI runs, xyz format by default", &
743 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
744 9582 : NULLIFY (keyword)
745 : CALL keyword_create(keyword, __LOCATION__, name="PINT_UNIT", &
746 : description="Specify the units of measurement for the velocities "// &
747 : "(currently works only for the path integral code). "// &
748 : "All available CP2K units can be used.", &
749 : usage="PINT_UNIT angstrom*au_t^-1", &
750 9582 : default_c_val="bohr*au_t^-1")
751 9582 : CALL section_add_keyword(section, keyword)
752 9582 : CALL keyword_release(keyword)
753 :
754 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
755 : description="The atomic velocities in the format: "// &
756 : "$ v_x \ v_y \ v_z$ "// &
757 : "The same order as for the atomic coordinates is assumed.", &
758 : repeats=.TRUE., usage="{Real} {Real} {Real}", &
759 9582 : type_of_var=real_t, n_var=3)
760 9582 : CALL section_add_keyword(section, keyword)
761 9582 : CALL keyword_release(keyword)
762 :
763 9582 : END SUBROUTINE create_velocity_section
764 :
765 : ! **************************************************************************************************
766 : !> \brief Creates the shell velocity section
767 : !> \param section the section to create
768 : !> \author teo
769 : ! **************************************************************************************************
770 9582 : SUBROUTINE create_shell_vel_section(section)
771 : TYPE(section_type), POINTER :: section
772 :
773 : TYPE(keyword_type), POINTER :: keyword
774 :
775 9582 : CPASSERT(.NOT. ASSOCIATED(section))
776 : CALL section_create(section, __LOCATION__, name="shell_velocity", &
777 : description="The velocities of shells for shell-model potentials, "// &
778 : "in xyz format ", &
779 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
780 9582 : NULLIFY (keyword)
781 :
782 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
783 : description="The shell particle velocities in the format: "// &
784 : "$v_x \ v_y \ v_z$ "// &
785 : "The same order as for the shell particle coordinates is assumed.", &
786 : repeats=.TRUE., usage="{Real} {Real} {Real}", &
787 9582 : type_of_var=real_t, n_var=3)
788 9582 : CALL section_add_keyword(section, keyword)
789 9582 : CALL keyword_release(keyword)
790 :
791 9582 : END SUBROUTINE create_shell_vel_section
792 :
793 : ! **************************************************************************************************
794 : !> \brief Creates the shell velocity section
795 : !> \param section the section to create
796 : !> \author teo
797 : ! **************************************************************************************************
798 9582 : SUBROUTINE create_core_vel_section(section)
799 : TYPE(section_type), POINTER :: section
800 :
801 : TYPE(keyword_type), POINTER :: keyword
802 :
803 9582 : CPASSERT(.NOT. ASSOCIATED(section))
804 : CALL section_create(section, __LOCATION__, name="core_velocity", &
805 : description="The velocities of cores for shell-model potentials, "// &
806 : "in xyz format ", &
807 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
808 9582 : NULLIFY (keyword)
809 :
810 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
811 : description="The core particle velocities in the format: "// &
812 : "$v_x \ v_y \ v_z$ "// &
813 : "The same order as for the core particle coordinates is assumed.", &
814 : repeats=.TRUE., usage="{Real} {Real} {Real}", &
815 9582 : type_of_var=real_t, n_var=3)
816 9582 : CALL section_add_keyword(section, keyword)
817 9582 : CALL keyword_release(keyword)
818 :
819 9582 : END SUBROUTINE create_core_vel_section
820 :
821 : ! **************************************************************************************************
822 : !> \brief Creates the &POTENTIAL section
823 : !> \param section the section to create
824 : !> \author teo
825 : ! **************************************************************************************************
826 9582 : SUBROUTINE create_potential_section(section)
827 : TYPE(section_type), POINTER :: section
828 :
829 : TYPE(keyword_type), POINTER :: keyword
830 :
831 : CALL section_create(section, __LOCATION__, name="potential", &
832 : description="Section used to specify Potentials.", &
833 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
834 9582 : NULLIFY (keyword)
835 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
836 : description="CP2K Pseudo Potential Standard Format (GTH, ALL)", &
837 9582 : repeats=.TRUE., type_of_var=lchar_t)
838 9582 : CALL section_add_keyword(section, keyword)
839 9582 : CALL keyword_release(keyword)
840 :
841 9582 : END SUBROUTINE create_potential_section
842 :
843 : ! **************************************************************************************************
844 : !> \brief Creates the &KG_POTENTIAL section
845 : !> \param section the section to create
846 : !> \author JGH
847 : ! **************************************************************************************************
848 9582 : SUBROUTINE create_kgpot_section(section)
849 : TYPE(section_type), POINTER :: section
850 :
851 : TYPE(keyword_type), POINTER :: keyword
852 :
853 : CALL section_create(section, __LOCATION__, name="kg_potential", &
854 : description="Section used to specify KG Potentials.", &
855 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
856 9582 : NULLIFY (keyword)
857 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
858 : description="CP2K KG TNADD Potential Standard Format (TNADD)", &
859 9582 : repeats=.TRUE., type_of_var=lchar_t)
860 9582 : CALL section_add_keyword(section, keyword)
861 9582 : CALL keyword_release(keyword)
862 :
863 9582 : END SUBROUTINE create_kgpot_section
864 :
865 : ! **************************************************************************************************
866 : !> \brief Creates the &BASIS section
867 : !> \param section the section to create
868 : !> \author teo
869 : ! **************************************************************************************************
870 19164 : SUBROUTINE create_basis_section(section)
871 : TYPE(section_type), POINTER :: section
872 :
873 : TYPE(keyword_type), POINTER :: keyword
874 :
875 : CALL section_create(section, __LOCATION__, name="BASIS", &
876 : description="Section used to specify a general basis set for QM calculations.", &
877 19164 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
878 :
879 19164 : NULLIFY (keyword)
880 :
881 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
882 : description="The type of basis set defined in this section.", &
883 : lone_keyword_c_val="Orbital", &
884 19164 : usage="Orbital", default_c_val="Orbital")
885 19164 : CALL section_add_keyword(section, keyword)
886 19164 : CALL keyword_release(keyword)
887 :
888 : CALL keyword_create( &
889 : keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
890 : repeats=.TRUE., type_of_var=lchar_t, &
891 : description="CP2K Basis Set Standard Format:"//newline//newline// &
892 : "```"//newline// &
893 : "Element symbol Name of the basis set Alias names"//newline// &
894 : "nset (repeat the following block of lines nset times)"//newline// &
895 : "n lmin lmax nexp nshell(lmin) nshell(lmin+1) ... nshell(lmax-1) nshell(lmax)"//newline// &
896 : "a(1) c(1,l,1) c(1,l,2) ... c(1,l,nshell(l)-1) c(1,l,nshell(l)), l=lmin,lmax"//newline// &
897 : "a(2) c(2,l,1) c(2,l,2) ... c(2,l,nshell(l)-1) c(2,l,nshell(l)), l=lmin,lmax"//newline// &
898 : " . . . . ."//newline// &
899 : " . . . . ."//newline// &
900 : " . . . . ."//newline// &
901 : "a(nexp-1) c(nexp-1,l,1) c(nexp-1,l,2) ... c(nexp-1,l,nshell(l)-1) c(nexp-1,l,nshell(l)), l=lmin,lmax"//newline// &
902 : "a(nexp) c(nexp,l,1) c(nexp,l,2) ... c(nexp,l,nshell(l)-1) c(nexp,l,nshell(l)), l=lmin,lmax"//newline// &
903 : newline// &
904 : newline// &
905 : "nset : Number of exponent sets"//newline// &
906 : "n : Principle quantum number (only for orbital label printing)"//newline// &
907 : "lmax : Maximum angular momentum quantum number l"//newline// &
908 : "lmin : Minimum angular momentum quantum number l"//newline// &
909 : "nshell(l): Number of shells for angular momentum quantum number l"//newline// &
910 : "a : Exponent"//newline// &
911 : "c : Contraction coefficient"//newline// &
912 19164 : "```")
913 19164 : CALL section_add_keyword(section, keyword)
914 19164 : CALL keyword_release(keyword)
915 :
916 19164 : END SUBROUTINE create_basis_section
917 :
918 : ! **************************************************************************************************
919 : !> \brief Creates the &COORD section
920 : !> \param section the section to create
921 : !> \author teo
922 : ! **************************************************************************************************
923 9582 : SUBROUTINE create_coord_section(section)
924 : TYPE(section_type), POINTER :: section
925 :
926 : TYPE(keyword_type), POINTER :: keyword
927 :
928 9582 : CPASSERT(.NOT. ASSOCIATED(section))
929 : CALL section_create(section, __LOCATION__, name="coord", &
930 : description="The coordinates for simple systems (like small QM cells) "// &
931 : "are specified here by default using explicit XYZ coordinates. "// &
932 : "Simple products and fractions combined with functions of a single "// &
933 : "number can be used like 2/3, 0.3*COS(60) or -SQRT(3)/2. "// &
934 : "More complex systems should be given via an external coordinate "// &
935 : "file in the SUBSYS%TOPOLOGY section.", &
936 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
937 9582 : NULLIFY (keyword)
938 : CALL keyword_create(keyword, __LOCATION__, name="UNIT", &
939 : description='Specify the unit of measurement for the coordinates in input'// &
940 : "All available CP2K units can be used.", &
941 9582 : usage="UNIT angstrom", default_c_val="angstrom")
942 9582 : CALL section_add_keyword(section, keyword)
943 9582 : CALL keyword_release(keyword)
944 :
945 : CALL keyword_create(keyword, __LOCATION__, name="SCALED", &
946 : description='Specify if the coordinates in input are scaled. '// &
947 : 'When true, the coordinates are given in multiples of the lattice vectors.', &
948 : usage="SCALED F", default_l_val=.FALSE., &
949 9582 : lone_keyword_l_val=.TRUE.)
950 9582 : CALL section_add_keyword(section, keyword)
951 9582 : CALL keyword_release(keyword)
952 :
953 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
954 : description="The atomic coordinates in the format:"//newline//newline// &
955 : "`ATOMIC_KIND X Y Z MOLNAME`"//newline//newline// &
956 : "The `MOLNAME` is optional. If not provided the molecule name "// &
957 : "is internally created. All other fields after `MOLNAME` are simply ignored.", &
958 : repeats=.TRUE., usage="{{String} {Real} {Real} {Real} {String}}", &
959 9582 : type_of_var=lchar_t)
960 9582 : CALL section_add_keyword(section, keyword)
961 9582 : CALL keyword_release(keyword)
962 9582 : END SUBROUTINE create_coord_section
963 :
964 : ! **************************************************************************************************
965 : !> \brief Creates the &SHELL_COORD section
966 : !> \param section the section to create
967 : !> \author teo
968 : ! **************************************************************************************************
969 9582 : SUBROUTINE create_shell_coord_section(section)
970 : TYPE(section_type), POINTER :: section
971 :
972 : TYPE(keyword_type), POINTER :: keyword
973 :
974 9582 : CPASSERT(.NOT. ASSOCIATED(section))
975 : CALL section_create(section, __LOCATION__, name="shell_coord", &
976 : description="The shell coordinates for the shell-model potentials"// &
977 : " xyz format with an additional column for the index of the corresponding particle", &
978 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
979 9582 : NULLIFY (keyword)
980 : CALL keyword_create(keyword, __LOCATION__, name="UNIT", &
981 : description='Specify the unit of measurement for the coordinates in input'// &
982 : "All available CP2K units can be used.", &
983 9582 : usage="UNIT angstrom", default_c_val="angstrom")
984 9582 : CALL section_add_keyword(section, keyword)
985 9582 : CALL keyword_release(keyword)
986 :
987 : CALL keyword_create(keyword, __LOCATION__, name="SCALED", &
988 : description='Specify if the coordinates in input are scaled. '// &
989 : 'When true, the coordinates are given in multiples of the lattice vectors.', &
990 : usage="SCALED F", default_l_val=.FALSE., &
991 9582 : lone_keyword_l_val=.TRUE.)
992 9582 : CALL section_add_keyword(section, keyword)
993 9582 : CALL keyword_release(keyword)
994 :
995 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
996 : description="The shell particle coordinates in the format:"//newline//newline// &
997 : "`ATOMIC_KIND X Y Z ATOMIC_INDEX`"//newline//newline// &
998 : "The `ATOMIC_INDEX` refers to the atom the shell particle belongs to.", &
999 : repeats=.TRUE., usage="{{String} {Real} {Real} {Real} {Integer}}", &
1000 9582 : type_of_var=lchar_t)
1001 9582 : CALL section_add_keyword(section, keyword)
1002 9582 : CALL keyword_release(keyword)
1003 :
1004 9582 : END SUBROUTINE create_shell_coord_section
1005 :
1006 : ! **************************************************************************************************
1007 : !> \brief Creates the &core_COORD section
1008 : !> \param section the section to create
1009 : !> \author teo
1010 : ! **************************************************************************************************
1011 9582 : SUBROUTINE create_core_coord_section(section)
1012 : TYPE(section_type), POINTER :: section
1013 :
1014 : TYPE(keyword_type), POINTER :: keyword
1015 :
1016 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1017 : CALL section_create(section, __LOCATION__, name="core_coord", &
1018 : description="The core coordinates for the shell-model potentials"// &
1019 : " xyz format with an additional column for the index of the corresponding particle", &
1020 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1021 9582 : NULLIFY (keyword)
1022 : CALL keyword_create(keyword, __LOCATION__, name="UNIT", &
1023 : description='Specify the unit of measurement for the coordinates in input'// &
1024 : "All available CP2K units can be used.", &
1025 9582 : usage="UNIT angstrom", default_c_val="angstrom")
1026 9582 : CALL section_add_keyword(section, keyword)
1027 9582 : CALL keyword_release(keyword)
1028 :
1029 : CALL keyword_create(keyword, __LOCATION__, name="SCALED", &
1030 : description='Specify if the coordinates in input are scaled. '// &
1031 : 'When true, the coordinates are given in multiples of the lattice vectors.', &
1032 : usage="SCALED F", default_l_val=.FALSE., &
1033 9582 : lone_keyword_l_val=.TRUE.)
1034 9582 : CALL section_add_keyword(section, keyword)
1035 9582 : CALL keyword_release(keyword)
1036 :
1037 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1038 : description="The core particle coordinates in the format:"//newline//newline// &
1039 : "`ATOMIC_KIND X Y Z ATOMIC_INDEX`"//newline//newline// &
1040 : "The `ATOMIC_INDEX` refers to the atom the core particle belongs to.", &
1041 : repeats=.TRUE., usage="{{String} {Real} {Real} {Real} {Integer}}", &
1042 9582 : type_of_var=lchar_t)
1043 9582 : CALL section_add_keyword(section, keyword)
1044 9582 : CALL keyword_release(keyword)
1045 :
1046 9582 : END SUBROUTINE create_core_coord_section
1047 :
1048 : ! **************************************************************************************************
1049 : !> \brief Creates the QM/MM section
1050 : !> \param section the section to create
1051 : !> \author teo
1052 : ! **************************************************************************************************
1053 9582 : SUBROUTINE create_kind_section(section)
1054 : TYPE(section_type), POINTER :: section
1055 :
1056 : TYPE(keyword_type), POINTER :: keyword
1057 : TYPE(section_type), POINTER :: subsection
1058 :
1059 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1060 :
1061 : CALL section_create(section, __LOCATION__, name="KIND", &
1062 : description="The description of the kind of the atoms (mostly for QM)", &
1063 9582 : n_keywords=20, n_subsections=1, repeats=.TRUE.)
1064 :
1065 9582 : NULLIFY (keyword)
1066 :
1067 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
1068 : description="The name of the kind described in this section.", &
1069 9582 : usage="H", default_c_val="DEFAULT")
1070 9582 : CALL section_add_keyword(section, keyword)
1071 9582 : CALL keyword_release(keyword)
1072 :
1073 : CALL keyword_create(keyword, __LOCATION__, name="BASIS_SET", &
1074 : description="The primary Gaussian basis set (NONE implies no basis used, meaningful with GHOST). "// &
1075 : "Defaults are set for TYPE {ORB} and FORM {GTO}. Possible values for TYPE are "// &
1076 : "{ORB, AUX, MIN, RI_AUX, LRI, ...}. Possible values for "// &
1077 : "FORM are {GTO, STO}. Where STO results in a GTO expansion of a Slater type basis. "// &
1078 : "If a value for FORM is given, also TYPE has to be set explicitly.", &
1079 : usage="BASIS_SET [type] [form] DZVP", type_of_var=char_t, default_c_vals=[" ", " ", " "], &
1080 : citations=[VandeVondele2005a, VandeVondele2007], &
1081 57492 : repeats=.TRUE., n_var=-1)
1082 9582 : CALL section_add_keyword(section, keyword)
1083 9582 : CALL keyword_release(keyword)
1084 :
1085 : ! old type basis set input keywords
1086 : ! kept for backward compatibility
1087 : CALL keyword_create( &
1088 : keyword, __LOCATION__, name="AUX_BASIS_SET", &
1089 : variants=s2a("AUXILIARY_BASIS_SET", "AUX_BASIS"), &
1090 : description="The auxiliary basis set (GTO type)", &
1091 : usage="AUX_BASIS_SET DZVP", default_c_val=" ", &
1092 : n_var=1, &
1093 : deprecation_notice="use 'BASIS_SET AUX ...' instead", &
1094 9582 : removed=.TRUE.)
1095 9582 : CALL section_add_keyword(section, keyword)
1096 9582 : CALL keyword_release(keyword)
1097 :
1098 : CALL keyword_create( &
1099 : keyword, __LOCATION__, name="RI_AUX_BASIS_SET", &
1100 : variants=s2a("RI_MP2_BASIS_SET", "RI_RPA_BASIS_SET", "RI_AUX_BASIS"), &
1101 : description="The RI auxiliary basis set used in WF_CORRELATION (GTO type)", &
1102 : usage="RI_AUX_BASIS_SET DZVP", default_c_val=" ", &
1103 : n_var=1, &
1104 : deprecation_notice="Use 'BASIS_SET RI_AUX ...' instead.", &
1105 9582 : removed=.TRUE.)
1106 9582 : CALL section_add_keyword(section, keyword)
1107 9582 : CALL keyword_release(keyword)
1108 :
1109 : CALL keyword_create( &
1110 : keyword, __LOCATION__, name="LRI_BASIS_SET", &
1111 : variants=s2a("LRI_BASIS"), &
1112 : description="The local resolution of identity basis set (GTO type)", &
1113 : usage="LRI_BASIS_SET", default_c_val=" ", &
1114 : n_var=1, &
1115 : deprecation_notice="Use 'BASIS_SET LRI ...' instead.", &
1116 9582 : removed=.TRUE.)
1117 9582 : CALL section_add_keyword(section, keyword)
1118 9582 : CALL keyword_release(keyword)
1119 :
1120 : CALL keyword_create( &
1121 : keyword, __LOCATION__, name="AUX_FIT_BASIS_SET", &
1122 : variants=s2a("AUXILIARY_FIT_BASIS_SET", "AUX_FIT_BASIS"), &
1123 : description="The auxiliary basis set (GTO type) for auxiliary density matrix method", &
1124 : usage="AUX_FIT_BASIS_SET DZVP", default_c_val=" ", &
1125 : citations=[Guidon2010], &
1126 : n_var=1, &
1127 : deprecation_notice="Use 'BASIS_SET AUX_FIT ...' instead.", &
1128 19164 : removed=.TRUE.)
1129 9582 : CALL section_add_keyword(section, keyword)
1130 9582 : CALL keyword_release(keyword)
1131 : ! end of old basis set keywords
1132 :
1133 : CALL keyword_create(keyword, __LOCATION__, name="ELEC_CONF", &
1134 : description="Specifies the electronic configuration used in construction the "// &
1135 : "atomic initial guess (see the pseudo potential file for the default values).", &
1136 : usage="ELEC_CONF n_elec(s) n_elec(p) n_elec(d) ... ", &
1137 9582 : n_var=-1, type_of_var=integer_t)
1138 9582 : CALL section_add_keyword(section, keyword)
1139 9582 : CALL keyword_release(keyword)
1140 :
1141 : CALL keyword_create(keyword, __LOCATION__, name="CORE_CORRECTION", &
1142 : description="Corrects the effective nuclear charge", &
1143 : usage="CORE_CORRECTION 1.0", n_var=1, &
1144 9582 : default_r_val=0.0_dp)
1145 9582 : CALL section_add_keyword(section, keyword)
1146 9582 : CALL keyword_release(keyword)
1147 :
1148 : CALL keyword_create(keyword, __LOCATION__, name="MAGNETIZATION", &
1149 : description="The magnetization used in the atomic initial guess. "// &
1150 : "Adds magnetization/2 spin-alpha electrons and removes magnetization/2 spin-beta electrons.", &
1151 : usage="MAGNETIZATION 0.5", n_var=1, &
1152 9582 : default_r_val=0.0_dp)
1153 9582 : CALL section_add_keyword(section, keyword)
1154 9582 : CALL keyword_release(keyword)
1155 :
1156 : CALL keyword_create(keyword, __LOCATION__, name="ELEMENT", &
1157 : variants=["ELEMENT_SYMBOL"], &
1158 : description="The element of the actual kind "// &
1159 : "(if not given it is inferred from the kind name)", &
1160 19164 : usage="ELEMENT O", type_of_var=char_t, n_var=1)
1161 9582 : CALL section_add_keyword(section, keyword)
1162 9582 : CALL keyword_release(keyword)
1163 :
1164 : CALL keyword_create(keyword, __LOCATION__, name="MASS", &
1165 : variants=s2a("ATOMIC_MASS", "ATOMIC_WEIGHT", "WEIGHT"), &
1166 : description="The mass of the atom "// &
1167 : "(if negative or non present it is inferred from the element symbol)", &
1168 9582 : usage="MASS 2.0", type_of_var=real_t, n_var=1)
1169 9582 : CALL section_add_keyword(section, keyword)
1170 9582 : CALL keyword_release(keyword)
1171 :
1172 : CALL keyword_create(keyword, __LOCATION__, name="POTENTIAL_FILE_NAME", &
1173 : description="The name of the file where to find this kinds pseudopotential."// &
1174 : " Default file is specified in DFT section.", &
1175 9582 : usage="POTENTIAL_FILE_NAME <PSEUDO-POTENTIAL-FILE-NAME>", default_c_val="-", n_var=1)
1176 9582 : CALL section_add_keyword(section, keyword)
1177 9582 : CALL keyword_release(keyword)
1178 :
1179 : CALL keyword_create(keyword, __LOCATION__, name="POTENTIAL_TYPE", &
1180 : description="The type of this kinds pseudopotential (ECP, ALL, GTH, UPS).", &
1181 : deprecation_notice="Use 'POTENTIAL <TYPE> ...' instead.", &
1182 9582 : usage="POTENTIAL_TYPE <TYPE>", default_c_val="", n_var=1)
1183 9582 : CALL section_add_keyword(section, keyword)
1184 9582 : CALL keyword_release(keyword)
1185 :
1186 : CALL keyword_create(keyword, __LOCATION__, name="POTENTIAL", &
1187 : variants=["POT"], &
1188 : description="The type (ECP, ALL, GTH, UPS) and name of the pseudopotential for the defined kind.", &
1189 : usage="POTENTIAL [type] <POTENTIAL-NAME>", type_of_var=char_t, default_c_vals=[" ", " "], &
1190 67074 : citations=[Goedecker1996, Hartwigsen1998, Krack2005], n_var=-1)
1191 9582 : CALL section_add_keyword(section, keyword)
1192 9582 : CALL keyword_release(keyword)
1193 :
1194 : CALL keyword_create(keyword, __LOCATION__, name="KG_POTENTIAL_FILE_NAME", &
1195 : description="The name of the file where to find this kinds KG potential."// &
1196 : " Default file is specified in DFT section.", &
1197 9582 : usage="KG_POTENTIAL_FILE_NAME <POTENTIAL-FILE-NAME>", default_c_val="-", n_var=1)
1198 9582 : CALL section_add_keyword(section, keyword)
1199 9582 : CALL keyword_release(keyword)
1200 :
1201 : CALL keyword_create(keyword, __LOCATION__, name="KG_POTENTIAL", &
1202 : variants=["KG_POT"], &
1203 : description="The name of the non-additive atomic kinetic energy potential.", &
1204 19164 : usage="KG_POTENTIAL <TNADD-POTENTIAL-NAME>", default_c_val="NONE", n_var=1)
1205 9582 : CALL section_add_keyword(section, keyword)
1206 9582 : CALL keyword_release(keyword)
1207 :
1208 : CALL keyword_create(keyword, __LOCATION__, name="ECP_SEMI_LOCAL", &
1209 : description="Use ECPs in the original semi-local form."// &
1210 : " This requires the availability of the corresponding integral library."// &
1211 : " If set to False, a fully nonlocal one-center expansion of the ECP is constructed.", &
1212 9582 : usage="ECP_SEMI_LOCAL {T,F}", default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1213 9582 : CALL section_add_keyword(section, keyword)
1214 9582 : CALL keyword_release(keyword)
1215 :
1216 : CALL keyword_create(keyword, __LOCATION__, name="COVALENT_RADIUS", &
1217 : description="Use this covalent radius (in Angstrom) for all atoms of "// &
1218 : "the atomic kind instead of the internally tabulated default value", &
1219 : usage="COVALENT_RADIUS 1.24", n_var=1, default_r_val=0.0_dp, &
1220 9582 : unit_str="angstrom")
1221 9582 : CALL section_add_keyword(section, keyword)
1222 9582 : CALL keyword_release(keyword)
1223 :
1224 : CALL keyword_create(keyword, __LOCATION__, name="VDW_RADIUS", &
1225 : description="Use this van der Waals radius (in Angstrom) for all atoms of "// &
1226 : "the atomic kind instead of the internally tabulated default value", &
1227 9582 : usage="VDW_RADIUS 1.85", n_var=1, default_r_val=0.0_dp, unit_str="angstrom")
1228 9582 : CALL section_add_keyword(section, keyword)
1229 9582 : CALL keyword_release(keyword)
1230 :
1231 : CALL keyword_create(keyword, __LOCATION__, name="HARD_EXP_RADIUS", &
1232 : description="The region where the hard density is supposed to be confined"// &
1233 : " (GAPW) (in Bohr, default is 1.2 for H and 1.512 otherwise)", &
1234 9582 : usage="HARD_EXP_RADIUS 0.9", type_of_var=real_t, n_var=1)
1235 9582 : CALL section_add_keyword(section, keyword)
1236 9582 : CALL keyword_release(keyword)
1237 :
1238 : CALL keyword_create(keyword, __LOCATION__, name="MAX_RAD_LOCAL", &
1239 : description="Max radius for the basis functions used to"// &
1240 : " generate the local projectors in GAPW [Bohr]", &
1241 9582 : usage="MAX_RAD_LOCAL 15.0", default_r_val=13.0_dp*bohr)
1242 9582 : CALL section_add_keyword(section, keyword)
1243 9582 : CALL keyword_release(keyword)
1244 :
1245 : CALL keyword_create(keyword, __LOCATION__, name="RHO0_EXP_RADIUS", &
1246 : description="the radius which defines the atomic region where "// &
1247 : "the hard compensation density is confined. "// &
1248 : "should be less than HARD_EXP_RADIUS (GAPW) (Bohr, default equals HARD_EXP_RADIUS)", &
1249 9582 : usage="RHO0_EXP_RADIUS 0.9", type_of_var=real_t, n_var=1)
1250 9582 : CALL section_add_keyword(section, keyword)
1251 9582 : CALL keyword_release(keyword)
1252 :
1253 : CALL keyword_create(keyword, __LOCATION__, name="LEBEDEV_GRID", &
1254 : description="The number of points for the angular part of "// &
1255 : "the local grid (GAPW)", &
1256 9582 : usage="LEBEDEV_GRID 40", default_i_val=50)
1257 9582 : CALL section_add_keyword(section, keyword)
1258 9582 : CALL keyword_release(keyword)
1259 :
1260 : CALL keyword_create(keyword, __LOCATION__, name="RADIAL_GRID", &
1261 : description="The number of points for the radial part of "// &
1262 : "the local grid (GAPW)", &
1263 9582 : usage="RADIAL_GRID 70", default_i_val=50)
1264 9582 : CALL section_add_keyword(section, keyword)
1265 9582 : CALL keyword_release(keyword)
1266 :
1267 : CALL keyword_create(keyword, __LOCATION__, name="MM_RADIUS", &
1268 : description="Defines the radius of the electrostatic multipole "// &
1269 : "of the atom in Fist. This radius applies to the charge, the "// &
1270 : "dipole and the quadrupole. When zero, the atom is treated as "// &
1271 : "a point multipole, otherwise it is treated as a Gaussian "// &
1272 : "charge distribution with the given radius: "// &
1273 : "p(x,y,z)*N*exp(-(x**2+y**2+z**2)/(2*MM_RADIUS**2)), where N is "// &
1274 : "a normalization constant. In the core-shell model, only the "// &
1275 : "shell is treated as a Gaussian and the core is always a point "// &
1276 : "charge.", &
1277 : usage="MM_RADIUS {real}", default_r_val=0.0_dp, type_of_var=real_t, &
1278 9582 : unit_str="angstrom", n_var=1)
1279 9582 : CALL section_add_keyword(section, keyword)
1280 9582 : CALL keyword_release(keyword)
1281 :
1282 : CALL keyword_create(keyword, __LOCATION__, name="DFTB3_PARAM", &
1283 : description="The third order parameter (derivative of hardness) used in "// &
1284 : "diagonal DFTB3 correction.", &
1285 9582 : usage="DFTB3_PARAM 0.2", default_r_val=0.0_dp)
1286 9582 : CALL section_add_keyword(section, keyword)
1287 9582 : CALL keyword_release(keyword)
1288 :
1289 : CALL keyword_create(keyword, __LOCATION__, name="LMAX_DFTB", &
1290 : description="The maximum l-quantum number of the DFTB basis for this kind.", &
1291 9582 : usage="LMAX_DFTB 1", default_i_val=-1)
1292 9582 : CALL section_add_keyword(section, keyword)
1293 9582 : CALL keyword_release(keyword)
1294 :
1295 : CALL keyword_create(keyword, __LOCATION__, name="MAO", &
1296 : description="The number of MAOs (Modified Atomic Orbitals) for this kind.", &
1297 9582 : usage="MAO 4", default_i_val=-1)
1298 9582 : CALL section_add_keyword(section, keyword)
1299 9582 : CALL keyword_release(keyword)
1300 :
1301 : ! Logicals
1302 : CALL keyword_create(keyword, __LOCATION__, name="SE_P_ORBITALS_ON_H", &
1303 : description="Forces the usage of p-orbitals on H for SEMI-EMPIRICAL calculations."// &
1304 : " This keyword applies only when the KIND is specifying an Hydrogen element."// &
1305 : " It is ignored in all other cases. ", &
1306 9582 : usage="SE_P_ORBITALS_ON_H", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1307 9582 : CALL section_add_keyword(section, keyword)
1308 9582 : CALL keyword_release(keyword)
1309 :
1310 : CALL keyword_create(keyword, __LOCATION__, name="GPW_TYPE", &
1311 : description="Force one type to be treated by the GPW scheme,"// &
1312 : " whatever are its primitives, even if the GAPW method is used", &
1313 9582 : usage="GPW_TYPE", default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1314 9582 : CALL section_add_keyword(section, keyword)
1315 9582 : CALL keyword_release(keyword)
1316 :
1317 : CALL keyword_create(keyword, __LOCATION__, &
1318 : name="GHOST", &
1319 : description="This keyword makes all atoms of this kind "// &
1320 : "ghost atoms, i.e. without pseudo or nuclear charge. "// &
1321 : "Useful to just have the basis set at that position (e.g. BSSE calculations), "// &
1322 : "or to have a non-interacting particle with BASIS_SET NONE", &
1323 : usage="GHOST", &
1324 : default_l_val=.FALSE., &
1325 9582 : lone_keyword_l_val=.TRUE.)
1326 9582 : CALL section_add_keyword(section, keyword)
1327 9582 : CALL keyword_release(keyword)
1328 :
1329 : CALL keyword_create(keyword, __LOCATION__, &
1330 : name="MONOVALENT", &
1331 : description="This keyword makes all atoms of this kind monovalent, i.e. with "// &
1332 : "a single electron and nuclear charge set to 1.0. Used to saturate dangling bonds, "// &
1333 : "ideally in conjunction with a monovalent pseudopotential. Currently GTH only.", &
1334 : usage="MONOVALENT", &
1335 : default_l_val=.FALSE., &
1336 9582 : lone_keyword_l_val=.TRUE.)
1337 9582 : CALL section_add_keyword(section, keyword)
1338 9582 : CALL keyword_release(keyword)
1339 :
1340 : CALL keyword_create(keyword, __LOCATION__, &
1341 : name="FLOATING_BASIS_CENTER", &
1342 : description="This keyword makes all atoms of this kind "// &
1343 : "floating functions, i.e. without pseudo or nuclear charge"// &
1344 : " which are subject to a geometry optimization in the outer SCF.", &
1345 : usage="FLOATING_BASIS_CENTER", &
1346 : default_l_val=.FALSE., &
1347 9582 : lone_keyword_l_val=.TRUE.)
1348 9582 : CALL section_add_keyword(section, keyword)
1349 9582 : CALL keyword_release(keyword)
1350 :
1351 : CALL keyword_create(keyword, __LOCATION__, &
1352 : name="NO_OPTIMIZE", &
1353 : description="Skip optimization of this type (used in specific basis set or"// &
1354 : " potential optimization schemes)", &
1355 : usage="NO_OPTIMIZE", &
1356 : default_l_val=.FALSE., &
1357 9582 : lone_keyword_l_val=.TRUE.)
1358 9582 : CALL section_add_keyword(section, keyword)
1359 9582 : CALL keyword_release(keyword)
1360 :
1361 : CALL keyword_create(keyword, __LOCATION__, name="PAO_BASIS_SIZE", &
1362 : description="The block size used for the polarized atomic orbital basis. "// &
1363 : "Setting PAO_BASIS_SIZE to the size of the primary basis or to a value "// &
1364 : "below one will disables the PAO method for the given atomic kind. "// &
1365 9582 : "By default PAO is disbabled.", default_i_val=0)
1366 9582 : CALL section_add_keyword(section, keyword)
1367 9582 : CALL keyword_release(keyword)
1368 :
1369 : CALL keyword_create(keyword, __LOCATION__, name="PAO_MODEL_FILE", type_of_var=lchar_t, &
1370 9582 : description="The filename of the PyTorch model for predicting PAO basis sets.")
1371 9582 : CALL section_add_keyword(section, keyword)
1372 9582 : CALL keyword_release(keyword)
1373 :
1374 9582 : NULLIFY (subsection)
1375 9582 : CALL create_pao_potential_section(subsection)
1376 9582 : CALL section_add_subsection(section, subsection)
1377 9582 : CALL section_release(subsection)
1378 :
1379 9582 : CALL create_pao_descriptor_section(subsection)
1380 9582 : CALL section_add_subsection(section, subsection)
1381 9582 : CALL section_release(subsection)
1382 :
1383 9582 : CALL create_basis_section(subsection)
1384 9582 : CALL section_add_subsection(section, subsection)
1385 9582 : CALL section_release(subsection)
1386 :
1387 9582 : CALL create_potential_section(subsection)
1388 9582 : CALL section_add_subsection(section, subsection)
1389 9582 : CALL section_release(subsection)
1390 :
1391 9582 : CALL create_kgpot_section(subsection)
1392 9582 : CALL section_add_subsection(section, subsection)
1393 9582 : CALL section_release(subsection)
1394 :
1395 9582 : CALL create_dft_plus_u_section(subsection)
1396 9582 : CALL section_add_subsection(section, subsection)
1397 9582 : CALL section_release(subsection)
1398 :
1399 9582 : CALL create_bs_section(subsection)
1400 9582 : CALL section_add_subsection(section, subsection)
1401 9582 : CALL section_release(subsection)
1402 :
1403 9582 : END SUBROUTINE create_kind_section
1404 :
1405 : ! **************************************************************************************************
1406 : !> \brief Creates the PAO_POTENTIAL section
1407 : !> \param section the section to create
1408 : !> \author Ole Schuett
1409 : ! **************************************************************************************************
1410 9582 : SUBROUTINE create_pao_potential_section(section)
1411 : TYPE(section_type), POINTER :: section
1412 :
1413 : TYPE(keyword_type), POINTER :: keyword
1414 :
1415 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1416 9582 : NULLIFY (keyword)
1417 :
1418 : CALL section_create(section, __LOCATION__, name="PAO_POTENTIAL", repeats=.TRUE., &
1419 9582 : description="Settings of the PAO potentials, which are atomic kind specific.")
1420 :
1421 : CALL keyword_create(keyword, __LOCATION__, name="MAXL", &
1422 : description="Maximum angular moment of the potential "// &
1423 9582 : "(must be an even number).", default_i_val=0)
1424 9582 : CALL section_add_keyword(section, keyword)
1425 9582 : CALL keyword_release(keyword)
1426 :
1427 : CALL keyword_create(keyword, __LOCATION__, name="BETA", &
1428 : description="Exponent of the Gaussian potential term.", &
1429 9582 : default_r_val=1.0_dp)
1430 9582 : CALL section_add_keyword(section, keyword)
1431 9582 : CALL keyword_release(keyword)
1432 :
1433 : CALL keyword_create(keyword, __LOCATION__, name="WEIGHT", &
1434 : description="Weight of Gaussian potential term.", &
1435 9582 : default_r_val=1.0_dp)
1436 9582 : CALL section_add_keyword(section, keyword)
1437 9582 : CALL keyword_release(keyword)
1438 :
1439 : CALL keyword_create(keyword, __LOCATION__, name="MAX_PROJECTOR", &
1440 : description="Maximum angular moment of the potential's projectors. "// &
1441 9582 : "Used only by the GTH parametrization", default_i_val=2)
1442 9582 : CALL section_add_keyword(section, keyword)
1443 9582 : CALL keyword_release(keyword)
1444 :
1445 9582 : END SUBROUTINE create_pao_potential_section
1446 :
1447 : ! **************************************************************************************************
1448 : !> \brief Creates the PAO_DESCRIPTOR section
1449 : !> \param section the section to create
1450 : !> \author Ole Schuett
1451 : ! **************************************************************************************************
1452 9582 : SUBROUTINE create_pao_descriptor_section(section)
1453 : TYPE(section_type), POINTER :: section
1454 :
1455 : TYPE(keyword_type), POINTER :: keyword
1456 :
1457 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1458 9582 : NULLIFY (keyword)
1459 :
1460 : CALL section_create(section, __LOCATION__, name="PAO_DESCRIPTOR", repeats=.TRUE., &
1461 9582 : description="Settings of the PAO descriptor, which are atomic kind specific.")
1462 :
1463 : CALL keyword_create(keyword, __LOCATION__, name="BETA", &
1464 : description="Exponent of the Gaussian potential term.", &
1465 9582 : default_r_val=1.0_dp)
1466 9582 : CALL section_add_keyword(section, keyword)
1467 9582 : CALL keyword_release(keyword)
1468 :
1469 : CALL keyword_create(keyword, __LOCATION__, name="SCREENING", &
1470 : description="Exponent of the Gaussian screening.", &
1471 9582 : default_r_val=0.2_dp)
1472 9582 : CALL section_add_keyword(section, keyword)
1473 9582 : CALL keyword_release(keyword)
1474 :
1475 : CALL keyword_create(keyword, __LOCATION__, name="WEIGHT", &
1476 : description="Weight of Gaussian potential term.", &
1477 9582 : default_r_val=1.0_dp)
1478 9582 : CALL section_add_keyword(section, keyword)
1479 9582 : CALL keyword_release(keyword)
1480 :
1481 9582 : END SUBROUTINE create_pao_descriptor_section
1482 :
1483 : ! **************************************************************************************************
1484 : !> \brief Create CP2K input section for BS method: imposing atomic orbital occupation
1485 : !> different from default in initialization of the density matrix
1486 : !> it works only with GUESS ATOMIC
1487 : !> \param section ...
1488 : !> \date 05.08.2009
1489 : !> \author MI
1490 : !> \version 1.0
1491 : ! **************************************************************************************************
1492 9582 : SUBROUTINE create_bs_section(section)
1493 :
1494 : TYPE(section_type), POINTER :: section
1495 :
1496 : TYPE(keyword_type), POINTER :: keyword
1497 : TYPE(section_type), POINTER :: subsection
1498 :
1499 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1500 :
1501 : CALL section_create(section, __LOCATION__, &
1502 : name="BS", &
1503 : description="Define the required atomic orbital occupation "// &
1504 : "assigned in initialization of the density matrix, by adding or "// &
1505 : "subtracting electrons from specific angular momentum channels. "// &
1506 : "It works only with GUESS ATOMIC.", &
1507 : n_keywords=0, &
1508 : n_subsections=2, &
1509 9582 : repeats=.FALSE.)
1510 :
1511 9582 : NULLIFY (keyword, subsection)
1512 :
1513 : CALL keyword_create(keyword, __LOCATION__, &
1514 : name="_SECTION_PARAMETERS_", &
1515 : description="controls the activation of the BS section", &
1516 : usage="&BS ON", &
1517 : default_l_val=.FALSE., &
1518 9582 : lone_keyword_l_val=.TRUE.)
1519 9582 : CALL section_add_keyword(section, keyword)
1520 9582 : CALL keyword_release(keyword)
1521 :
1522 : CALL section_create(subsection, __LOCATION__, name="ALPHA", description="alpha spin", &
1523 : n_keywords=3, &
1524 : n_subsections=0, &
1525 9582 : repeats=.FALSE.)
1526 :
1527 : CALL keyword_create(keyword, __LOCATION__, &
1528 : name="NEL", &
1529 : description="Orbital ccupation change per angular momentum quantum number. "// &
1530 : "In unrestricted calculations applied to spin alpha.", &
1531 : repeats=.FALSE., &
1532 : n_var=-1, &
1533 : default_i_val=-1, &
1534 9582 : usage="NEL 2")
1535 9582 : CALL section_add_keyword(subsection, keyword)
1536 9582 : CALL keyword_release(keyword)
1537 :
1538 : CALL keyword_create(keyword, __LOCATION__, &
1539 : name="L", &
1540 : variants=["L"], &
1541 : description="Angular momentum quantum number of the "// &
1542 : "orbitals whose occupation is changed", &
1543 : repeats=.FALSE., &
1544 : n_var=-1, &
1545 : default_i_val=-1, &
1546 19164 : usage="L 2")
1547 9582 : CALL section_add_keyword(subsection, keyword)
1548 9582 : CALL keyword_release(keyword)
1549 :
1550 : CALL keyword_create(keyword, __LOCATION__, &
1551 : name="N", &
1552 : variants=["N"], &
1553 : description="Principal quantum number of the "// &
1554 : "orbitals whose occupation is changed. "// &
1555 : "Default is the first not occupied", &
1556 : repeats=.FALSE., &
1557 : n_var=-1, &
1558 : default_i_val=0, &
1559 19164 : usage="N 2")
1560 9582 : CALL section_add_keyword(subsection, keyword)
1561 9582 : CALL keyword_release(keyword)
1562 9582 : CALL section_add_subsection(section, subsection)
1563 9582 : CALL section_release(subsection)
1564 :
1565 : CALL section_create(subsection, __LOCATION__, name="BETA", description="beta spin", &
1566 : n_keywords=3, &
1567 : n_subsections=0, &
1568 9582 : repeats=.FALSE.)
1569 :
1570 : CALL keyword_create(keyword, __LOCATION__, &
1571 : name="NEL", &
1572 : description="Orbital ccupation change per angular momentum quantum number. "// &
1573 : "Applied to spin beta and active only in unrestricted calculations.", &
1574 : repeats=.FALSE., &
1575 : n_var=-1, &
1576 : default_i_val=-1, &
1577 9582 : usage="NEL 2")
1578 9582 : CALL section_add_keyword(subsection, keyword)
1579 9582 : CALL keyword_release(keyword)
1580 :
1581 : CALL keyword_create(keyword, __LOCATION__, &
1582 : name="L", &
1583 : description="Angular momentum quantum number of the "// &
1584 : "orbitals of beta spin whose occupation is changed. "// &
1585 : "Active only for unrestricted calculations", &
1586 : repeats=.FALSE., &
1587 : n_var=-1, &
1588 : default_i_val=-1, &
1589 9582 : usage="L 2")
1590 9582 : CALL section_add_keyword(subsection, keyword)
1591 9582 : CALL keyword_release(keyword)
1592 :
1593 : CALL keyword_create(keyword, __LOCATION__, &
1594 : name="N", &
1595 : description="Principal quantum number of the "// &
1596 : "orbitals of beta spin whose occupation is changed. "// &
1597 : "Default is the first not occupied. "// &
1598 : "Active only for unrestricted calculations", &
1599 : repeats=.FALSE., &
1600 : n_var=-1, &
1601 : default_i_val=0, &
1602 9582 : usage="N 2")
1603 9582 : CALL section_add_keyword(subsection, keyword)
1604 9582 : CALL keyword_release(keyword)
1605 :
1606 9582 : CALL section_add_subsection(section, subsection)
1607 9582 : CALL section_release(subsection)
1608 :
1609 9582 : END SUBROUTINE create_bs_section
1610 :
1611 : ! **************************************************************************************************
1612 : !> \brief Create the topology section for FIST.. and the base is running running...
1613 : !> Contains all information regarding topology to be read in input file..
1614 : !> \param section the section to create
1615 : !> \author teo
1616 : ! **************************************************************************************************
1617 9582 : SUBROUTINE create_topology_section(section)
1618 : TYPE(section_type), POINTER :: section
1619 :
1620 : TYPE(keyword_type), POINTER :: keyword
1621 : TYPE(section_type), POINTER :: print_key, subsection
1622 :
1623 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1624 : CALL section_create(section, __LOCATION__, name="TOPOLOGY", &
1625 : description="Section specifying information regarding how to handle the topology"// &
1626 : " for classical runs.", &
1627 9582 : n_keywords=5, n_subsections=0, repeats=.FALSE.)
1628 :
1629 9582 : NULLIFY (keyword, print_key)
1630 : ! Logical
1631 : CALL keyword_create(keyword, __LOCATION__, name="USE_ELEMENT_AS_KIND", &
1632 : description="Kinds are generated according to the element name."// &
1633 : " Default=True for SE and TB methods.", &
1634 : usage="USE_ELEMENT_AS_KIND logical", &
1635 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1636 9582 : CALL section_add_keyword(section, keyword)
1637 9582 : CALL keyword_release(keyword)
1638 :
1639 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_OCCUP", &
1640 : variants=["CHARGE_O"], &
1641 : description="Read MM charges from the OCCUP field of PDB file.", &
1642 : usage="CHARGE_OCCUP logical", &
1643 19164 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1644 9582 : CALL section_add_keyword(section, keyword)
1645 9582 : CALL keyword_release(keyword)
1646 :
1647 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_BETA", &
1648 : variants=["CHARGE_B"], &
1649 : description="Read MM charges from the BETA field of PDB file.", &
1650 : usage="CHARGE_BETA logical", &
1651 19164 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1652 9582 : CALL section_add_keyword(section, keyword)
1653 9582 : CALL keyword_release(keyword)
1654 :
1655 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_EXTENDED", &
1656 : description="Read MM charges from the very last field of PDB file (starting from column 81)."// &
1657 : " No limitations of number of digits.", &
1658 : usage="CHARGE_EXTENDED logical", &
1659 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1660 9582 : CALL section_add_keyword(section, keyword)
1661 9582 : CALL keyword_release(keyword)
1662 :
1663 : CALL keyword_create(keyword, __LOCATION__, name="PARA_RES", &
1664 : description="For a protein, each residue is now considered a molecule", &
1665 : usage="PARA_RES logical", &
1666 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1667 9582 : CALL section_add_keyword(section, keyword)
1668 9582 : CALL keyword_release(keyword)
1669 :
1670 : CALL keyword_create(keyword, __LOCATION__, name="MOL_CHECK", &
1671 : description="Check molecules have the same number of atom and names.", &
1672 : usage="MOL_CHECK logical", &
1673 9582 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1674 9582 : CALL section_add_keyword(section, keyword)
1675 9582 : CALL keyword_release(keyword)
1676 :
1677 : CALL keyword_create(keyword, __LOCATION__, name="USE_G96_VELOCITY", &
1678 : description="Use the velocities in the G96 coordinate files as the starting velocity", &
1679 : usage="USE_G96_VELOCITY logical", &
1680 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1681 9582 : CALL section_add_keyword(section, keyword)
1682 9582 : CALL keyword_release(keyword)
1683 :
1684 : ! Character
1685 : CALL keyword_create(keyword, __LOCATION__, name="COORD_FILE_NAME", &
1686 : variants=s2a("COORD_FILE"), &
1687 : description="Specifies the filename that contains coordinates. "// &
1688 : "In case the CELL section is not set explicitly but this file "// &
1689 : "contains cell information, including CIF, PDB and (Extended) XYZ "// &
1690 : "formats, this file is also parsed for setting up the simulation cell.", &
1691 9582 : usage="COORD_FILE_NAME <FILENAME>", type_of_var=lchar_t)
1692 9582 : CALL section_add_keyword(section, keyword)
1693 9582 : CALL keyword_release(keyword)
1694 :
1695 : CALL keyword_create(keyword, __LOCATION__, name="COORD_FILE_FORMAT", &
1696 : variants=s2a("COORDINATE"), &
1697 : description="Set up the way in which coordinates will be read.", &
1698 : usage="COORD_FILE_FORMAT (OFF|PDB|XYZ|G96|CRD|CIF|XTL|CP2K)", &
1699 : enum_c_vals=s2a("OFF", "PDB", "XYZ", "G96", "CRD", "CIF", "XTL", "CP2K"), &
1700 : enum_i_vals=[do_coord_off, do_coord_pdb, do_coord_xyz, do_coord_g96, do_coord_crd, &
1701 : do_coord_cif, do_coord_xtl, do_coord_cp2k], &
1702 : enum_desc=s2a( &
1703 : "Coordinates read in the &COORD section of the input file", &
1704 : "Coordinates provided through a PDB file format", &
1705 : "Coordinates provided through an XYZ file format", &
1706 : "Coordinates provided through a GROMOS96 file format", &
1707 : "Coordinates provided through an AMBER file format", &
1708 : "Coordinates provided through a CIF (Crystallographic Information File) file format", &
1709 : "Coordinates provided through a XTL (MSI native) file format", &
1710 : "Read the coordinates in CP2K &COORD section format from an external file. "// &
1711 : "NOTE: This file will be overwritten with the latest coordinates."), &
1712 9582 : default_i_val=do_coord_off)
1713 9582 : CALL section_add_keyword(section, keyword)
1714 9582 : CALL keyword_release(keyword)
1715 :
1716 : CALL keyword_create(keyword, __LOCATION__, name="NUMBER_OF_ATOMS", &
1717 : variants=s2a("NATOMS", "NATOM"), &
1718 : description="Optionally define the number of atoms read from an external file "// &
1719 : "(see COORD_FILE_NAME) if the COORD_FILE_FORMAT CP2K is used", &
1720 : repeats=.FALSE., &
1721 : n_var=1, &
1722 : type_of_var=integer_t, &
1723 : default_i_val=-1, &
1724 9582 : usage="NATOMS 768000")
1725 9582 : CALL section_add_keyword(section, keyword)
1726 9582 : CALL keyword_release(keyword)
1727 :
1728 9582 : CALL connectivity_framework(section, do_conn_generate)
1729 :
1730 : CALL keyword_create(keyword, __LOCATION__, name="DISABLE_EXCLUSION_LISTS", &
1731 : description="Do not build any exclusion lists.", &
1732 : usage="DISABLE_EXCLUSION_LISTS", &
1733 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1734 9582 : CALL section_add_keyword(section, keyword)
1735 9582 : CALL keyword_release(keyword)
1736 :
1737 : CALL keyword_create(keyword, __LOCATION__, name="EXCLUDE_VDW", &
1738 : description="Specifies which kind of Van der Waals interaction to skip.", &
1739 : usage="EXCLUDE_VDW (1-1||1-2||1-3||1-4)", &
1740 : enum_c_vals=s2a("1-1", "1-2", "1-3", "1-4"), &
1741 : enum_i_vals=[do_skip_11, do_skip_12, do_skip_13, do_skip_14], &
1742 9582 : default_i_val=do_skip_13)
1743 9582 : CALL section_add_keyword(section, keyword)
1744 9582 : CALL keyword_release(keyword)
1745 :
1746 : CALL keyword_create(keyword, __LOCATION__, name="EXCLUDE_EI", &
1747 : description="Specifies which kind of Electrostatic interaction to skip.", &
1748 : usage="EXCLUDE_EI (1-1||1-2||1-3||1-4)", &
1749 : enum_c_vals=s2a("1-1", "1-2", "1-3", "1-4"), &
1750 : enum_i_vals=[do_skip_11, do_skip_12, do_skip_13, do_skip_14], &
1751 9582 : default_i_val=do_skip_13)
1752 9582 : CALL section_add_keyword(section, keyword)
1753 9582 : CALL keyword_release(keyword)
1754 :
1755 : CALL keyword_create(keyword, __LOCATION__, name="AUTOGEN_EXCLUDE_LISTS", &
1756 : description="When True, the exclude lists are solely based on"// &
1757 : " the bond data in the topology. The (minimal)"// &
1758 : " number of bonds between two atoms is used to"// &
1759 : " determine if the atom pair is added to an"// &
1760 : " exclusion list. When False, 1-2 exclusion is based"// &
1761 : " on bonds in the topology, 1-3 exclusion is based"// &
1762 : " on bonds and bends in the topology, 1-4 exclusion"// &
1763 : " is based on bonds, bends and dihedrals in the"// &
1764 : " topology. This implies that a missing dihedral in"// &
1765 : " the topology will cause the corresponding 1-4 pair"// &
1766 : " not to be in the exclusion list, in case 1-4"// &
1767 : " exclusion is requested for VDW or EI interactions.", &
1768 : usage="AUTOGEN_EXCLUDE_LISTS logical", &
1769 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1770 9582 : CALL section_add_keyword(section, keyword)
1771 9582 : CALL keyword_release(keyword)
1772 :
1773 : CALL keyword_create( &
1774 : keyword, __LOCATION__, name="MULTIPLE_UNIT_CELL", &
1775 : description="Specifies the numbers of repetition in space (X, Y, Z) of the defined cell, "// &
1776 : "assuming it as a unit cell. This keyword affects only the coordinates specification. The same keyword "// &
1777 : "in SUBSYS%CELL%MULTIPLE_UNIT_CELL should be modified in order to affect the cell "// &
1778 : "specification.", usage="MULTIPLE_UNIT_CELL 1 1 1", &
1779 9582 : n_var=3, default_i_vals=[1, 1, 1], repeats=.FALSE.)
1780 9582 : CALL section_add_keyword(section, keyword)
1781 9582 : CALL keyword_release(keyword)
1782 :
1783 : CALL keyword_create(keyword, __LOCATION__, name="MEMORY_PROGRESSION_FACTOR", &
1784 : description="This keyword is quite technical and should normally not be changed by the user. It "// &
1785 : "affects the memory allocation during the construction of the topology. It does NOT affect the "// &
1786 : "memory used once the topology is built.", &
1787 9582 : n_var=1, default_r_val=1.2_dp, repeats=.FALSE.)
1788 9582 : CALL section_add_keyword(section, keyword)
1789 9582 : CALL keyword_release(keyword)
1790 :
1791 : CALL cp_print_key_section_create(print_key, __LOCATION__, "DUMP_PDB", &
1792 : description="controls the dumping of the PDB at the starting geometry", &
1793 9582 : print_level=debug_print_level, filename="dump")
1794 9582 : CALL section_add_subsection(section, print_key)
1795 :
1796 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_OCCUP", &
1797 : variants=["CHARGE_O"], &
1798 : description="Write the MM charges to the OCCUP field of the PDB file", &
1799 : usage="CHARGE_OCCUP logical", &
1800 19164 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1801 9582 : CALL section_add_keyword(print_key, keyword)
1802 9582 : CALL keyword_release(keyword)
1803 :
1804 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_BETA", &
1805 : variants=["CHARGE_B"], &
1806 : description="Write the MM charges to the BETA field of the PDB file", &
1807 : usage="CHARGE_BETA logical", &
1808 19164 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1809 9582 : CALL section_add_keyword(print_key, keyword)
1810 9582 : CALL keyword_release(keyword)
1811 :
1812 : CALL keyword_create(keyword, __LOCATION__, name="CHARGE_EXTENDED", &
1813 : description="Write the MM charges to the very last field of the PDB file (starting from column 81)", &
1814 : usage="CHARGE_EXTENDED logical", &
1815 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1816 9582 : CALL section_add_keyword(print_key, keyword)
1817 9582 : CALL keyword_release(keyword)
1818 :
1819 9582 : CALL section_release(print_key)
1820 :
1821 : CALL cp_print_key_section_create(print_key, __LOCATION__, "DUMP_PSF", &
1822 : description="controls the dumping of the PSF connectivity", &
1823 9582 : print_level=debug_print_level, filename="dump")
1824 9582 : CALL section_add_subsection(section, print_key)
1825 9582 : CALL section_release(print_key)
1826 :
1827 9582 : NULLIFY (subsection)
1828 9582 : CALL create_exclude_list_section(subsection, "EXCLUDE_VDW_LIST")
1829 9582 : CALL section_add_subsection(section, subsection)
1830 9582 : CALL section_release(subsection)
1831 :
1832 9582 : CALL create_exclude_list_section(subsection, "EXCLUDE_EI_LIST")
1833 9582 : CALL section_add_subsection(section, subsection)
1834 9582 : CALL section_release(subsection)
1835 :
1836 9582 : CALL create_center_section(subsection)
1837 9582 : CALL section_add_subsection(section, subsection)
1838 9582 : CALL section_release(subsection)
1839 :
1840 9582 : CALL create_generate_section(subsection)
1841 9582 : CALL section_add_subsection(section, subsection)
1842 9582 : CALL section_release(subsection)
1843 :
1844 9582 : CALL create_molset_section(subsection)
1845 9582 : CALL section_add_subsection(section, subsection)
1846 9582 : CALL section_release(subsection)
1847 :
1848 9582 : END SUBROUTINE create_topology_section
1849 :
1850 : ! **************************************************************************************************
1851 : !> \brief Setup a list of fine exclusion elements
1852 : !> \param section the section to create
1853 : !> \param header ...
1854 : !> \author Teodoro Laino [tlaino] - 12.2009
1855 : ! **************************************************************************************************
1856 19164 : SUBROUTINE create_exclude_list_section(section, header)
1857 : TYPE(section_type), POINTER :: section
1858 : CHARACTER(LEN=*), INTENT(IN) :: header
1859 :
1860 : TYPE(keyword_type), POINTER :: keyword
1861 :
1862 19164 : CPASSERT(.NOT. ASSOCIATED(section))
1863 19164 : NULLIFY (keyword)
1864 : CALL section_create(section, __LOCATION__, TRIM(header), &
1865 : description="Speficy bonds (via atom kinds) for fine tuning of 1-2 "// &
1866 : "exclusion lists. If this section is not present the 1-2 exclusion is "// &
1867 : "applied to all bond kinds. When this section is present the 1-2 exclusion "// &
1868 : "is applied ONLY to the bonds defined herein. This section allows ONLY fine tuning of 1-2 "// &
1869 : "interactions. ", &
1870 19164 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1871 :
1872 : CALL keyword_create(keyword, __LOCATION__, name="BOND", &
1873 : description="Specify the atom kinds involved in the bond for which 1-2 exclusion holds.", &
1874 : usage="BOND {KIND1} {KIND2}", type_of_var=char_t, &
1875 19164 : n_var=2)
1876 19164 : CALL section_add_keyword(section, keyword)
1877 19164 : CALL keyword_release(keyword)
1878 19164 : END SUBROUTINE create_exclude_list_section
1879 :
1880 : ! **************************************************************************************************
1881 : !> \brief Specify keywords used to center molecule in the box
1882 : !> \param section the section to create
1883 : !> \author Teodoro Laino [tlaino] - University of Zurich - 06.2009
1884 : ! **************************************************************************************************
1885 9582 : SUBROUTINE create_center_section(section)
1886 : TYPE(section_type), POINTER :: section
1887 :
1888 : TYPE(keyword_type), POINTER :: keyword
1889 :
1890 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1891 9582 : NULLIFY (keyword)
1892 : CALL section_create(section, __LOCATION__, "CENTER_COORDINATES", &
1893 : description="Allows centering the coordinates of the system in the box. "// &
1894 : "The centering point can be defined by the user.", &
1895 9582 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1896 :
1897 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
1898 : description="Controls the activation of the centering method", &
1899 : usage="&CENTER_COORDINATES T", &
1900 : default_l_val=.FALSE., &
1901 9582 : lone_keyword_l_val=.TRUE.)
1902 9582 : CALL section_add_keyword(section, keyword)
1903 9582 : CALL keyword_release(keyword)
1904 :
1905 : CALL keyword_create(keyword, __LOCATION__, name="CENTER_POINT", &
1906 : description="Specify the point used for centering the coordinates. Default is to "// &
1907 : "center the system in cell/2. ", type_of_var=real_t, n_var=3, &
1908 9582 : repeats=.FALSE.)
1909 9582 : CALL section_add_keyword(section, keyword)
1910 9582 : CALL keyword_release(keyword)
1911 9582 : END SUBROUTINE create_center_section
1912 :
1913 : ! **************************************************************************************************
1914 : !> \brief Specify keywords used to setup several molecules with few connectivity files
1915 : !> \param section the section to create
1916 : !> \author Teodoro Laino [tlaino] - University of Zurich - 08.2008
1917 : ! **************************************************************************************************
1918 9582 : SUBROUTINE create_molset_section(section)
1919 : TYPE(section_type), POINTER :: section
1920 :
1921 : TYPE(keyword_type), POINTER :: keyword
1922 : TYPE(section_type), POINTER :: subsection, subsubsection
1923 :
1924 9582 : CPASSERT(.NOT. ASSOCIATED(section))
1925 9582 : NULLIFY (keyword, subsection, subsubsection)
1926 : CALL section_create(section, __LOCATION__, name="MOL_SET", &
1927 : description="Specify the connectivity of a full system specifying the connectivity"// &
1928 : " of the fragments of the system.", &
1929 9582 : n_keywords=2, n_subsections=0, repeats=.FALSE.)
1930 :
1931 : ! MOLECULES
1932 : CALL section_create(subsection, __LOCATION__, name="MOLECULE", &
1933 : description="Specify information about the connectivity of single molecules", &
1934 9582 : n_keywords=2, n_subsections=0, repeats=.TRUE.)
1935 :
1936 : CALL keyword_create(keyword, __LOCATION__, name="NMOL", &
1937 : description="number of molecules ", &
1938 9582 : usage="NMOL {integer}", default_i_val=1)
1939 9582 : CALL section_add_keyword(subsection, keyword)
1940 9582 : CALL keyword_release(keyword)
1941 :
1942 9582 : CALL connectivity_framework(subsection, do_conn_psf)
1943 9582 : CALL section_add_subsection(section, subsection)
1944 9582 : CALL section_release(subsection)
1945 :
1946 : ! MERGE MOLECULES
1947 : CALL section_create(subsection, __LOCATION__, name="MERGE_MOLECULES", &
1948 : description="Enables the creation of connecting bridges (bonds, angles, torsions, impropers)"// &
1949 : " between the two or more molecules defined with independent connectivity.", &
1950 9582 : n_keywords=2, n_subsections=0, repeats=.FALSE.)
1951 :
1952 : CALL section_create(subsubsection, __LOCATION__, name="bonds", &
1953 9582 : description="Defines new bonds", n_keywords=2, n_subsections=0, repeats=.FALSE.)
1954 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1955 : description="Two integer indexes per line defining the new bond."// &
1956 : " Indexes must be relative to the full system and not to the single molecules", &
1957 : repeats=.TRUE., &
1958 9582 : usage="{Integer} {Integer}", type_of_var=integer_t, n_var=2)
1959 9582 : CALL section_add_keyword(subsubsection, keyword)
1960 9582 : CALL keyword_release(keyword)
1961 9582 : CALL section_add_subsection(subsection, subsubsection)
1962 9582 : CALL section_release(subsubsection)
1963 :
1964 : CALL section_create(subsubsection, __LOCATION__, name="angles", &
1965 : description="Defines new angles", n_keywords=2, n_subsections=0, &
1966 9582 : repeats=.FALSE.)
1967 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1968 : description="Three integer indexes per line defining the new angle"// &
1969 : " Indexes must be relative to the full system and not to the single molecules", repeats=.TRUE., &
1970 9582 : usage="{Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=3)
1971 9582 : CALL section_add_keyword(subsubsection, keyword)
1972 9582 : CALL keyword_release(keyword)
1973 9582 : CALL section_add_subsection(subsection, subsubsection)
1974 9582 : CALL section_release(subsubsection)
1975 :
1976 : CALL section_create(subsubsection, __LOCATION__, name="torsions", &
1977 : description="Defines new torsions", n_keywords=2, n_subsections=0, &
1978 9582 : repeats=.FALSE.)
1979 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1980 : description="Four integer indexes per line defining the new torsion"// &
1981 : " Indexes must be relative to the full system and not to the single molecules", repeats=.TRUE., &
1982 9582 : usage="{Integer} {Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=4)
1983 9582 : CALL section_add_keyword(subsubsection, keyword)
1984 9582 : CALL keyword_release(keyword)
1985 9582 : CALL section_add_subsection(subsection, subsubsection)
1986 9582 : CALL section_release(subsubsection)
1987 :
1988 : CALL section_create(subsubsection, __LOCATION__, name="impropers", &
1989 : description="Defines new impropers", n_keywords=2, n_subsections=0, &
1990 9582 : repeats=.FALSE.)
1991 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
1992 : description="Four integer indexes per line defining the new improper"// &
1993 : " Indexes must be relative to the full system and not to the single molecules", repeats=.TRUE., &
1994 9582 : usage="{Integer} {Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=4)
1995 9582 : CALL section_add_keyword(subsubsection, keyword)
1996 9582 : CALL keyword_release(keyword)
1997 9582 : CALL section_add_subsection(subsection, subsubsection)
1998 9582 : CALL section_release(subsubsection)
1999 :
2000 9582 : CALL section_add_subsection(section, subsection)
2001 9582 : CALL section_release(subsection)
2002 :
2003 9582 : END SUBROUTINE create_molset_section
2004 :
2005 : ! **************************************************************************************************
2006 : !> \brief Specify keywords used to generate connectivity
2007 : !> \param section the section to create
2008 : !> \author Teodoro Laino [tlaino] - University of Zurich - 08.2008
2009 : ! **************************************************************************************************
2010 9582 : SUBROUTINE create_generate_section(section)
2011 : TYPE(section_type), POINTER :: section
2012 :
2013 : TYPE(keyword_type), POINTER :: keyword
2014 : TYPE(section_type), POINTER :: subsection
2015 :
2016 9582 : CPASSERT(.NOT. ASSOCIATED(section))
2017 9582 : NULLIFY (keyword, subsection)
2018 : CALL section_create(section, __LOCATION__, name="GENERATE", &
2019 : description="Setup of keywords controlling the generation of the connectivity", &
2020 9582 : n_keywords=2, n_subsections=0, repeats=.TRUE.)
2021 :
2022 : CALL keyword_create(keyword, __LOCATION__, name="REORDER", &
2023 : description="Reorder a list of atomic coordinates into order so it can be packed correctly.", &
2024 : usage="REORDER <LOGICAL>", &
2025 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
2026 9582 : CALL section_add_keyword(section, keyword)
2027 9582 : CALL keyword_release(keyword)
2028 :
2029 : CALL keyword_create(keyword, __LOCATION__, name="CREATE_MOLECULES", &
2030 : description="Create molecules names and definition. Can be used to override the"// &
2031 : " molecules specifications of a possible input connectivity or to create molecules"// &
2032 : " specifications for file types as XYZ, missing of molecules definitions.", &
2033 : usage="CREATE_MOLECULES <LOGICAL>", &
2034 9582 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
2035 9582 : CALL section_add_keyword(section, keyword)
2036 9582 : CALL keyword_release(keyword)
2037 :
2038 : CALL keyword_create(keyword, __LOCATION__, name="BONDPARM", &
2039 : description="Used in conjunction with BONDPARM_FACTOR to "// &
2040 : "help determine wheather there is bonding "// &
2041 : "between two atoms based on a distance criteria. "// &
2042 : "Can use covalent radii information or VDW radii information", &
2043 : usage="BONDPARM (COVALENT||VDW)", &
2044 : enum_c_vals=s2a("COVALENT", "VDW"), &
2045 : enum_i_vals=[do_bondparm_covalent, do_bondparm_vdw], &
2046 9582 : default_i_val=do_bondparm_covalent)
2047 9582 : CALL section_add_keyword(section, keyword)
2048 9582 : CALL keyword_release(keyword)
2049 :
2050 : CALL keyword_create(keyword, __LOCATION__, name="BONDPARM_FACTOR", &
2051 : description="Used in conjunction with BONDPARM to help "// &
2052 : "determine wheather there is bonding between "// &
2053 : "two atoms based on a distance criteria.", &
2054 9582 : usage="bondparm_factor {real}", default_r_val=1.1_dp)
2055 9582 : CALL section_add_keyword(section, keyword)
2056 9582 : CALL keyword_release(keyword)
2057 :
2058 : CALL keyword_create(keyword, __LOCATION__, name="BONDLENGTH_MAX", &
2059 : description="Maximum distance to generate neighbor lists to build connectivity", &
2060 : usage="BONDLENGTH_MAX <real>", &
2061 : default_r_val=cp_unit_to_cp2k(value=3.0_dp, unit_str="angstrom"), &
2062 9582 : unit_str="angstrom")
2063 9582 : CALL section_add_keyword(section, keyword)
2064 9582 : CALL keyword_release(keyword)
2065 :
2066 : CALL keyword_create(keyword, __LOCATION__, name="BONDLENGTH_MIN", &
2067 : description="Minimum distance to generate neighbor lists to build connectivity", &
2068 : usage="BONDLENGTH_MIN <real>", &
2069 : default_r_val=cp_unit_to_cp2k(value=0.01_dp, unit_str="angstrom"), &
2070 9582 : unit_str="angstrom")
2071 9582 : CALL section_add_keyword(section, keyword)
2072 9582 : CALL keyword_release(keyword)
2073 :
2074 : ! BONDS
2075 : CALL section_create(subsection, __LOCATION__, name="BOND", &
2076 : description="Section used to add/remove bonds in the connectivity."// &
2077 : " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2078 9582 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
2079 :
2080 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
2081 : description="controls the activation of the bond", &
2082 : usage="&BOND (ADD|REMOVE)", &
2083 : enum_c_vals=s2a("ADD", "REMOVE"), &
2084 : enum_i_vals=[do_add, do_remove], &
2085 9582 : default_i_val=do_add)
2086 9582 : CALL section_add_keyword(subsection, keyword)
2087 9582 : CALL keyword_release(keyword)
2088 :
2089 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2090 : description="Specifies two atomic index united by a covalent bond", &
2091 : usage="ATOMS {integer} {integer}", type_of_var=integer_t, n_var=2, &
2092 9582 : repeats=.TRUE.)
2093 9582 : CALL section_add_keyword(subsection, keyword)
2094 9582 : CALL keyword_release(keyword)
2095 :
2096 9582 : CALL section_add_subsection(section, subsection)
2097 9582 : CALL section_release(subsection)
2098 :
2099 : ! ANGLES
2100 : CALL section_create(subsection, __LOCATION__, name="ANGLE", &
2101 : description="Section used to add/remove angles in the connectivity."// &
2102 : " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2103 9582 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
2104 :
2105 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
2106 : description="controls the activation of the bond", &
2107 : usage="&ANGLE (ADD|REMOVE)", &
2108 : enum_c_vals=s2a("ADD", "REMOVE"), &
2109 : enum_i_vals=[do_add, do_remove], &
2110 9582 : default_i_val=do_add)
2111 9582 : CALL section_add_keyword(subsection, keyword)
2112 9582 : CALL keyword_release(keyword)
2113 :
2114 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2115 : description="Specifies two atomic index united by a covalent bond", &
2116 : usage="ATOMS {integer} {integer} {integer} ", type_of_var=integer_t, n_var=3, &
2117 9582 : repeats=.TRUE.)
2118 9582 : CALL section_add_keyword(subsection, keyword)
2119 9582 : CALL keyword_release(keyword)
2120 :
2121 9582 : CALL section_add_subsection(section, subsection)
2122 9582 : CALL section_release(subsection)
2123 :
2124 : ! TORSIONS
2125 : CALL section_create(subsection, __LOCATION__, name="TORSION", &
2126 : description="Section used to add/remove torsion in the connectivity."// &
2127 : " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2128 9582 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
2129 :
2130 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
2131 : description="controls the activation of the bond", &
2132 : usage="&TORSION (ADD|REMOVE)", &
2133 : enum_c_vals=s2a("ADD", "REMOVE"), &
2134 : enum_i_vals=[do_add, do_remove], &
2135 9582 : default_i_val=do_add)
2136 9582 : CALL section_add_keyword(subsection, keyword)
2137 9582 : CALL keyword_release(keyword)
2138 :
2139 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2140 : description="Specifies two atomic index united by a covalent bond", &
2141 : usage="ATOMS {integer} {integer} {integer} {integer} ", type_of_var=integer_t, n_var=4, &
2142 9582 : repeats=.TRUE.)
2143 9582 : CALL section_add_keyword(subsection, keyword)
2144 9582 : CALL keyword_release(keyword)
2145 :
2146 9582 : CALL section_add_subsection(section, subsection)
2147 9582 : CALL section_release(subsection)
2148 :
2149 : ! IMPROPERS
2150 : CALL section_create(subsection, __LOCATION__, name="IMPROPER", &
2151 : description="Section used to add/remove improper in the connectivity."// &
2152 : " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2153 9582 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
2154 :
2155 : CALL keyword_create(keyword, __LOCATION__, name="_SECTION_PARAMETERS_", &
2156 : description="controls the activation of the bond", &
2157 : usage="&IMPROPER (ADD|REMOVE)", &
2158 : enum_c_vals=s2a("ADD", "REMOVE"), &
2159 : enum_i_vals=[do_add, do_remove], &
2160 9582 : default_i_val=do_add)
2161 9582 : CALL section_add_keyword(subsection, keyword)
2162 9582 : CALL keyword_release(keyword)
2163 :
2164 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2165 : description="Specifies two atomic index united by a covalent bond", &
2166 : usage="ATOMS {integer} {integer} {integer} {integer} ", type_of_var=integer_t, n_var=4, &
2167 9582 : repeats=.TRUE.)
2168 9582 : CALL section_add_keyword(subsection, keyword)
2169 9582 : CALL keyword_release(keyword)
2170 :
2171 9582 : CALL section_add_subsection(section, subsection)
2172 9582 : CALL section_release(subsection)
2173 :
2174 : ! ISOLATED ATOMS
2175 : CALL section_create(subsection, __LOCATION__, name="ISOLATED_ATOMS", &
2176 : description=" This section specifies the atoms that one considers isolated. Useful when present "// &
2177 9582 : "ions in solution.", n_keywords=1, n_subsections=0, repeats=.FALSE.)
2178 : CALL keyword_create(keyword, __LOCATION__, name="LIST", &
2179 : description="Specifies a list of atomic indexes of the isolated ion", &
2180 : usage="LIST {integer}", type_of_var=integer_t, n_var=-1, &
2181 9582 : repeats=.TRUE.)
2182 9582 : CALL section_add_keyword(subsection, keyword)
2183 9582 : CALL keyword_release(keyword)
2184 :
2185 9582 : CALL section_add_subsection(section, subsection)
2186 9582 : CALL section_release(subsection)
2187 :
2188 : ! Neighbor lists keys and printing handling the construction of NL for the connectivity
2189 9582 : CALL create_neighbor_lists_section(subsection)
2190 9582 : CALL section_add_subsection(section, subsection)
2191 9582 : CALL section_release(subsection)
2192 :
2193 9582 : CALL create_gen_print_section(subsection)
2194 9582 : CALL section_add_subsection(section, subsection)
2195 9582 : CALL section_release(subsection)
2196 :
2197 9582 : END SUBROUTINE create_generate_section
2198 :
2199 : ! **************************************************************************************************
2200 : !> \brief Create the print gen section
2201 : !> \param section the section to create
2202 : !> \author teo
2203 : ! **************************************************************************************************
2204 9582 : SUBROUTINE create_gen_print_section(section)
2205 : TYPE(section_type), POINTER :: section
2206 :
2207 : TYPE(section_type), POINTER :: print_key
2208 :
2209 9582 : CPASSERT(.NOT. ASSOCIATED(section))
2210 : CALL section_create(section, __LOCATION__, name="print", &
2211 : description="Section of possible print options in GENERATE code.", &
2212 9582 : n_keywords=0, n_subsections=1, repeats=.FALSE.)
2213 :
2214 9582 : NULLIFY (print_key)
2215 : CALL cp_print_key_section_create(print_key, __LOCATION__, "NEIGHBOR_LISTS", &
2216 : description="Activates the printing of the neighbor lists used"// &
2217 : " for generating the connectivity.", print_level=high_print_level, &
2218 9582 : filename="", unit_str="angstrom")
2219 9582 : CALL section_add_subsection(section, print_key)
2220 9582 : CALL section_release(print_key)
2221 :
2222 : CALL cp_print_key_section_create(print_key, __LOCATION__, "SUBCELL", &
2223 : description="Activates the printing of the subcells used for the "// &
2224 : "generation of neighbor lists for connectivity.", &
2225 9582 : print_level=high_print_level, filename="__STD_OUT__")
2226 9582 : CALL section_add_subsection(section, print_key)
2227 9582 : CALL section_release(print_key)
2228 :
2229 9582 : END SUBROUTINE create_gen_print_section
2230 :
2231 : ! **************************************************************************************************
2232 : !> \brief Specify keywords used to define connectivity
2233 : !> \param section the section to create
2234 : !> \param default ...
2235 : !> \author teo
2236 : ! **************************************************************************************************
2237 19164 : SUBROUTINE connectivity_framework(section, default)
2238 : TYPE(section_type), POINTER :: section
2239 : INTEGER, INTENT(IN) :: default
2240 :
2241 : TYPE(keyword_type), POINTER :: keyword
2242 :
2243 19164 : CPASSERT(ASSOCIATED(section))
2244 19164 : NULLIFY (keyword)
2245 : CALL keyword_create(keyword, __LOCATION__, name="CONN_FILE_NAME", &
2246 : variants=["CONN_FILE"], &
2247 : description="Specifies the filename that contains the molecular connectivity.", &
2248 38328 : usage="CONN_FILE_NAME <FILENAME>", type_of_var=lchar_t)
2249 19164 : CALL section_add_keyword(section, keyword)
2250 19164 : CALL keyword_release(keyword)
2251 :
2252 : CALL keyword_create(keyword, __LOCATION__, name="CONN_FILE_FORMAT", &
2253 : variants=["CONNECTIVITY"], &
2254 : description="Ways to determine and generate a molecules. "// &
2255 : "Default is to use GENERATE", &
2256 : usage="CONN_FILE_FORMAT (PSF|UPSF|MOL_SET|GENERATE|OFF|G87|G96|AMBER|USER)", &
2257 : enum_c_vals=s2a("PSF", "UPSF", "MOL_SET", "GENERATE", "OFF", "G87", "G96", "AMBER", "USER"), &
2258 : enum_i_vals=[do_conn_psf, &
2259 : do_conn_psf_u, &
2260 : do_conn_mol_set, &
2261 : do_conn_generate, &
2262 : do_conn_off, &
2263 : do_conn_g87, &
2264 : do_conn_g96, &
2265 : do_conn_amb7, &
2266 : do_conn_user], &
2267 : enum_desc=s2a("Use a PSF file to determine the connectivity."// &
2268 : " (support standard CHARMM/XPLOR and EXT CHARMM)", &
2269 : "Read a PSF file in an unformatted way (useful for not so standard PSF).", &
2270 : "Use multiple PSF (for now...) files to generate the whole system.", &
2271 : "Use a simple distance criteria. (Look at keyword BONDPARM)", &
2272 : "Do not generate molecules. (e.g. for QS or ill defined systems)", &
2273 : "Use GROMOS G87 topology file.", &
2274 : "Use GROMOS G96 topology file.", &
2275 : "Use AMBER topology file for reading connectivity (compatible starting from AMBER V.7)", &
2276 : "Allows the definition of molecules and residues based on the 5th and 6th column of "// &
2277 : "the COORD section. This option can be handy for the definition of molecules with QS "// &
2278 : "or to save memory in the case of very large systems (use PARA_RES off)."), &
2279 38328 : default_i_val=default)
2280 19164 : CALL section_add_keyword(section, keyword)
2281 19164 : CALL keyword_release(keyword)
2282 19164 : END SUBROUTINE connectivity_framework
2283 :
2284 : ! **************************************************************************************************
2285 : !> \brief Create CP2K input section for the DFT+U method parameters
2286 : !> \param section ...
2287 : !> \date 01.11.2007
2288 : !> \author Matthias Krack (MK)
2289 : !> \version 1.0
2290 : ! **************************************************************************************************
2291 9582 : SUBROUTINE create_dft_plus_u_section(section)
2292 :
2293 : TYPE(section_type), POINTER :: section
2294 :
2295 : TYPE(keyword_type), POINTER :: keyword
2296 : TYPE(section_type), POINTER :: subsection
2297 :
2298 9582 : CPASSERT(.NOT. ASSOCIATED(section))
2299 :
2300 : CALL section_create(section, __LOCATION__, &
2301 : name="DFT_PLUS_U", &
2302 : description="Define the parameters for a DFT+U run", &
2303 : n_keywords=3, &
2304 : n_subsections=1, &
2305 9582 : repeats=.FALSE.)
2306 9582 : NULLIFY (keyword)
2307 :
2308 : CALL keyword_create(keyword, __LOCATION__, &
2309 : name="_SECTION_PARAMETERS_", &
2310 : description="Controls the activation of the DFT+U section", &
2311 : usage="&DFT_PLUS_U ON", &
2312 : default_l_val=.FALSE., &
2313 9582 : lone_keyword_l_val=.TRUE.)
2314 9582 : CALL section_add_keyword(section, keyword)
2315 9582 : CALL keyword_release(keyword)
2316 :
2317 : CALL keyword_create(keyword, __LOCATION__, &
2318 : name="L", &
2319 : description="Angular momentum quantum number of the "// &
2320 : "orbitals to which the correction is applied", &
2321 : repeats=.FALSE., &
2322 : n_var=1, &
2323 : type_of_var=integer_t, &
2324 : default_i_val=-1, &
2325 9582 : usage="L 2")
2326 9582 : CALL section_add_keyword(section, keyword)
2327 9582 : CALL keyword_release(keyword)
2328 :
2329 : CALL keyword_create(keyword, __LOCATION__, &
2330 : name="U_MINUS_J", &
2331 : variants=["U_EFF"], &
2332 : description="Effective parameter U(eff) = U - J", &
2333 : repeats=.FALSE., &
2334 : n_var=1, &
2335 : type_of_var=real_t, &
2336 : default_r_val=0.0_dp, &
2337 : unit_str="au_e", &
2338 19164 : usage="U_MINUS_J [eV] 1.4")
2339 9582 : CALL section_add_keyword(section, keyword)
2340 9582 : CALL keyword_release(keyword)
2341 :
2342 : CALL keyword_create(keyword, __LOCATION__, &
2343 : name="N", &
2344 : description="principal quantum number of the "// &
2345 : "orbitals to which the correction is applied. Ignored unless pwdft is used for the calculations", &
2346 : repeats=.FALSE., &
2347 : n_var=1, &
2348 : type_of_var=integer_t, &
2349 : default_i_val=-1, &
2350 9582 : usage="N 2")
2351 9582 : CALL section_add_keyword(section, keyword)
2352 9582 : CALL keyword_release(keyword)
2353 :
2354 : CALL keyword_create(keyword, __LOCATION__, &
2355 : name="U", &
2356 : description="U parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2357 : repeats=.FALSE., &
2358 : n_var=1, &
2359 : type_of_var=real_t, &
2360 : default_r_val=0.0_dp, &
2361 : unit_str="au_e", &
2362 9582 : usage="U [eV] 1.4")
2363 9582 : CALL section_add_keyword(section, keyword)
2364 9582 : CALL keyword_release(keyword)
2365 :
2366 : CALL keyword_create(keyword, __LOCATION__, &
2367 : name="J", &
2368 : description="J parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2369 : repeats=.FALSE., &
2370 : n_var=1, &
2371 : type_of_var=real_t, &
2372 : default_r_val=0.0_dp, &
2373 : unit_str="au_e", &
2374 9582 : usage="J [eV] 1.4")
2375 9582 : CALL section_add_keyword(section, keyword)
2376 9582 : CALL keyword_release(keyword)
2377 :
2378 : CALL keyword_create(keyword, __LOCATION__, &
2379 : name="alpha", &
2380 : description="alpha parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2381 : repeats=.FALSE., &
2382 : n_var=1, &
2383 : type_of_var=real_t, &
2384 : default_r_val=0.0_dp, &
2385 : unit_str="au_e", &
2386 9582 : usage="alpha [eV] 1.4")
2387 9582 : CALL section_add_keyword(section, keyword)
2388 9582 : CALL keyword_release(keyword)
2389 :
2390 : CALL keyword_create(keyword, __LOCATION__, &
2391 : name="beta", &
2392 : description="beta parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2393 : repeats=.FALSE., &
2394 : n_var=1, &
2395 : type_of_var=real_t, &
2396 : default_r_val=0.0_dp, &
2397 : unit_str="au_e", &
2398 9582 : usage="beta [eV] 1.4")
2399 9582 : CALL section_add_keyword(section, keyword)
2400 9582 : CALL keyword_release(keyword)
2401 :
2402 : CALL keyword_create(keyword, __LOCATION__, &
2403 : name="J0", &
2404 : description="J0 parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2405 : repeats=.FALSE., &
2406 : n_var=1, &
2407 : type_of_var=real_t, &
2408 : default_r_val=0.0_dp, &
2409 : unit_str="au_e", &
2410 9582 : usage="J0 [eV] 1.4")
2411 9582 : CALL section_add_keyword(section, keyword)
2412 9582 : CALL keyword_release(keyword)
2413 :
2414 : CALL keyword_create(keyword, __LOCATION__, &
2415 : name="occupation", &
2416 : description="number of electrons in the hubbard shell. Ignored unless pwdft is used", &
2417 : repeats=.FALSE., &
2418 : n_var=1, &
2419 : type_of_var=real_t, &
2420 : default_r_val=0.0_dp, &
2421 9582 : usage="occupation 6")
2422 9582 : CALL section_add_keyword(section, keyword)
2423 9582 : CALL keyword_release(keyword)
2424 :
2425 : CALL keyword_create(keyword, __LOCATION__, &
2426 : name="U_RAMPING", &
2427 : description="Increase the effective U parameter stepwise using the specified "// &
2428 : "increment until the target value given by U_MINUS_J is reached.", &
2429 : repeats=.FALSE., &
2430 : n_var=1, &
2431 : type_of_var=real_t, &
2432 : default_r_val=0.0_dp, &
2433 : unit_str="au_e", &
2434 9582 : usage="U_RAMPING [eV] 0.1")
2435 9582 : CALL section_add_keyword(section, keyword)
2436 9582 : CALL keyword_release(keyword)
2437 :
2438 : CALL keyword_create(keyword, __LOCATION__, &
2439 : name="EPS_U_RAMPING", &
2440 : description="Threshold value (SCF convergence) for incrementing the effective "// &
2441 : "U value when U ramping is active.", &
2442 : repeats=.FALSE., &
2443 : n_var=1, &
2444 : type_of_var=real_t, &
2445 : default_r_val=1.0E-5_dp, &
2446 9582 : usage="EPS_U_RAMPING 1.0E-6")
2447 9582 : CALL section_add_keyword(section, keyword)
2448 9582 : CALL keyword_release(keyword)
2449 :
2450 : CALL keyword_create(keyword, __LOCATION__, &
2451 : name="INIT_U_RAMPING_EACH_SCF", &
2452 : description="Set the initial U ramping value to zero before each wavefunction optimisation. "// &
2453 : "The default is to apply U ramping only for the initial wavefunction optimisation.", &
2454 : repeats=.FALSE., &
2455 : default_l_val=.FALSE., &
2456 : lone_keyword_l_val=.TRUE., &
2457 9582 : usage="INIT_U_RAMPING_EACH_SCF on")
2458 9582 : CALL section_add_keyword(section, keyword)
2459 9582 : CALL keyword_release(keyword)
2460 :
2461 9582 : NULLIFY (subsection)
2462 :
2463 : CALL section_create(subsection, __LOCATION__, &
2464 : name="ENFORCE_OCCUPATION", &
2465 : description="Enforce and control a special (initial) orbital occupation. "// &
2466 : "Note, this feature works only for the methods MULLIKEN and LOWDIN. "// &
2467 : "It should only be used to prepare an initial configuration. An "// &
2468 : "inadequate parameter choice can easily inhibit SCF convergence.", &
2469 : n_keywords=5, &
2470 : n_subsections=0, &
2471 9582 : repeats=.FALSE.)
2472 :
2473 : CALL keyword_create(keyword, __LOCATION__, &
2474 : name="_SECTION_PARAMETERS_", &
2475 : description="Controls the activation of the ENFORCE_OCCUPATION section", &
2476 : usage="&ENFORCE_OCCUPATION ON", &
2477 : default_l_val=.FALSE., &
2478 9582 : lone_keyword_l_val=.TRUE.)
2479 9582 : CALL section_add_keyword(subsection, keyword)
2480 9582 : CALL keyword_release(keyword)
2481 :
2482 : CALL keyword_create(keyword, __LOCATION__, name="NELEC", &
2483 : variants=["N_ELECTRONS"], &
2484 : description="Number of alpha and beta electrons. An occupation (per spin) smaller than 0.5 is ignored.", &
2485 : repeats=.FALSE., &
2486 : n_var=-1, &
2487 : type_of_var=real_t, &
2488 : default_r_val=0.0_dp, &
2489 19164 : usage="NELEC 5.0 4.0")
2490 9582 : CALL section_add_keyword(subsection, keyword)
2491 9582 : CALL keyword_release(keyword)
2492 :
2493 : CALL keyword_create(keyword, __LOCATION__, &
2494 : name="ORBITALS", &
2495 : variants=["M"], &
2496 : description="Select orbitals and occupation order. An input of 1 to 2*L+1 integer values in "// &
2497 : "the range -L to L defining the M values of the spherical orbitals is expected.", &
2498 : repeats=.FALSE., &
2499 : n_var=-1, &
2500 : type_of_var=integer_t, &
2501 : default_i_val=0, &
2502 19164 : usage="ORBITALS 0 +1 -1")
2503 9582 : CALL section_add_keyword(subsection, keyword)
2504 9582 : CALL keyword_release(keyword)
2505 :
2506 : CALL keyword_create(keyword, __LOCATION__, &
2507 : name="EPS_SCF", &
2508 : description="The occupation constraint is enforced until this threshold value "// &
2509 : "for the SCF convergence criterion is reached", &
2510 : repeats=.FALSE., &
2511 : n_var=1, &
2512 : type_of_var=real_t, &
2513 : default_r_val=1.0E30_dp, &
2514 9582 : usage="EPS_SCF 0.001")
2515 9582 : CALL section_add_keyword(subsection, keyword)
2516 9582 : CALL keyword_release(keyword)
2517 :
2518 : CALL keyword_create(keyword, __LOCATION__, &
2519 : name="MAX_SCF", &
2520 : description="The occupation constraint is applied for this number of initial SCF iterations", &
2521 : repeats=.FALSE., &
2522 : n_var=1, &
2523 : type_of_var=integer_t, &
2524 : default_i_val=-1, &
2525 9582 : usage="MAX_SCF 5")
2526 9582 : CALL section_add_keyword(subsection, keyword)
2527 9582 : CALL keyword_release(keyword)
2528 :
2529 : CALL keyword_create(keyword, __LOCATION__, &
2530 : name="SMEAR", &
2531 : description="The occupation constraint is applied with smearing", &
2532 : repeats=.FALSE., &
2533 : default_l_val=.FALSE., &
2534 : lone_keyword_l_val=.TRUE., &
2535 9582 : usage="SMEAR ON")
2536 9582 : CALL section_add_keyword(subsection, keyword)
2537 9582 : CALL keyword_release(keyword)
2538 :
2539 9582 : CALL section_add_subsection(section, subsection)
2540 9582 : CALL section_release(subsection)
2541 :
2542 9582 : END SUBROUTINE create_dft_plus_u_section
2543 :
2544 : END MODULE input_cp2k_subsys
|