(git:374b731)
Loading...
Searching...
No Matches
input_cp2k_subsys.F
Go to the documentation of this file.
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! **************************************************************************************************
15
16 USE bibliography, ONLY: goedecker1996, &
17 guidon2010, &
19 krack2005, &
22 USE cell_types, ONLY: &
32 USE cp_units, ONLY: cp_unit_to_cp2k
33 USE input_constants, ONLY: &
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, &
70
71CONTAINS
72
73! **************************************************************************************************
74!> \brief creates the cell section
75!> \param section ...
76!> \param periodic ...
77!> \author Ole Schuett
78! **************************************************************************************************
79 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 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 "COS, EXP, LOG, LOG10, SIN, SQRT, and TAN are available.")
91 CALL create_cell_section_low(section, periodic)
92
93 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 "Check the &CELL section for further details.")
99 CALL create_cell_section_low(subsection, periodic)
100 CALL section_add_subsection(section, subsection)
101 CALL section_release(subsection)
102
103 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 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 my_periodic = use_perd_xyz
119 IF (PRESENT(periodic)) my_periodic = periodic
120
121 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 n_var=3, type_of_var=real_t, repeats=.false.)
127 CALL section_add_keyword(section, keyword)
128 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 n_var=3, type_of_var=real_t, repeats=.false.)
135 CALL section_add_keyword(section, keyword)
136 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 n_var=3, type_of_var=real_t, repeats=.false.)
143 CALL section_add_keyword(section, keyword)
144 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 n_var=3, type_of_var=real_t, repeats=.false.)
154 CALL section_add_keyword(section, keyword)
155 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 repeats=.false.)
168 CALL section_add_keyword(section, keyword)
169 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 type_of_var=lchar_t)
175 CALL section_add_keyword(section, keyword)
176 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 default_i_val=do_cell_cp2k)
187 CALL section_add_keyword(section, keyword)
188 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, &
201 default_i_val=my_periodic)
202 CALL section_add_keyword(section, keyword)
203 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 n_var=3, default_i_vals=(/1, 1, 1/), repeats=.false.)
211 CALL section_add_keyword(section, keyword)
212 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 &ne; b &ne; c &ne; a, &alpha; &ne; &beta; &ne; &gamma; &ne; &alpha; &ne; 90&deg;)", &
220 "Monoclinic (a &ne; b &ne; c, &alpha; = &gamma; = 90&deg;, &beta; &ne; 90&deg;)", &
221 "Monoclinic (a = b &ne; c, &alpha; = &beta; = 90&deg;, &gamma; &ne; 90&deg;)", &
222 "Orthorhombic (a &ne; b &ne; c, &alpha; = &beta; = &gamma; = 90&deg;)", &
223 "Tetragonal (a = b &ne; c, &alpha; = &beta; = &gamma; = 90&deg;)", &
224 "Tetragonal (a = c &ne; b, &alpha; = &beta; = &gamma; = 90&deg;)", &
225 "Tetragonal (a &ne; b = c, &alpha; = &beta; = &gamma; = 90&deg;)", &
226 "Tetragonal (alias for TETRAGONAL_AB)", &
227 "Rhombohedral (a = b = c, &alpha; = &beta; = &gamma; &ne; 90&deg;)", &
228 "Hexagonal (alias for HEXAGONAL_GAMMA_60)", &
229 "Hexagonal (a = b &ne; c, &alpha; = &beta; = 90&deg;, &gamma; = 60&deg;)", &
230 "Hexagonal (a = b &ne; c, &alpha; = &beta; = 90&deg;, &gamma; = 120&deg;)", &
231 "Cubic (a = b = c, &alpha; = &beta; = &gamma; = 90&deg;)"), &
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"), &
239 default_i_val=cell_sym_none)
240 CALL section_add_keyword(section, keyword)
241 CALL keyword_release(keyword)
242
243 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 SUBROUTINE create_rng_section(section)
251 TYPE(section_type), POINTER :: section
252
253 TYPE(keyword_type), POINTER :: keyword
254
255 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 n_keywords=1, n_subsections=0, repeats=.false.)
259 NULLIFY (keyword)
260
261 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
262 description="Specify an initial RNG stream record", repeats=.true., &
263 usage="{RNG record string}", type_of_var=lchar_t)
264 CALL section_add_keyword(section, keyword)
265 CALL keyword_release(keyword)
266
267 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 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 cpassert(.NOT. ASSOCIATED(section))
282 CALL section_create(section, __location__, name="subsys", &
283 description="a subsystem: coordinates, topology, molecules and cell", &
284 n_keywords=1, n_subsections=9, repeats=.false.)
285
286 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 default_i_vals=(/12345/))
296 CALL section_add_keyword(section, keyword)
297 CALL keyword_release(keyword)
298
299 NULLIFY (subsection)
300
301 CALL create_rng_section(subsection)
302 CALL section_add_subsection(section, subsection)
303 CALL section_release(subsection)
304
305 CALL create_cell_section(subsection)
306 CALL section_add_subsection(section, subsection)
307 CALL section_release(subsection)
308
309 CALL create_coord_section(subsection)
310 CALL section_add_subsection(section, subsection)
311 CALL section_release(subsection)
312
313 CALL create_velocity_section(subsection)
314 CALL section_add_subsection(section, subsection)
315 CALL section_release(subsection)
316
317 CALL create_kind_section(subsection)
318 CALL section_add_subsection(section, subsection)
319 CALL section_release(subsection)
320
321 CALL create_topology_section(subsection)
322 CALL section_add_subsection(section, subsection)
323 CALL section_release(subsection)
324
325 CALL create_colvar_section(section=subsection)
326 CALL section_add_subsection(section, subsection)
327 CALL section_release(subsection)
328
329 CALL create_multipole_section(subsection)
330 CALL section_add_subsection(section, subsection)
331 CALL section_release(subsection)
332
333 CALL create_shell_coord_section(subsection)
334 CALL section_add_subsection(section, subsection)
335 CALL section_release(subsection)
336
337 CALL create_shell_vel_section(subsection)
338 CALL section_add_subsection(section, subsection)
339 CALL section_release(subsection)
340 CALL create_core_coord_section(subsection)
341 CALL section_add_subsection(section, subsection)
342 CALL section_release(subsection)
343
344 CALL create_core_vel_section(subsection)
345 CALL section_add_subsection(section, subsection)
346 CALL section_release(subsection)
347
348 CALL create_subsys_print_section(subsection)
349 CALL section_add_subsection(section, subsection)
350 CALL section_release(subsection)
351
352 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 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 NULLIFY (print_key, keyword)
366 cpassert(.NOT. ASSOCIATED(section))
367 CALL section_create(section, __location__, name="print", &
368 description="Controls printings related to the subsys", &
369 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 print_level=medium_print_level, filename="__STD_OUT__")
376 CALL section_add_subsection(section, print_key)
377 CALL section_release(print_key)
378
379 CALL create_structure_data_section(print_key)
380 CALL section_add_subsection(section, print_key)
381 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 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 default_r_val=0.5_dp*bohr, unit_str="angstrom")
395 CALL section_add_keyword(print_key, keyword)
396 CALL keyword_release(keyword)
397 CALL section_add_subsection(section, print_key)
398 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 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 default_l_val=.false., lone_keyword_l_val=.true.)
406 CALL section_add_keyword(print_key, keyword)
407 CALL keyword_release(keyword)
408 CALL keyword_create(keyword, __location__, name="cif_info", &
409 description="Prints information when parsing CIF files.", &
410 default_l_val=.false., lone_keyword_l_val=.true.)
411 CALL section_add_keyword(print_key, keyword)
412 CALL keyword_release(keyword)
413 CALL keyword_create(keyword, __location__, name="pdb_info", &
414 description="Prints information when parsing PDB files.", &
415 default_l_val=.false., lone_keyword_l_val=.true.)
416 CALL section_add_keyword(print_key, keyword)
417 CALL keyword_release(keyword)
418 CALL keyword_create(keyword, __location__, name="xyz_info", &
419 description="Prints information when parsing XYZ files.", &
420 default_l_val=.false., lone_keyword_l_val=.true.)
421 CALL section_add_keyword(print_key, keyword)
422 CALL keyword_release(keyword)
423 CALL keyword_create(keyword, __location__, name="psf_info", &
424 description="Prints information when parsing PSF files.", &
425 default_l_val=.false., lone_keyword_l_val=.true.)
426 CALL section_add_keyword(print_key, keyword)
427 CALL keyword_release(keyword)
428 CALL keyword_create(keyword, __location__, name="amber_info", &
429 description="Prints information when parsing ABER topology files.", &
430 default_l_val=.false., lone_keyword_l_val=.true.)
431 CALL section_add_keyword(print_key, keyword)
432 CALL keyword_release(keyword)
433 CALL keyword_create(keyword, __location__, name="g96_info", &
434 description="Prints information when parsing G96 files.", &
435 default_l_val=.false., lone_keyword_l_val=.true.)
436 CALL section_add_keyword(print_key, keyword)
437 CALL keyword_release(keyword)
438 CALL keyword_create(keyword, __location__, name="crd_info", &
439 description="Prints information when parsing CRD files.", &
440 default_l_val=.false., lone_keyword_l_val=.true.)
441 CALL section_add_keyword(print_key, keyword)
442 CALL keyword_release(keyword)
443 CALL keyword_create(keyword, __location__, name="gtop_info", &
444 description="Prints information when parsing GROMOS topology files.", &
445 default_l_val=.false., lone_keyword_l_val=.true.)
446 CALL section_add_keyword(print_key, keyword)
447 CALL keyword_release(keyword)
448 CALL keyword_create(keyword, __location__, name="util_info", &
449 description="Prints information regarding topology utilities", &
450 default_l_val=.false., lone_keyword_l_val=.true.)
451 CALL section_add_keyword(print_key, keyword)
452 CALL keyword_release(keyword)
453 CALL keyword_create(keyword, __location__, name="generate_info", &
454 description="Prints information regarding topology generation", &
455 default_l_val=.false., lone_keyword_l_val=.true.)
456 CALL section_add_keyword(print_key, keyword)
457 CALL keyword_release(keyword)
458 CALL section_add_subsection(section, print_key)
459 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 unit_str="angstrom")
465 CALL section_add_subsection(section, print_key)
466 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 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 default_l_val=.false., lone_keyword_l_val=.true.)
476 CALL section_add_keyword(print_key, keyword)
477 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 default_l_val=.false., lone_keyword_l_val=.true.)
481 CALL section_add_keyword(print_key, keyword)
482 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 default_l_val=.false., lone_keyword_l_val=.true.)
486 CALL section_add_keyword(print_key, keyword)
487 CALL keyword_release(keyword)
488 CALL section_add_subsection(section, print_key)
489 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 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 default_l_val=.false., lone_keyword_l_val=.true.)
497 CALL section_add_keyword(print_key, keyword)
498 CALL keyword_release(keyword)
499 CALL keyword_create(keyword, __location__, name="EPS_GEO", &
500 description="Accuracy required for symmetry detection", &
501 default_r_val=1.0e-4_dp)
502 CALL section_add_keyword(print_key, keyword)
503 CALL keyword_release(keyword)
504 CALL keyword_create(keyword, __location__, name="STANDARD_ORIENTATION", &
505 description="Print molecular coordinates in standard orientation", &
506 default_l_val=.false., lone_keyword_l_val=.true.)
507 CALL section_add_keyword(print_key, keyword)
508 CALL keyword_release(keyword)
509 CALL keyword_create(keyword, __location__, name="INERTIA", &
510 description="Print molecular inertia tensor", &
511 default_l_val=.false., lone_keyword_l_val=.true.)
512 CALL section_add_keyword(print_key, keyword)
513 CALL keyword_release(keyword)
514 CALL keyword_create(keyword, __location__, name="SYMMETRY_ELEMENTS", &
515 description="Print symmetry elements", &
516 default_l_val=.false., lone_keyword_l_val=.true.)
517 CALL section_add_keyword(print_key, keyword)
518 CALL keyword_release(keyword)
519 CALL keyword_create(keyword, __location__, name="ALL", &
520 description="Print all symmetry information", &
521 default_l_val=.false., lone_keyword_l_val=.true.)
522 CALL section_add_keyword(print_key, keyword)
523 CALL keyword_release(keyword)
524 CALL keyword_create(keyword, __location__, name="ROTATION_MATRICES", &
525 description="All the rotation matrices of the point group", &
526 default_l_val=.false.)
527 CALL section_add_keyword(print_key, keyword)
528 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 default_c_val="NONE")
533 CALL section_add_keyword(print_key, keyword)
534 CALL keyword_release(keyword)
535 CALL section_add_subsection(section, print_key)
536 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 print_level=medium_print_level, filename="__STD_OUT__")
541 CALL section_add_subsection(section, print_key)
542 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 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 default_l_val=.true., lone_keyword_l_val=.true.)
551 CALL section_add_keyword(print_key, keyword)
552 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 default_l_val=.true., lone_keyword_l_val=.true.)
557 CALL section_add_keyword(print_key, keyword)
558 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 default_l_val=.true., lone_keyword_l_val=.true.)
563 CALL section_add_keyword(print_key, keyword)
564 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 default_l_val=.true., lone_keyword_l_val=.true.)
569 CALL section_add_keyword(print_key, keyword)
570 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 default_l_val=.true., lone_keyword_l_val=.true.)
575 CALL section_add_keyword(print_key, keyword)
576 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 default_l_val=.true., lone_keyword_l_val=.true.)
582 CALL section_add_keyword(print_key, keyword)
583 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 default_l_val=.true., lone_keyword_l_val=.true.)
589 CALL section_add_keyword(print_key, keyword)
590 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 default_l_val=.true., lone_keyword_l_val=.true.)
595 CALL section_add_keyword(print_key, keyword)
596 CALL keyword_release(keyword)
597
598 CALL section_add_subsection(section, print_key)
599 CALL section_release(print_key)
600
601 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 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 cpassert(.NOT. ASSOCIATED(section))
615 CALL section_create(section, __location__, name="multipoles", &
616 description="Specifies the dipoles and quadrupoles for particles.", &
617 n_keywords=1, n_subsections=0, repeats=.false.)
618
619 NULLIFY (keyword, subsection)
620 CALL section_create(subsection, __location__, name="dipoles", &
621 description="Specifies the dipoles of the particles.", &
622 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 type_of_var=real_t, n_var=3)
628 CALL section_add_keyword(subsection, keyword)
629 CALL keyword_release(keyword)
630 CALL section_add_subsection(section, subsection)
631 CALL section_release(subsection)
632
633 CALL section_create(subsection, __location__, name="quadrupoles", &
634 description="Specifies the quadrupoles of the particles.", &
635 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 type_of_var=real_t, n_var=6)
641 CALL section_add_keyword(subsection, keyword)
642 CALL keyword_release(keyword)
643 CALL section_add_subsection(section, subsection)
644 CALL section_release(subsection)
645
646 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 SUBROUTINE create_structure_data_section(print_key)
654 TYPE(section_type), POINTER :: print_key
655
656 TYPE(keyword_type), POINTER :: keyword
657
658 cpassert(.NOT. ASSOCIATED(print_key))
659
660 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 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 type_of_var=integer_t)
672 CALL section_add_keyword(print_key, keyword)
673 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 type_of_var=integer_t)
680 CALL section_add_keyword(print_key, keyword)
681 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 type_of_var=integer_t)
687 CALL section_add_keyword(print_key, keyword)
688 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 type_of_var=integer_t)
694 CALL section_add_keyword(print_key, keyword)
695 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 repeats=.true., type_of_var=integer_t)
702 CALL section_add_keyword(print_key, keyword)
703 CALL keyword_release(keyword)
704
705 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 SUBROUTINE create_velocity_section(section)
713 TYPE(section_type), POINTER :: section
714
715 TYPE(keyword_type), POINTER :: keyword
716
717 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 n_keywords=1, n_subsections=0, repeats=.false.)
722 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 default_c_val="bohr*au_t^-1")
729 CALL section_add_keyword(section, keyword)
730 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 type_of_var=real_t, n_var=3)
738 CALL section_add_keyword(section, keyword)
739 CALL keyword_release(keyword)
740
741 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 SUBROUTINE create_shell_vel_section(section)
749 TYPE(section_type), POINTER :: section
750
751 TYPE(keyword_type), POINTER :: keyword
752
753 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 n_keywords=1, n_subsections=0, repeats=.false.)
758 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 type_of_var=real_t, n_var=3)
766 CALL section_add_keyword(section, keyword)
767 CALL keyword_release(keyword)
768
769 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 SUBROUTINE create_core_vel_section(section)
777 TYPE(section_type), POINTER :: section
778
779 TYPE(keyword_type), POINTER :: keyword
780
781 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 n_keywords=1, n_subsections=0, repeats=.false.)
786 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 type_of_var=real_t, n_var=3)
794 CALL section_add_keyword(section, keyword)
795 CALL keyword_release(keyword)
796
797 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 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 n_keywords=1, n_subsections=0, repeats=.false.)
812 NULLIFY (keyword)
813 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
814 description="CP2K Pseudo Potential Standard Format (GTH, ALL)", &
815 repeats=.true., type_of_var=lchar_t)
816 CALL section_add_keyword(section, keyword)
817 CALL keyword_release(keyword)
818
819 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 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 n_keywords=1, n_subsections=0, repeats=.false.)
834 NULLIFY (keyword)
835 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
836 description="CP2K KG TNADD Potential Standard Format (TNADD)", &
837 repeats=.true., type_of_var=lchar_t)
838 CALL section_add_keyword(section, keyword)
839 CALL keyword_release(keyword)
840
841 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 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 n_keywords=1, n_subsections=0, repeats=.true.)
856
857 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 usage="Orbital", default_c_val="Orbital")
863 CALL section_add_keyword(section, keyword)
864 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 "```")
891 CALL section_add_keyword(section, keyword)
892 CALL keyword_release(keyword)
893
894 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 SUBROUTINE create_coord_section(section)
902 TYPE(section_type), POINTER :: section
903
904 TYPE(keyword_type), POINTER :: keyword
905
906 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 n_keywords=1, n_subsections=0, repeats=.false.)
915 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 usage="UNIT angstrom", default_c_val="angstrom")
920 CALL section_add_keyword(section, keyword)
921 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 lone_keyword_l_val=.true.)
928 CALL section_add_keyword(section, keyword)
929 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 type_of_var=lchar_t)
938 CALL section_add_keyword(section, keyword)
939 CALL keyword_release(keyword)
940 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 SUBROUTINE create_shell_coord_section(section)
948 TYPE(section_type), POINTER :: section
949
950 TYPE(keyword_type), POINTER :: keyword
951
952 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 n_keywords=1, n_subsections=0, repeats=.false.)
957 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 usage="UNIT angstrom", default_c_val="angstrom")
962 CALL section_add_keyword(section, keyword)
963 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 lone_keyword_l_val=.true.)
970 CALL section_add_keyword(section, keyword)
971 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 type_of_var=lchar_t)
979 CALL section_add_keyword(section, keyword)
980 CALL keyword_release(keyword)
981
982 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 SUBROUTINE create_core_coord_section(section)
990 TYPE(section_type), POINTER :: section
991
992 TYPE(keyword_type), POINTER :: keyword
993
994 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 n_keywords=1, n_subsections=0, repeats=.false.)
999 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 usage="UNIT angstrom", default_c_val="angstrom")
1004 CALL section_add_keyword(section, keyword)
1005 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 lone_keyword_l_val=.true.)
1012 CALL section_add_keyword(section, keyword)
1013 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 type_of_var=lchar_t)
1021 CALL section_add_keyword(section, keyword)
1022 CALL keyword_release(keyword)
1023
1024 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 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 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 n_keywords=19, n_subsections=1, repeats=.true.)
1042
1043 NULLIFY (keyword)
1044
1045 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1046 description="The name of the kind described in this section.", &
1047 usage="H", default_c_val="DEFAULT")
1048 CALL section_add_keyword(section, keyword)
1049 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 repeats=.true., n_var=-1)
1060 CALL section_add_keyword(section, keyword)
1061 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 removed=.true.)
1073 CALL section_add_keyword(section, keyword)
1074 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 removed=.true.)
1084 CALL section_add_keyword(section, keyword)
1085 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 removed=.true.)
1095 CALL section_add_keyword(section, keyword)
1096 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 removed=.true.)
1107 CALL section_add_keyword(section, keyword)
1108 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 n_var=-1, type_of_var=integer_t)
1116 CALL section_add_keyword(section, keyword)
1117 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 default_r_val=0.0_dp)
1123 CALL section_add_keyword(section, keyword)
1124 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 default_r_val=0.0_dp)
1131 CALL section_add_keyword(section, keyword)
1132 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 usage="ELEMENT O", type_of_var=char_t, n_var=1)
1139 CALL section_add_keyword(section, keyword)
1140 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 usage="MASS 2.0", type_of_var=real_t, n_var=1)
1147 CALL section_add_keyword(section, keyword)
1148 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 usage="POTENTIAL_FILE_NAME <PSEUDO-POTENTIAL-FILE-NAME>", default_c_val="-", n_var=1)
1154 CALL section_add_keyword(section, keyword)
1155 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 usage="POTENTIAL_TYPE <TYPE>", default_c_val="", n_var=1)
1161 CALL section_add_keyword(section, keyword)
1162 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 citations=(/goedecker1996, hartwigsen1998, krack2005/), n_var=-1)
1169 CALL section_add_keyword(section, keyword)
1170 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 usage="KG_POTENTIAL_FILE_NAME <POTENTIAL-FILE-NAME>", default_c_val="-", n_var=1)
1176 CALL section_add_keyword(section, keyword)
1177 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 usage="KG_POTENTIAL <TNADD-POTENTIAL-NAME>", default_c_val="NONE", n_var=1)
1183 CALL section_add_keyword(section, keyword)
1184 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 usage="ECP_SEMI_LOCAL {T,F}", default_l_val=.true., lone_keyword_l_val=.true.)
1191 CALL section_add_keyword(section, keyword)
1192 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 unit_str="angstrom")
1199 CALL section_add_keyword(section, keyword)
1200 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 usage="VDW_RADIUS 1.85", n_var=1, default_r_val=0.0_dp, unit_str="angstrom")
1206 CALL section_add_keyword(section, keyword)
1207 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 usage="HARD_EXP_RADIUS 0.9", type_of_var=real_t, n_var=1)
1213 CALL section_add_keyword(section, keyword)
1214 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 usage="MAX_RAD_LOCAL 15.0", default_r_val=13.0_dp*bohr)
1220 CALL section_add_keyword(section, keyword)
1221 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 usage="RHO_EXP_RADIUS 0.9", type_of_var=real_t, n_var=1)
1228 CALL section_add_keyword(section, keyword)
1229 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 usage="LEBEDEV_GRID 40", default_i_val=50)
1235 CALL section_add_keyword(section, keyword)
1236 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 usage="RADIAL_GRID 70", default_i_val=50)
1242 CALL section_add_keyword(section, keyword)
1243 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 unit_str="angstrom", n_var=1)
1257 CALL section_add_keyword(section, keyword)
1258 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 usage="DFTB3_PARAM 0.2", default_r_val=0.0_dp)
1264 CALL section_add_keyword(section, keyword)
1265 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 usage="LMAX_DFTB 1", default_i_val=-1)
1270 CALL section_add_keyword(section, keyword)
1271 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 usage="MAO 4", default_i_val=-1)
1276 CALL section_add_keyword(section, keyword)
1277 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 usage="SE_P_ORBITALS_ON_H", default_l_val=.false., lone_keyword_l_val=.true.)
1285 CALL section_add_keyword(section, keyword)
1286 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 usage="GPW_TYPE", default_l_val=.false., lone_keyword_l_val=.true.)
1292 CALL section_add_keyword(section, keyword)
1293 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 lone_keyword_l_val=.true.)
1304 CALL section_add_keyword(section, keyword)
1305 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 lone_keyword_l_val=.true.)
1315 CALL section_add_keyword(section, keyword)
1316 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 lone_keyword_l_val=.true.)
1325 CALL section_add_keyword(section, keyword)
1326 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 "By default PAO is disbabled.", default_i_val=0)
1333 CALL section_add_keyword(section, keyword)
1334 CALL keyword_release(keyword)
1335
1336 NULLIFY (subsection)
1337 CALL create_pao_potential_section(subsection)
1338 CALL section_add_subsection(section, subsection)
1339 CALL section_release(subsection)
1340
1341 CALL create_pao_descriptor_section(subsection)
1342 CALL section_add_subsection(section, subsection)
1343 CALL section_release(subsection)
1344
1345 CALL create_basis_section(subsection)
1346 CALL section_add_subsection(section, subsection)
1347 CALL section_release(subsection)
1348
1349 CALL create_potential_section(subsection)
1350 CALL section_add_subsection(section, subsection)
1351 CALL section_release(subsection)
1352
1353 CALL create_kgpot_section(subsection)
1354 CALL section_add_subsection(section, subsection)
1355 CALL section_release(subsection)
1356
1357 CALL create_dft_plus_u_section(subsection)
1358 CALL section_add_subsection(section, subsection)
1359 CALL section_release(subsection)
1360
1361 CALL create_bs_section(subsection)
1362 CALL section_add_subsection(section, subsection)
1363 CALL section_release(subsection)
1364
1365 END SUBROUTINE create_kind_section
1366
1367! **************************************************************************************************
1368!> \brief Creates the PAO_POTENTIAL section
1369!> \param section the section to create
1370!> \author Ole Schuett
1371! **************************************************************************************************
1372 SUBROUTINE create_pao_potential_section(section)
1373 TYPE(section_type), POINTER :: section
1374
1375 TYPE(keyword_type), POINTER :: keyword
1376
1377 cpassert(.NOT. ASSOCIATED(section))
1378 NULLIFY (keyword)
1379
1380 CALL section_create(section, __location__, name="PAO_POTENTIAL", repeats=.true., &
1381 description="Settings of the PAO potentials, which are atomic kind specific.")
1382
1383 CALL keyword_create(keyword, __location__, name="MAXL", &
1384 description="Maximum angular moment of the potential "// &
1385 "(must be an even number).", default_i_val=0)
1386 CALL section_add_keyword(section, keyword)
1387 CALL keyword_release(keyword)
1388
1389 CALL keyword_create(keyword, __location__, name="BETA", &
1390 description="Exponent of the Gaussian potential term.", &
1391 default_r_val=1.0_dp)
1392 CALL section_add_keyword(section, keyword)
1393 CALL keyword_release(keyword)
1394
1395 CALL keyword_create(keyword, __location__, name="WEIGHT", &
1396 description="Weight of Gaussian potential term.", &
1397 default_r_val=1.0_dp)
1398 CALL section_add_keyword(section, keyword)
1399 CALL keyword_release(keyword)
1400
1401 CALL keyword_create(keyword, __location__, name="MAX_PROJECTOR", &
1402 description="Maximum angular moment of the potential's projectors. "// &
1403 "Used only by the GTH parametrization", default_i_val=2)
1404 CALL section_add_keyword(section, keyword)
1405 CALL keyword_release(keyword)
1406
1407 END SUBROUTINE create_pao_potential_section
1408
1409! **************************************************************************************************
1410!> \brief Creates the PAO_DESCRIPTOR section
1411!> \param section the section to create
1412!> \author Ole Schuett
1413! **************************************************************************************************
1414 SUBROUTINE create_pao_descriptor_section(section)
1415 TYPE(section_type), POINTER :: section
1416
1417 TYPE(keyword_type), POINTER :: keyword
1418
1419 cpassert(.NOT. ASSOCIATED(section))
1420 NULLIFY (keyword)
1421
1422 CALL section_create(section, __location__, name="PAO_DESCRIPTOR", repeats=.true., &
1423 description="Settings of the PAO descriptor, which are atomic kind specific.")
1424
1425 CALL keyword_create(keyword, __location__, name="BETA", &
1426 description="Exponent of the Gaussian potential term.", &
1427 default_r_val=1.0_dp)
1428 CALL section_add_keyword(section, keyword)
1429 CALL keyword_release(keyword)
1430
1431 CALL keyword_create(keyword, __location__, name="SCREENING", &
1432 description="Exponent of the Gaussian screening.", &
1433 default_r_val=0.2_dp)
1434 CALL section_add_keyword(section, keyword)
1435 CALL keyword_release(keyword)
1436
1437 CALL keyword_create(keyword, __location__, name="WEIGHT", &
1438 description="Weight of Gaussian potential term.", &
1439 default_r_val=1.0_dp)
1440 CALL section_add_keyword(section, keyword)
1441 CALL keyword_release(keyword)
1442
1443 END SUBROUTINE create_pao_descriptor_section
1444
1445! **************************************************************************************************
1446!> \brief Create CP2K input section for BS method: imposing atomic orbital occupation
1447!> different from default in initialization of the density matrix
1448!> it works only with GUESS ATOMIC
1449!> \param section ...
1450!> \date 05.08.2009
1451!> \author MI
1452!> \version 1.0
1453! **************************************************************************************************
1454 SUBROUTINE create_bs_section(section)
1455
1456 TYPE(section_type), POINTER :: section
1457
1458 TYPE(keyword_type), POINTER :: keyword
1459 TYPE(section_type), POINTER :: subsection
1460
1461 cpassert(.NOT. ASSOCIATED(section))
1462
1463 CALL section_create(section, __location__, &
1464 name="BS", &
1465 description="Define the required atomic orbital occupation "// &
1466 "assigned in initialization of the density matrix, by adding or "// &
1467 "subtracting electrons from specific angular momentum channels. "// &
1468 "It works only with GUESS ATOMIC.", &
1469 n_keywords=0, &
1470 n_subsections=2, &
1471 repeats=.false.)
1472
1473 NULLIFY (keyword, subsection)
1474
1475 CALL keyword_create(keyword, __location__, &
1476 name="_SECTION_PARAMETERS_", &
1477 description="controls the activation of the BS section", &
1478 usage="&BS ON", &
1479 default_l_val=.false., &
1480 lone_keyword_l_val=.true.)
1481 CALL section_add_keyword(section, keyword)
1482 CALL keyword_release(keyword)
1483
1484 CALL section_create(subsection, __location__, name="ALPHA", description="alpha spin", &
1485 n_keywords=3, &
1486 n_subsections=0, &
1487 repeats=.false.)
1488
1489 CALL keyword_create(keyword, __location__, &
1490 name="NEL", &
1491 description="Orbital ccupation change per angular momentum quantum number. "// &
1492 "In unrestricted calculations applied to spin alpha.", &
1493 repeats=.false., &
1494 n_var=-1, &
1495 default_i_val=-1, &
1496 usage="NEL 2")
1497 CALL section_add_keyword(subsection, keyword)
1498 CALL keyword_release(keyword)
1499
1500 CALL keyword_create(keyword, __location__, &
1501 name="L", &
1502 variants=(/"L"/), &
1503 description="Angular momentum quantum number of the "// &
1504 "orbitals whose occupation is changed", &
1505 repeats=.false., &
1506 n_var=-1, &
1507 default_i_val=-1, &
1508 usage="L 2")
1509 CALL section_add_keyword(subsection, keyword)
1510 CALL keyword_release(keyword)
1511
1512 CALL keyword_create(keyword, __location__, &
1513 name="N", &
1514 variants=(/"N"/), &
1515 description="Principal quantum number of the "// &
1516 "orbitals whose occupation is changed. "// &
1517 "Default is the first not occupied", &
1518 repeats=.false., &
1519 n_var=-1, &
1520 default_i_val=0, &
1521 usage="N 2")
1522 CALL section_add_keyword(subsection, keyword)
1523 CALL keyword_release(keyword)
1524 CALL section_add_subsection(section, subsection)
1525 CALL section_release(subsection)
1526
1527 CALL section_create(subsection, __location__, name="BETA", description="beta spin", &
1528 n_keywords=3, &
1529 n_subsections=0, &
1530 repeats=.false.)
1531
1532 CALL keyword_create(keyword, __location__, &
1533 name="NEL", &
1534 description="Orbital ccupation change per angular momentum quantum number. "// &
1535 "Applied to spin beta and active only in unrestricted calculations.", &
1536 repeats=.false., &
1537 n_var=-1, &
1538 default_i_val=-1, &
1539 usage="NEL 2")
1540 CALL section_add_keyword(subsection, keyword)
1541 CALL keyword_release(keyword)
1542
1543 CALL keyword_create(keyword, __location__, &
1544 name="L", &
1545 description="Angular momentum quantum number of the "// &
1546 "orbitals of beta spin whose occupation is changed. "// &
1547 "Active only for unrestricted calculations", &
1548 repeats=.false., &
1549 n_var=-1, &
1550 default_i_val=-1, &
1551 usage="L 2")
1552 CALL section_add_keyword(subsection, keyword)
1553 CALL keyword_release(keyword)
1554
1555 CALL keyword_create(keyword, __location__, &
1556 name="N", &
1557 description="Principal quantum number of the "// &
1558 "orbitals of beta spin whose occupation is changed. "// &
1559 "Default is the first not occupied. "// &
1560 "Active only for unrestricted calculations", &
1561 repeats=.false., &
1562 n_var=-1, &
1563 default_i_val=0, &
1564 usage="N 2")
1565 CALL section_add_keyword(subsection, keyword)
1566 CALL keyword_release(keyword)
1567
1568 CALL section_add_subsection(section, subsection)
1569 CALL section_release(subsection)
1570
1571 END SUBROUTINE create_bs_section
1572
1573! **************************************************************************************************
1574!> \brief Create the topology section for FIST.. and the base is running running...
1575!> Contains all information regarding topology to be read in input file..
1576!> \param section the section to create
1577!> \author teo
1578! **************************************************************************************************
1579 SUBROUTINE create_topology_section(section)
1580 TYPE(section_type), POINTER :: section
1581
1582 TYPE(keyword_type), POINTER :: keyword
1583 TYPE(section_type), POINTER :: print_key, subsection
1584
1585 cpassert(.NOT. ASSOCIATED(section))
1586 CALL section_create(section, __location__, name="TOPOLOGY", &
1587 description="Section specifying information regarding how to handle the topology"// &
1588 " for classical runs.", &
1589 n_keywords=5, n_subsections=0, repeats=.false.)
1590
1591 NULLIFY (keyword, print_key)
1592 ! Logical
1593 CALL keyword_create(keyword, __location__, name="USE_ELEMENT_AS_KIND", &
1594 description="Kinds are generated according to the element name."// &
1595 " Default=True for SE and TB methods.", &
1596 usage="USE_ELEMENT_AS_KIND logical", &
1597 default_l_val=.false., lone_keyword_l_val=.true.)
1598 CALL section_add_keyword(section, keyword)
1599 CALL keyword_release(keyword)
1600
1601 CALL keyword_create(keyword, __location__, name="CHARGE_OCCUP", &
1602 variants=(/"CHARGE_O"/), &
1603 description="Read MM charges from the OCCUP field of PDB file.", &
1604 usage="CHARGE_OCCUP logical", &
1605 default_l_val=.false., lone_keyword_l_val=.true.)
1606 CALL section_add_keyword(section, keyword)
1607 CALL keyword_release(keyword)
1608
1609 CALL keyword_create(keyword, __location__, name="CHARGE_BETA", &
1610 variants=(/"CHARGE_B"/), &
1611 description="Read MM charges from the BETA field of PDB file.", &
1612 usage="CHARGE_BETA logical", &
1613 default_l_val=.false., lone_keyword_l_val=.true.)
1614 CALL section_add_keyword(section, keyword)
1615 CALL keyword_release(keyword)
1616
1617 CALL keyword_create(keyword, __location__, name="CHARGE_EXTENDED", &
1618 description="Read MM charges from the very last field of PDB file (starting from column 81)."// &
1619 " No limitations of number of digits.", &
1620 usage="CHARGE_EXTENDED logical", &
1621 default_l_val=.false., lone_keyword_l_val=.true.)
1622 CALL section_add_keyword(section, keyword)
1623 CALL keyword_release(keyword)
1624
1625 CALL keyword_create(keyword, __location__, name="PARA_RES", &
1626 description="For a protein, each residue is now considered a molecule", &
1627 usage="PARA_RES logical", &
1628 default_l_val=.true., lone_keyword_l_val=.true.)
1629 CALL section_add_keyword(section, keyword)
1630 CALL keyword_release(keyword)
1631
1632 CALL keyword_create(keyword, __location__, name="MOL_CHECK", &
1633 description="Check molecules have the same number of atom and names.", &
1634 usage="MOL_CHECK logical", &
1635 default_l_val=.true., lone_keyword_l_val=.true.)
1636 CALL section_add_keyword(section, keyword)
1637 CALL keyword_release(keyword)
1638
1639 CALL keyword_create(keyword, __location__, name="USE_G96_VELOCITY", &
1640 description="Use the velocities in the G96 coordinate files as the starting velocity", &
1641 usage="USE_G96_VELOCITY logical", &
1642 default_l_val=.false., lone_keyword_l_val=.true.)
1643 CALL section_add_keyword(section, keyword)
1644 CALL keyword_release(keyword)
1645
1646 ! Character
1647 CALL keyword_create(keyword, __location__, name="COORD_FILE_NAME", &
1648 variants=s2a("COORD_FILE"), &
1649 description="Specifies the filename that contains coordinates.", &
1650 usage="COORD_FILE_NAME <FILENAME>", type_of_var=lchar_t)
1651 CALL section_add_keyword(section, keyword)
1652 CALL keyword_release(keyword)
1653
1654 CALL keyword_create(keyword, __location__, name="COORD_FILE_FORMAT", &
1655 variants=s2a("COORDINATE"), &
1656 description="Set up the way in which coordinates will be read.", &
1657 usage="COORD_FILE_FORMAT (OFF|PDB|XYZ|G96|CRD|CIF|XTL|CP2K)", &
1658 enum_c_vals=s2a("OFF", "PDB", "XYZ", "G96", "CRD", "CIF", "XTL", "CP2K"), &
1661 enum_desc=s2a( &
1662 "Coordinates read in the &COORD section of the input file", &
1663 "Coordinates provided through a PDB file format", &
1664 "Coordinates provided through an XYZ file format", &
1665 "Coordinates provided through a GROMOS96 file format", &
1666 "Coordinates provided through an AMBER file format", &
1667 "Coordinates provided through a CIF (Crystallographic Information File) file format", &
1668 "Coordinates provided through a XTL (MSI native) file format", &
1669 "Read the coordinates in CP2K &COORD section format from an external file. "// &
1670 "NOTE: This file will be overwritten with the latest coordinates."), &
1671 default_i_val=do_coord_off)
1672 CALL section_add_keyword(section, keyword)
1673 CALL keyword_release(keyword)
1674
1675 CALL keyword_create(keyword, __location__, name="NUMBER_OF_ATOMS", &
1676 variants=s2a("NATOMS", "NATOM"), &
1677 description="Optionally define the number of atoms read from an external file "// &
1678 "(see COORD_FILE_NAME) if the COORD_FILE_FORMAT CP2K is used", &
1679 repeats=.false., &
1680 n_var=1, &
1681 type_of_var=integer_t, &
1682 default_i_val=-1, &
1683 usage="NATOMS 768000")
1684 CALL section_add_keyword(section, keyword)
1685 CALL keyword_release(keyword)
1686
1687 CALL connectivity_framework(section, do_conn_generate)
1688
1689 CALL keyword_create(keyword, __location__, name="DISABLE_EXCLUSION_LISTS", &
1690 description="Do not build any exclusion lists.", &
1691 usage="DISABLE_EXCLUSION_LISTS", &
1692 default_l_val=.false., lone_keyword_l_val=.true.)
1693 CALL section_add_keyword(section, keyword)
1694 CALL keyword_release(keyword)
1695
1696 CALL keyword_create(keyword, __location__, name="EXCLUDE_VDW", &
1697 description="Specifies which kind of Van der Waals interaction to skip.", &
1698 usage="EXCLUDE_VDW (1-1||1-2||1-3||1-4)", &
1699 enum_c_vals=s2a("1-1", "1-2", "1-3", "1-4"), &
1700 enum_i_vals=(/do_skip_11, do_skip_12, do_skip_13, do_skip_14/), &
1701 default_i_val=do_skip_13)
1702 CALL section_add_keyword(section, keyword)
1703 CALL keyword_release(keyword)
1704
1705 CALL keyword_create(keyword, __location__, name="EXCLUDE_EI", &
1706 description="Specifies which kind of Electrostatic interaction to skip.", &
1707 usage="EXCLUDE_EI (1-1||1-2||1-3||1-4)", &
1708 enum_c_vals=s2a("1-1", "1-2", "1-3", "1-4"), &
1709 enum_i_vals=(/do_skip_11, do_skip_12, do_skip_13, do_skip_14/), &
1710 default_i_val=do_skip_13)
1711 CALL section_add_keyword(section, keyword)
1712 CALL keyword_release(keyword)
1713
1714 CALL keyword_create(keyword, __location__, name="AUTOGEN_EXCLUDE_LISTS", &
1715 description="When True, the exclude lists are solely based on"// &
1716 " the bond data in the topology. The (minimal)"// &
1717 " number of bonds between two atoms is used to"// &
1718 " determine if the atom pair is added to an"// &
1719 " exclusion list. When False, 1-2 exclusion is based"// &
1720 " on bonds in the topology, 1-3 exclusion is based"// &
1721 " on bonds and bends in the topology, 1-4 exclusion"// &
1722 " is based on bonds, bends and dihedrals in the"// &
1723 " topology. This implies that a missing dihedral in"// &
1724 " the topology will cause the corresponding 1-4 pair"// &
1725 " not to be in the exclusion list, in case 1-4"// &
1726 " exclusion is requested for VDW or EI interactions.", &
1727 usage="AUTOGEN_EXCLUDE_LISTS logical", &
1728 default_l_val=.false., lone_keyword_l_val=.true.)
1729 CALL section_add_keyword(section, keyword)
1730 CALL keyword_release(keyword)
1731
1732 CALL keyword_create( &
1733 keyword, __location__, name="MULTIPLE_UNIT_CELL", &
1734 description="Specifies the numbers of repetition in space (X, Y, Z) of the defined cell, "// &
1735 "assuming it as a unit cell. This keyword affects only the coordinates specification. The same keyword "// &
1736 "in SUBSYS%CELL%MULTIPLE_UNIT_CELL should be modified in order to affect the cell "// &
1737 "specification.", usage="MULTIPLE_UNIT_CELL 1 1 1", &
1738 n_var=3, default_i_vals=(/1, 1, 1/), repeats=.false.)
1739 CALL section_add_keyword(section, keyword)
1740 CALL keyword_release(keyword)
1741
1742 CALL keyword_create(keyword, __location__, name="MEMORY_PROGRESSION_FACTOR", &
1743 description="This keyword is quite technical and should normally not be changed by the user. It "// &
1744 "affects the memory allocation during the construction of the topology. It does NOT affect the "// &
1745 "memory used once the topology is built.", &
1746 n_var=1, default_r_val=1.2_dp, repeats=.false.)
1747 CALL section_add_keyword(section, keyword)
1748 CALL keyword_release(keyword)
1749
1750 CALL cp_print_key_section_create(print_key, __location__, "DUMP_PDB", &
1751 description="controls the dumping of the PDB at the starting geometry", &
1752 print_level=debug_print_level, filename="dump")
1753 CALL section_add_subsection(section, print_key)
1754
1755 CALL keyword_create(keyword, __location__, name="CHARGE_OCCUP", &
1756 variants=(/"CHARGE_O"/), &
1757 description="Write the MM charges to the OCCUP field of the PDB file", &
1758 usage="CHARGE_OCCUP logical", &
1759 default_l_val=.false., lone_keyword_l_val=.true.)
1760 CALL section_add_keyword(print_key, keyword)
1761 CALL keyword_release(keyword)
1762
1763 CALL keyword_create(keyword, __location__, name="CHARGE_BETA", &
1764 variants=(/"CHARGE_B"/), &
1765 description="Write the MM charges to the BETA field of the PDB file", &
1766 usage="CHARGE_BETA logical", &
1767 default_l_val=.false., lone_keyword_l_val=.true.)
1768 CALL section_add_keyword(print_key, keyword)
1769 CALL keyword_release(keyword)
1770
1771 CALL keyword_create(keyword, __location__, name="CHARGE_EXTENDED", &
1772 description="Write the MM charges to the very last field of the PDB file (starting from column 81)", &
1773 usage="CHARGE_EXTENDED logical", &
1774 default_l_val=.false., lone_keyword_l_val=.true.)
1775 CALL section_add_keyword(print_key, keyword)
1776 CALL keyword_release(keyword)
1777
1778 CALL section_release(print_key)
1779
1780 CALL cp_print_key_section_create(print_key, __location__, "DUMP_PSF", &
1781 description="controls the dumping of the PSF connectivity", &
1782 print_level=debug_print_level, filename="dump")
1783 CALL section_add_subsection(section, print_key)
1784 CALL section_release(print_key)
1785
1786 NULLIFY (subsection)
1787 CALL create_exclude_list_section(subsection, "EXCLUDE_VDW_LIST")
1788 CALL section_add_subsection(section, subsection)
1789 CALL section_release(subsection)
1790
1791 CALL create_exclude_list_section(subsection, "EXCLUDE_EI_LIST")
1792 CALL section_add_subsection(section, subsection)
1793 CALL section_release(subsection)
1794
1795 CALL create_center_section(subsection)
1796 CALL section_add_subsection(section, subsection)
1797 CALL section_release(subsection)
1798
1799 CALL create_generate_section(subsection)
1800 CALL section_add_subsection(section, subsection)
1801 CALL section_release(subsection)
1802
1803 CALL create_molset_section(subsection)
1804 CALL section_add_subsection(section, subsection)
1805 CALL section_release(subsection)
1806
1807 END SUBROUTINE create_topology_section
1808
1809! **************************************************************************************************
1810!> \brief Setup a list of fine exclusion elements
1811!> \param section the section to create
1812!> \param header ...
1813!> \author Teodoro Laino [tlaino] - 12.2009
1814! **************************************************************************************************
1815 SUBROUTINE create_exclude_list_section(section, header)
1816 TYPE(section_type), POINTER :: section
1817 CHARACTER(LEN=*), INTENT(IN) :: header
1818
1819 TYPE(keyword_type), POINTER :: keyword
1820
1821 cpassert(.NOT. ASSOCIATED(section))
1822 NULLIFY (keyword)
1823 CALL section_create(section, __location__, trim(header), &
1824 description="Speficy bonds (via atom kinds) for fine tuning of 1-2 "// &
1825 "exclusion lists. If this section is not present the 1-2 exclusion is "// &
1826 "applied to all bond kinds. When this section is present the 1-2 exclusion "// &
1827 "is applied ONLY to the bonds defined herein. This section allows ONLY fine tuning of 1-2 "// &
1828 "interactions. ", &
1829 n_keywords=1, n_subsections=0, repeats=.false.)
1830
1831 CALL keyword_create(keyword, __location__, name="BOND", &
1832 description="Specify the atom kinds involved in the bond for which 1-2 exclusion holds.", &
1833 usage="BOND {KIND1} {KIND2}", type_of_var=char_t, &
1834 n_var=2)
1835 CALL section_add_keyword(section, keyword)
1836 CALL keyword_release(keyword)
1837 END SUBROUTINE create_exclude_list_section
1838
1839! **************************************************************************************************
1840!> \brief Specify keywords used to center molecule in the box
1841!> \param section the section to create
1842!> \author Teodoro Laino [tlaino] - University of Zurich - 06.2009
1843! **************************************************************************************************
1844 SUBROUTINE create_center_section(section)
1845 TYPE(section_type), POINTER :: section
1846
1847 TYPE(keyword_type), POINTER :: keyword
1848
1849 cpassert(.NOT. ASSOCIATED(section))
1850 NULLIFY (keyword)
1851 CALL section_create(section, __location__, "CENTER_COORDINATES", &
1852 description="Allows centering the coordinates of the system in the box. "// &
1853 "The centering point can be defined by the user.", &
1854 n_keywords=1, n_subsections=0, repeats=.false.)
1855
1856 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
1857 description="Controls the activation of the centering method", &
1858 usage="&CENTER_COORDINATES T", &
1859 default_l_val=.false., &
1860 lone_keyword_l_val=.true.)
1861 CALL section_add_keyword(section, keyword)
1862 CALL keyword_release(keyword)
1863
1864 CALL keyword_create(keyword, __location__, name="CENTER_POINT", &
1865 description="Specify the point used for centering the coordinates. Default is to "// &
1866 "center the system in cell/2. ", type_of_var=real_t, n_var=3, &
1867 repeats=.false.)
1868 CALL section_add_keyword(section, keyword)
1869 CALL keyword_release(keyword)
1870 END SUBROUTINE create_center_section
1871
1872! **************************************************************************************************
1873!> \brief Specify keywords used to setup several molecules with few connectivity files
1874!> \param section the section to create
1875!> \author Teodoro Laino [tlaino] - University of Zurich - 08.2008
1876! **************************************************************************************************
1877 SUBROUTINE create_molset_section(section)
1878 TYPE(section_type), POINTER :: section
1879
1880 TYPE(keyword_type), POINTER :: keyword
1881 TYPE(section_type), POINTER :: subsection, subsubsection
1882
1883 cpassert(.NOT. ASSOCIATED(section))
1884 NULLIFY (keyword, subsection, subsubsection)
1885 CALL section_create(section, __location__, name="MOL_SET", &
1886 description="Specify the connectivity of a full system specifying the connectivity"// &
1887 " of the fragments of the system.", &
1888 n_keywords=2, n_subsections=0, repeats=.false.)
1889
1890 ! MOLECULES
1891 CALL section_create(subsection, __location__, name="MOLECULE", &
1892 description="Specify information about the connectivity of single molecules", &
1893 n_keywords=2, n_subsections=0, repeats=.true.)
1894
1895 CALL keyword_create(keyword, __location__, name="NMOL", &
1896 description="number of molecules ", &
1897 usage="NMOL {integer}", default_i_val=1)
1898 CALL section_add_keyword(subsection, keyword)
1899 CALL keyword_release(keyword)
1900
1901 CALL connectivity_framework(subsection, do_conn_psf)
1902 CALL section_add_subsection(section, subsection)
1903 CALL section_release(subsection)
1904
1905 ! MERGE MOLECULES
1906 CALL section_create(subsection, __location__, name="MERGE_MOLECULES", &
1907 description="Enables the creation of connecting bridges (bonds, angles, torsions, impropers)"// &
1908 " between the two or more molecules defined with independent connectivity.", &
1909 n_keywords=2, n_subsections=0, repeats=.false.)
1910
1911 CALL section_create(subsubsection, __location__, name="bonds", &
1912 description="Defines new bonds", n_keywords=2, n_subsections=0, repeats=.false.)
1913 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
1914 description="Two integer indexes per line defining the new bond."// &
1915 " Indexes must be relative to the full system and not to the single molecules", &
1916 repeats=.true., &
1917 usage="{Integer} {Integer}", type_of_var=integer_t, n_var=2)
1918 CALL section_add_keyword(subsubsection, keyword)
1919 CALL keyword_release(keyword)
1920 CALL section_add_subsection(subsection, subsubsection)
1921 CALL section_release(subsubsection)
1922
1923 CALL section_create(subsubsection, __location__, name="angles", &
1924 description="Defines new angles", n_keywords=2, n_subsections=0, &
1925 repeats=.false.)
1926 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
1927 description="Three integer indexes per line defining the new angle"// &
1928 " Indexes must be relative to the full system and not to the single molecules", repeats=.true., &
1929 usage="{Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=3)
1930 CALL section_add_keyword(subsubsection, keyword)
1931 CALL keyword_release(keyword)
1932 CALL section_add_subsection(subsection, subsubsection)
1933 CALL section_release(subsubsection)
1934
1935 CALL section_create(subsubsection, __location__, name="torsions", &
1936 description="Defines new torsions", n_keywords=2, n_subsections=0, &
1937 repeats=.false.)
1938 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
1939 description="Four integer indexes per line defining the new torsion"// &
1940 " Indexes must be relative to the full system and not to the single molecules", repeats=.true., &
1941 usage="{Integer} {Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=4)
1942 CALL section_add_keyword(subsubsection, keyword)
1943 CALL keyword_release(keyword)
1944 CALL section_add_subsection(subsection, subsubsection)
1945 CALL section_release(subsubsection)
1946
1947 CALL section_create(subsubsection, __location__, name="impropers", &
1948 description="Defines new impropers", n_keywords=2, n_subsections=0, &
1949 repeats=.false.)
1950 CALL keyword_create(keyword, __location__, name="_DEFAULT_KEYWORD_", &
1951 description="Four integer indexes per line defining the new improper"// &
1952 " Indexes must be relative to the full system and not to the single molecules", repeats=.true., &
1953 usage="{Integer} {Integer} {Integer} {Integer}", type_of_var=integer_t, n_var=4)
1954 CALL section_add_keyword(subsubsection, keyword)
1955 CALL keyword_release(keyword)
1956 CALL section_add_subsection(subsection, subsubsection)
1957 CALL section_release(subsubsection)
1958
1959 CALL section_add_subsection(section, subsection)
1960 CALL section_release(subsection)
1961
1962 END SUBROUTINE create_molset_section
1963
1964! **************************************************************************************************
1965!> \brief Specify keywords used to generate connectivity
1966!> \param section the section to create
1967!> \author Teodoro Laino [tlaino] - University of Zurich - 08.2008
1968! **************************************************************************************************
1969 SUBROUTINE create_generate_section(section)
1970 TYPE(section_type), POINTER :: section
1971
1972 TYPE(keyword_type), POINTER :: keyword
1973 TYPE(section_type), POINTER :: subsection
1974
1975 cpassert(.NOT. ASSOCIATED(section))
1976 NULLIFY (keyword, subsection)
1977 CALL section_create(section, __location__, name="GENERATE", &
1978 description="Setup of keywords controlling the generation of the connectivity", &
1979 n_keywords=2, n_subsections=0, repeats=.true.)
1980
1981 CALL keyword_create(keyword, __location__, name="REORDER", &
1982 description="Reorder a list of atomic coordinates into order so it can be packed correctly.", &
1983 usage="REORDER <LOGICAL>", &
1984 default_l_val=.false., lone_keyword_l_val=.true.)
1985 CALL section_add_keyword(section, keyword)
1986 CALL keyword_release(keyword)
1987
1988 CALL keyword_create(keyword, __location__, name="CREATE_MOLECULES", &
1989 description="Create molecules names and definition. Can be used to override the"// &
1990 " molecules specifications of a possible input connectivity or to create molecules"// &
1991 " specifications for file types as XYZ, missing of molecules definitions.", &
1992 usage="CREATE_MOLECULES <LOGICAL>", &
1993 default_l_val=.false., lone_keyword_l_val=.true.)
1994 CALL section_add_keyword(section, keyword)
1995 CALL keyword_release(keyword)
1996
1997 CALL keyword_create(keyword, __location__, name="BONDPARM", &
1998 description="Used in conjunction with BONDPARM_FACTOR to "// &
1999 "help determine wheather there is bonding "// &
2000 "between two atoms based on a distance criteria. "// &
2001 "Can use covalent radii information or VDW radii information", &
2002 usage="BONDPARM (COVALENT||VDW)", &
2003 enum_c_vals=s2a("COVALENT", "VDW"), &
2004 enum_i_vals=(/do_bondparm_covalent, do_bondparm_vdw/), &
2005 default_i_val=do_bondparm_covalent)
2006 CALL section_add_keyword(section, keyword)
2007 CALL keyword_release(keyword)
2008
2009 CALL keyword_create(keyword, __location__, name="BONDPARM_FACTOR", &
2010 description="Used in conjunction with BONDPARM to help "// &
2011 "determine wheather there is bonding between "// &
2012 "two atoms based on a distance criteria.", &
2013 usage="bondparm_factor {real}", default_r_val=1.1_dp)
2014 CALL section_add_keyword(section, keyword)
2015 CALL keyword_release(keyword)
2016
2017 CALL keyword_create(keyword, __location__, name="BONDLENGTH_MAX", &
2018 description="Maximum distance to generate neighbor lists to build connectivity", &
2019 usage="BONDLENGTH_MAX <real>", &
2020 default_r_val=cp_unit_to_cp2k(value=3.0_dp, unit_str="angstrom"), &
2021 unit_str="angstrom")
2022 CALL section_add_keyword(section, keyword)
2023 CALL keyword_release(keyword)
2024
2025 CALL keyword_create(keyword, __location__, name="BONDLENGTH_MIN", &
2026 description="Minimum distance to generate neighbor lists to build connectivity", &
2027 usage="BONDLENGTH_MIN <real>", &
2028 default_r_val=cp_unit_to_cp2k(value=0.01_dp, unit_str="angstrom"), &
2029 unit_str="angstrom")
2030 CALL section_add_keyword(section, keyword)
2031 CALL keyword_release(keyword)
2032
2033 ! BONDS
2034 CALL section_create(subsection, __location__, name="BOND", &
2035 description="Section used to add/remove bonds in the connectivity."// &
2036 " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2037 n_keywords=1, n_subsections=0, repeats=.true.)
2038
2039 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
2040 description="controls the activation of the bond", &
2041 usage="&BOND (ADD|REMOVE)", &
2042 enum_c_vals=s2a("ADD", "REMOVE"), &
2043 enum_i_vals=(/do_add, do_remove/), &
2044 default_i_val=do_add)
2045 CALL section_add_keyword(subsection, keyword)
2046 CALL keyword_release(keyword)
2047
2048 CALL keyword_create(keyword, __location__, name="ATOMS", &
2049 description="Specifies two atomic index united by a covalent bond", &
2050 usage="ATOMS {integer} {integer}", type_of_var=integer_t, n_var=2, &
2051 repeats=.true.)
2052 CALL section_add_keyword(subsection, keyword)
2053 CALL keyword_release(keyword)
2054
2055 CALL section_add_subsection(section, subsection)
2056 CALL section_release(subsection)
2057
2058 ! ANGLES
2059 CALL section_create(subsection, __location__, name="ANGLE", &
2060 description="Section used to add/remove angles in the connectivity."// &
2061 " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2062 n_keywords=1, n_subsections=0, repeats=.true.)
2063
2064 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
2065 description="controls the activation of the bond", &
2066 usage="&ANGLE (ADD|REMOVE)", &
2067 enum_c_vals=s2a("ADD", "REMOVE"), &
2068 enum_i_vals=(/do_add, do_remove/), &
2069 default_i_val=do_add)
2070 CALL section_add_keyword(subsection, keyword)
2071 CALL keyword_release(keyword)
2072
2073 CALL keyword_create(keyword, __location__, name="ATOMS", &
2074 description="Specifies two atomic index united by a covalent bond", &
2075 usage="ATOMS {integer} {integer} {integer} ", type_of_var=integer_t, n_var=3, &
2076 repeats=.true.)
2077 CALL section_add_keyword(subsection, keyword)
2078 CALL keyword_release(keyword)
2079
2080 CALL section_add_subsection(section, subsection)
2081 CALL section_release(subsection)
2082
2083 ! TORSIONS
2084 CALL section_create(subsection, __location__, name="TORSION", &
2085 description="Section used to add/remove torsion in the connectivity."// &
2086 " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2087 n_keywords=1, n_subsections=0, repeats=.true.)
2088
2089 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
2090 description="controls the activation of the bond", &
2091 usage="&TORSION (ADD|REMOVE)", &
2092 enum_c_vals=s2a("ADD", "REMOVE"), &
2093 enum_i_vals=(/do_add, do_remove/), &
2094 default_i_val=do_add)
2095 CALL section_add_keyword(subsection, keyword)
2096 CALL keyword_release(keyword)
2097
2098 CALL keyword_create(keyword, __location__, name="ATOMS", &
2099 description="Specifies two atomic index united by a covalent bond", &
2100 usage="ATOMS {integer} {integer} {integer} {integer} ", type_of_var=integer_t, n_var=4, &
2101 repeats=.true.)
2102 CALL section_add_keyword(subsection, keyword)
2103 CALL keyword_release(keyword)
2104
2105 CALL section_add_subsection(section, subsection)
2106 CALL section_release(subsection)
2107
2108 ! IMPROPERS
2109 CALL section_create(subsection, __location__, name="IMPROPER", &
2110 description="Section used to add/remove improper in the connectivity."// &
2111 " Useful for systems with a complex connectivity, difficult to find out automatically.", &
2112 n_keywords=1, n_subsections=0, repeats=.true.)
2113
2114 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
2115 description="controls the activation of the bond", &
2116 usage="&IMPROPER (ADD|REMOVE)", &
2117 enum_c_vals=s2a("ADD", "REMOVE"), &
2118 enum_i_vals=(/do_add, do_remove/), &
2119 default_i_val=do_add)
2120 CALL section_add_keyword(subsection, keyword)
2121 CALL keyword_release(keyword)
2122
2123 CALL keyword_create(keyword, __location__, name="ATOMS", &
2124 description="Specifies two atomic index united by a covalent bond", &
2125 usage="ATOMS {integer} {integer} {integer} {integer} ", type_of_var=integer_t, n_var=4, &
2126 repeats=.true.)
2127 CALL section_add_keyword(subsection, keyword)
2128 CALL keyword_release(keyword)
2129
2130 CALL section_add_subsection(section, subsection)
2131 CALL section_release(subsection)
2132
2133 ! ISOLATED ATOMS
2134 CALL section_create(subsection, __location__, name="ISOLATED_ATOMS", &
2135 description=" This section specifies the atoms that one considers isolated. Useful when present "// &
2136 "ions in solution.", n_keywords=1, n_subsections=0, repeats=.false.)
2137 CALL keyword_create(keyword, __location__, name="LIST", &
2138 description="Specifies a list of atomic indexes of the isolated ion", &
2139 usage="LIST {integer}", type_of_var=integer_t, n_var=-1, &
2140 repeats=.true.)
2141 CALL section_add_keyword(subsection, keyword)
2142 CALL keyword_release(keyword)
2143
2144 CALL section_add_subsection(section, subsection)
2145 CALL section_release(subsection)
2146
2147 ! Neighbor lists keys and printing handling the construction of NL for the connectivity
2148 CALL create_neighbor_lists_section(subsection)
2149 CALL section_add_subsection(section, subsection)
2150 CALL section_release(subsection)
2151
2152 CALL create_gen_print_section(subsection)
2153 CALL section_add_subsection(section, subsection)
2154 CALL section_release(subsection)
2155
2156 END SUBROUTINE create_generate_section
2157
2158! **************************************************************************************************
2159!> \brief Create the print gen section
2160!> \param section the section to create
2161!> \author teo
2162! **************************************************************************************************
2163 SUBROUTINE create_gen_print_section(section)
2164 TYPE(section_type), POINTER :: section
2165
2166 TYPE(section_type), POINTER :: print_key
2167
2168 cpassert(.NOT. ASSOCIATED(section))
2169 CALL section_create(section, __location__, name="print", &
2170 description="Section of possible print options in GENERATE code.", &
2171 n_keywords=0, n_subsections=1, repeats=.false.)
2172
2173 NULLIFY (print_key)
2174 CALL cp_print_key_section_create(print_key, __location__, "NEIGHBOR_LISTS", &
2175 description="Activates the printing of the neighbor lists used"// &
2176 " for generating the connectivity.", print_level=high_print_level, &
2177 filename="", unit_str="angstrom")
2178 CALL section_add_subsection(section, print_key)
2179 CALL section_release(print_key)
2180
2181 CALL cp_print_key_section_create(print_key, __location__, "SUBCELL", &
2182 description="Activates the printing of the subcells used for the "// &
2183 "generation of neighbor lists for connectivity.", &
2184 print_level=high_print_level, filename="__STD_OUT__")
2185 CALL section_add_subsection(section, print_key)
2186 CALL section_release(print_key)
2187
2188 END SUBROUTINE create_gen_print_section
2189
2190! **************************************************************************************************
2191!> \brief Specify keywords used to define connectivity
2192!> \param section the section to create
2193!> \param default ...
2194!> \author teo
2195! **************************************************************************************************
2196 SUBROUTINE connectivity_framework(section, default)
2197 TYPE(section_type), POINTER :: section
2198 INTEGER, INTENT(IN) :: default
2199
2200 TYPE(keyword_type), POINTER :: keyword
2201
2202 cpassert(ASSOCIATED(section))
2203 NULLIFY (keyword)
2204 CALL keyword_create(keyword, __location__, name="CONN_FILE_NAME", &
2205 variants=(/"CONN_FILE"/), &
2206 description="Specifies the filename that contains the molecular connectivity.", &
2207 usage="CONN_FILE_NAME <FILENAME>", type_of_var=lchar_t)
2208 CALL section_add_keyword(section, keyword)
2209 CALL keyword_release(keyword)
2210
2211 CALL keyword_create(keyword, __location__, name="CONN_FILE_FORMAT", &
2212 variants=(/"CONNECTIVITY"/), &
2213 description="Ways to determine and generate a molecules. "// &
2214 "Default is to use GENERATE", &
2215 usage="CONN_FILE_FORMAT (PSF|UPSF|MOL_SET|GENERATE|OFF|G87|G96|AMBER|USER)", &
2216 enum_c_vals=s2a("PSF", "UPSF", "MOL_SET", "GENERATE", "OFF", "G87", "G96", "AMBER", "USER"), &
2217 enum_i_vals=(/do_conn_psf, &
2218 do_conn_psf_u, &
2221 do_conn_off, &
2222 do_conn_g87, &
2223 do_conn_g96, &
2224 do_conn_amb7, &
2225 do_conn_user/), &
2226 enum_desc=s2a("Use a PSF file to determine the connectivity."// &
2227 " (support standard CHARMM/XPLOR and EXT CHARMM)", &
2228 "Read a PSF file in an unformatted way (useful for not so standard PSF).", &
2229 "Use multiple PSF (for now...) files to generate the whole system.", &
2230 "Use a simple distance criteria. (Look at keyword BONDPARM)", &
2231 "Do not generate molecules. (e.g. for QS or ill defined systems)", &
2232 "Use GROMOS G87 topology file.", &
2233 "Use GROMOS G96 topology file.", &
2234 "Use AMBER topology file for reading connectivity (compatible starting from AMBER V.7)", &
2235 "Allows the definition of molecules and residues based on the 5th and 6th column of "// &
2236 "the COORD section. This option can be handy for the definition of molecules with QS "// &
2237 "or to save memory in the case of very large systems (use PARA_RES off)."), &
2238 default_i_val=default)
2239 CALL section_add_keyword(section, keyword)
2240 CALL keyword_release(keyword)
2241 END SUBROUTINE connectivity_framework
2242
2243! **************************************************************************************************
2244!> \brief Create CP2K input section for the DFT+U method parameters
2245!> \param section ...
2246!> \date 01.11.2007
2247!> \author Matthias Krack (MK)
2248!> \version 1.0
2249! **************************************************************************************************
2250 SUBROUTINE create_dft_plus_u_section(section)
2251
2252 TYPE(section_type), POINTER :: section
2253
2254 TYPE(keyword_type), POINTER :: keyword
2255 TYPE(section_type), POINTER :: subsection
2256
2257 cpassert(.NOT. ASSOCIATED(section))
2258
2259 CALL section_create(section, __location__, &
2260 name="DFT_PLUS_U", &
2261 description="Define the parameters for a DFT+U run", &
2262 n_keywords=3, &
2263 n_subsections=1, &
2264 repeats=.false.)
2265 NULLIFY (keyword)
2266
2267 CALL keyword_create(keyword, __location__, &
2268 name="_SECTION_PARAMETERS_", &
2269 description="Controls the activation of the DFT+U section", &
2270 usage="&DFT_PLUS_U ON", &
2271 default_l_val=.false., &
2272 lone_keyword_l_val=.true.)
2273 CALL section_add_keyword(section, keyword)
2274 CALL keyword_release(keyword)
2275
2276 CALL keyword_create(keyword, __location__, &
2277 name="L", &
2278 description="Angular momentum quantum number of the "// &
2279 "orbitals to which the correction is applied", &
2280 repeats=.false., &
2281 n_var=1, &
2282 type_of_var=integer_t, &
2283 default_i_val=-1, &
2284 usage="L 2")
2285 CALL section_add_keyword(section, keyword)
2286 CALL keyword_release(keyword)
2287
2288 CALL keyword_create(keyword, __location__, &
2289 name="U_MINUS_J", &
2290 variants=(/"U_EFF"/), &
2291 description="Effective parameter U(eff) = U - J", &
2292 repeats=.false., &
2293 n_var=1, &
2294 type_of_var=real_t, &
2295 default_r_val=0.0_dp, &
2296 unit_str="au_e", &
2297 usage="U_MINUS_J [eV] 1.4")
2298 CALL section_add_keyword(section, keyword)
2299 CALL keyword_release(keyword)
2300
2301 CALL keyword_create(keyword, __location__, &
2302 name="N", &
2303 description="principal quantum number of the "// &
2304 "orbitals to which the correction is applied. Ignored unless pwdft is used for the calculations", &
2305 repeats=.false., &
2306 n_var=1, &
2307 type_of_var=integer_t, &
2308 default_i_val=-1, &
2309 usage="N 2")
2310 CALL section_add_keyword(section, keyword)
2311 CALL keyword_release(keyword)
2312
2313 CALL keyword_create(keyword, __location__, &
2314 name="U", &
2315 description="U parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2316 repeats=.false., &
2317 n_var=1, &
2318 type_of_var=real_t, &
2319 default_r_val=0.0_dp, &
2320 unit_str="au_e", &
2321 usage="U [eV] 1.4")
2322 CALL section_add_keyword(section, keyword)
2323 CALL keyword_release(keyword)
2324
2325 CALL keyword_create(keyword, __location__, &
2326 name="J", &
2327 description="J parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2328 repeats=.false., &
2329 n_var=1, &
2330 type_of_var=real_t, &
2331 default_r_val=0.0_dp, &
2332 unit_str="au_e", &
2333 usage="J [eV] 1.4")
2334 CALL section_add_keyword(section, keyword)
2335 CALL keyword_release(keyword)
2336
2337 CALL keyword_create(keyword, __location__, &
2338 name="alpha", &
2339 description="alpha parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2340 repeats=.false., &
2341 n_var=1, &
2342 type_of_var=real_t, &
2343 default_r_val=0.0_dp, &
2344 unit_str="au_e", &
2345 usage="alpha [eV] 1.4")
2346 CALL section_add_keyword(section, keyword)
2347 CALL keyword_release(keyword)
2348
2349 CALL keyword_create(keyword, __location__, &
2350 name="beta", &
2351 description="beta parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2352 repeats=.false., &
2353 n_var=1, &
2354 type_of_var=real_t, &
2355 default_r_val=0.0_dp, &
2356 unit_str="au_e", &
2357 usage="beta [eV] 1.4")
2358 CALL section_add_keyword(section, keyword)
2359 CALL keyword_release(keyword)
2360
2361 CALL keyword_create(keyword, __location__, &
2362 name="J0", &
2363 description="J0 parameter in the theory of Dudarev et al. Ignored unless pwdft is used", &
2364 repeats=.false., &
2365 n_var=1, &
2366 type_of_var=real_t, &
2367 default_r_val=0.0_dp, &
2368 unit_str="au_e", &
2369 usage="J0 [eV] 1.4")
2370 CALL section_add_keyword(section, keyword)
2371 CALL keyword_release(keyword)
2372
2373 CALL keyword_create(keyword, __location__, &
2374 name="occupation", &
2375 description="number of electrons in the hubbard shell. Ignored unless pwdft is used", &
2376 repeats=.false., &
2377 n_var=1, &
2378 type_of_var=real_t, &
2379 default_r_val=0.0_dp, &
2380 usage="occupation 6")
2381 CALL section_add_keyword(section, keyword)
2382 CALL keyword_release(keyword)
2383
2384 CALL keyword_create(keyword, __location__, &
2385 name="U_RAMPING", &
2386 description="Increase the effective U parameter stepwise using the specified "// &
2387 "increment until the target value given by U_MINUS_J is reached.", &
2388 repeats=.false., &
2389 n_var=1, &
2390 type_of_var=real_t, &
2391 default_r_val=0.0_dp, &
2392 unit_str="au_e", &
2393 usage="U_RAMPING [eV] 0.1")
2394 CALL section_add_keyword(section, keyword)
2395 CALL keyword_release(keyword)
2396
2397 CALL keyword_create(keyword, __location__, &
2398 name="EPS_U_RAMPING", &
2399 description="Threshold value (SCF convergence) for incrementing the effective "// &
2400 "U value when U ramping is active.", &
2401 repeats=.false., &
2402 n_var=1, &
2403 type_of_var=real_t, &
2404 default_r_val=1.0e-5_dp, &
2405 usage="EPS_U_RAMPING 1.0E-6")
2406 CALL section_add_keyword(section, keyword)
2407 CALL keyword_release(keyword)
2408
2409 CALL keyword_create(keyword, __location__, &
2410 name="INIT_U_RAMPING_EACH_SCF", &
2411 description="Set the initial U ramping value to zero before each wavefunction optimisation. "// &
2412 "The default is to apply U ramping only for the initial wavefunction optimisation.", &
2413 repeats=.false., &
2414 default_l_val=.false., &
2415 lone_keyword_l_val=.true., &
2416 usage="INIT_U_RAMPING_EACH_SCF on")
2417 CALL section_add_keyword(section, keyword)
2418 CALL keyword_release(keyword)
2419
2420 NULLIFY (subsection)
2421
2422 CALL section_create(subsection, __location__, &
2423 name="ENFORCE_OCCUPATION", &
2424 description="Enforce and control a special (initial) orbital occupation. "// &
2425 "Note, this feature works only for the methods MULLIKEN and LOWDIN. "// &
2426 "It should only be used to prepare an initial configuration. An "// &
2427 "inadequate parameter choice can easily inhibit SCF convergence.", &
2428 n_keywords=5, &
2429 n_subsections=0, &
2430 repeats=.false.)
2431
2432 CALL keyword_create(keyword, __location__, &
2433 name="_SECTION_PARAMETERS_", &
2434 description="Controls the activation of the ENFORCE_OCCUPATION section", &
2435 usage="&ENFORCE_OCCUPATION ON", &
2436 default_l_val=.false., &
2437 lone_keyword_l_val=.true.)
2438 CALL section_add_keyword(subsection, keyword)
2439 CALL keyword_release(keyword)
2440
2441 CALL keyword_create(keyword, __location__, name="NELEC", &
2442 variants=(/"N_ELECTRONS"/), &
2443 description="Number of alpha and beta electrons. An occupation (per spin) smaller than 0.5 is ignored.", &
2444 repeats=.false., &
2445 n_var=-1, &
2446 type_of_var=real_t, &
2447 default_r_val=0.0_dp, &
2448 usage="NELEC 5.0 4.0")
2449 CALL section_add_keyword(subsection, keyword)
2450 CALL keyword_release(keyword)
2451
2452 CALL keyword_create(keyword, __location__, &
2453 name="ORBITALS", &
2454 variants=(/"M"/), &
2455 description="Select orbitals and occupation order. An input of 1 to 2*L+1 integer values in "// &
2456 "the range -L to L defining the M values of the spherical orbitals is expected.", &
2457 repeats=.false., &
2458 n_var=-1, &
2459 type_of_var=integer_t, &
2460 default_i_val=0, &
2461 usage="ORBITALS 0 +1 -1")
2462 CALL section_add_keyword(subsection, keyword)
2463 CALL keyword_release(keyword)
2464
2465 CALL keyword_create(keyword, __location__, &
2466 name="EPS_SCF", &
2467 description="The occupation constraint is enforced until this threshold value "// &
2468 "for the SCF convergence criterion is reached", &
2469 repeats=.false., &
2470 n_var=1, &
2471 type_of_var=real_t, &
2472 default_r_val=1.0e30_dp, &
2473 usage="EPS_SCF 0.001")
2474 CALL section_add_keyword(subsection, keyword)
2475 CALL keyword_release(keyword)
2476
2477 CALL keyword_create(keyword, __location__, &
2478 name="MAX_SCF", &
2479 description="The occupation constraint is applied for this number of initial SCF iterations", &
2480 repeats=.false., &
2481 n_var=1, &
2482 type_of_var=integer_t, &
2483 default_i_val=-1, &
2484 usage="MAX_SCF 5")
2485 CALL section_add_keyword(subsection, keyword)
2486 CALL keyword_release(keyword)
2487
2488 CALL keyword_create(keyword, __location__, &
2489 name="SMEAR", &
2490 description="The occupation constraint is applied with smearing", &
2491 repeats=.false., &
2492 default_l_val=.false., &
2493 lone_keyword_l_val=.true., &
2494 usage="SMEAR ON")
2495 CALL section_add_keyword(subsection, keyword)
2496 CALL keyword_release(keyword)
2497
2498 CALL section_add_subsection(section, subsection)
2499 CALL section_release(subsection)
2500
2501 END SUBROUTINE create_dft_plus_u_section
2502
2503END MODULE input_cp2k_subsys
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public guidon2010
integer, save, public goedecker1996
integer, save, public vandevondele2005a
integer, save, public vandevondele2007
integer, save, public hartwigsen1998
integer, save, public krack2005
Handles all functions related to the CELL.
Definition cell_types.F:15
integer, parameter, public use_perd_xyz
Definition cell_types.F:42
integer, parameter, public cell_sym_monoclinic
Definition cell_types.F:29
integer, parameter, public use_perd_y
Definition cell_types.F:42
integer, parameter, public cell_sym_triclinic
Definition cell_types.F:29
integer, parameter, public cell_sym_tetragonal_ab
Definition cell_types.F:29
integer, parameter, public use_perd_xz
Definition cell_types.F:42
integer, parameter, public cell_sym_rhombohedral
Definition cell_types.F:29
integer, parameter, public use_perd_x
Definition cell_types.F:42
integer, parameter, public cell_sym_tetragonal_ac
Definition cell_types.F:29
integer, parameter, public use_perd_z
Definition cell_types.F:42
integer, parameter, public use_perd_yz
Definition cell_types.F:42
integer, parameter, public use_perd_none
Definition cell_types.F:42
integer, parameter, public cell_sym_hexagonal_gamma_60
Definition cell_types.F:29
integer, parameter, public cell_sym_orthorhombic
Definition cell_types.F:29
integer, parameter, public cell_sym_none
Definition cell_types.F:29
integer, parameter, public cell_sym_hexagonal_gamma_120
Definition cell_types.F:29
integer, parameter, public cell_sym_monoclinic_gamma_ab
Definition cell_types.F:29
integer, parameter, public cell_sym_cubic
Definition cell_types.F:29
integer, parameter, public use_perd_xy
Definition cell_types.F:42
integer, parameter, public cell_sym_tetragonal_bc
Definition cell_types.F:29
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
integer, parameter, public medium_print_level
integer, parameter, public high_print_level
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition cp_units.F:1150
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_remove
integer, parameter, public do_conn_g87
integer, parameter, public do_coord_crd
integer, parameter, public do_skip_13
integer, parameter, public do_coord_xtl
integer, parameter, public do_conn_psf_u
integer, parameter, public do_bondparm_covalent
integer, parameter, public dump_pdb
integer, parameter, public do_conn_off
integer, parameter, public do_bondparm_vdw
integer, parameter, public do_conn_user
integer, parameter, public do_coord_cif
integer, parameter, public do_conn_amb7
integer, parameter, public do_cell_cif
integer, parameter, public do_cell_xsc
integer, parameter, public do_conn_psf
integer, parameter, public do_skip_12
integer, parameter, public do_add
integer, parameter, public do_conn_generate
integer, parameter, public do_conn_mol_set
integer, parameter, public do_conn_g96
integer, parameter, public do_coord_pdb
integer, parameter, public do_skip_11
integer, parameter, public do_cell_cp2k
integer, parameter, public do_skip_14
integer, parameter, public gaussian
integer, parameter, public do_coord_xyz
integer, parameter, public do_coord_off
integer, parameter, public do_coord_g96
integer, parameter, public do_coord_cp2k
recursive subroutine, public create_colvar_section(section, skip_recursive_colvar)
creates the colvar section
creates the mm section of the input
subroutine, public create_neighbor_lists_section(section)
This section specifies the input parameters for generation of neighbor lists.
builds the subsystem section of the input
subroutine, public create_subsys_section(section)
creates the structure of a subsys, i.e. a full set of atoms+mol+bounds+cell
subroutine, public create_rng_section(section)
Creates the random number restart section.
subroutine, public create_basis_section(section)
Creates the &BASIS section.
subroutine, public create_structure_data_section(print_key)
creates structure data section for output.. both subsys (for initialization) and motion section....
subroutine, public create_cell_section(section, periodic)
creates the cell section
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations)
creates a list of keywords
subroutine, public section_add_keyword(section, keyword)
adds a keyword to the given section
subroutine, public section_add_subsection(section, subsection)
adds a subsection to the given section
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
a wrapper for basic fortran types.
integer, parameter, public real_t
integer, parameter, public lchar_t
integer, parameter, public char_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public bohr
Definition physcon.F:147
Utilities for string manipulations.
character(len=1), parameter, public newline
represent a keyword in the input
represent a section of the input file