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