Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \par History
10 : !> - taken out of input_cp2k_motion
11 : !> \author teo & fawzi
12 : ! **************************************************************************************************
13 : MODULE input_cp2k_colvar
14 : USE bibliography, ONLY: Branduardi2007
15 : USE colvar_types, ONLY: &
16 : do_clv_fix_point, do_clv_geo_center, do_clv_x, do_clv_xy, do_clv_xyz, do_clv_xz, do_clv_y, &
17 : do_clv_yz, do_clv_z, plane_def_atoms, plane_def_vec
18 : USE cp_output_handling, ONLY: add_last_numeric,&
19 : cp_print_key_section_create,&
20 : high_print_level,&
21 : low_print_level
22 : USE cp_units, ONLY: cp_unit_to_cp2k
23 : USE input_constants, ONLY: gaussian,&
24 : numerical,&
25 : rmsd_all,&
26 : rmsd_list,&
27 : rmsd_weightlist
28 : USE input_keyword_types, ONLY: keyword_create,&
29 : keyword_release,&
30 : keyword_type
31 : USE input_section_types, ONLY: section_add_keyword,&
32 : section_add_subsection,&
33 : section_create,&
34 : section_release,&
35 : section_type
36 : USE input_val_types, ONLY: char_t,&
37 : integer_t,&
38 : lchar_t,&
39 : real_t
40 : USE kinds, ONLY: dp
41 : USE string_utilities, ONLY: s2a
42 : #include "./base/base_uses.f90"
43 :
44 : IMPLICIT NONE
45 : PRIVATE
46 :
47 : LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .FALSE.
48 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_colvar'
49 :
50 : PUBLIC :: create_colvar_section, &
51 : create_colvar_xyz_d_section, &
52 : create_colvar_xyz_od_section
53 :
54 : CONTAINS
55 :
56 : ! **************************************************************************************************
57 : !> \brief creates the colvar section
58 : !> \param section the section to be created
59 : !> \param skip_recursive_colvar ...
60 : !> \author teo
61 : ! **************************************************************************************************
62 36696 : RECURSIVE SUBROUTINE create_colvar_section(section, skip_recursive_colvar)
63 : TYPE(section_type), POINTER :: section
64 : LOGICAL, OPTIONAL :: skip_recursive_colvar
65 :
66 : LOGICAL :: skip
67 : TYPE(section_type), POINTER :: print_key, subsection
68 :
69 36696 : skip = .FALSE.
70 36696 : IF (PRESENT(skip_recursive_colvar)) skip = skip_recursive_colvar
71 36696 : CPASSERT(.NOT. ASSOCIATED(section))
72 : CALL section_create(section, __LOCATION__, name="COLVAR", &
73 : description="This section specifies the nature of the collective variables.", &
74 36696 : n_keywords=1, n_subsections=1, repeats=.TRUE.)
75 36696 : NULLIFY (subsection, print_key)
76 :
77 : CALL create_colvar_var_section(subsection=subsection, &
78 36696 : section=section, skip_recursive_colvar=skip)
79 :
80 : CALL section_create(subsection, __LOCATION__, name="PRINT", &
81 : description="Controls the printing of the colvar specifications", &
82 36696 : n_keywords=0, n_subsections=1, repeats=.TRUE.)
83 36696 : NULLIFY (print_key)
84 : CALL cp_print_key_section_create(print_key, __LOCATION__, "PROGRAM_RUN_INFO", &
85 : description="Controls the printing of basic information during colvar setup.", &
86 36696 : print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
87 36696 : CALL section_add_subsection(subsection, print_key)
88 36696 : CALL section_release(print_key)
89 36696 : CALL section_add_subsection(section, subsection)
90 36696 : CALL section_release(subsection)
91 :
92 36696 : CALL create_clv_info_section(subsection)
93 36696 : CALL section_add_subsection(section, subsection)
94 36696 : CALL section_release(subsection)
95 :
96 36696 : END SUBROUTINE create_colvar_section
97 :
98 : ! **************************************************************************************************
99 : !> \brief Create the restart section for colvar restraints
100 : !> This section will be only used for restraint restarts.
101 : !> Constraints are handled automatically
102 : !> \param section the section to create
103 : !> \author Teodoro Laino 08.2006
104 : ! **************************************************************************************************
105 36696 : SUBROUTINE create_clv_info_section(section)
106 : TYPE(section_type), POINTER :: section
107 :
108 : TYPE(keyword_type), POINTER :: keyword
109 :
110 36696 : CPASSERT(.NOT. ASSOCIATED(section))
111 36696 : NULLIFY (keyword)
112 : CALL section_create(section, __LOCATION__, name="COLVAR_FUNC_INFO", &
113 : description="Specify further data possibly used by colvars, depending "// &
114 : "on the starting geometry, for computing the functions value.", &
115 36696 : n_subsections=0, repeats=.FALSE.)
116 :
117 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
118 : description="Colvar function data."// &
119 : " The order is an internal order. So if you decide to edit/modify/add these values by hand"// &
120 : " you should know very well what you are doing.!", repeats=.TRUE., &
121 36696 : usage="{Real} ...", type_of_var=real_t, n_var=-1)
122 36696 : CALL section_add_keyword(section, keyword)
123 36696 : CALL keyword_release(keyword)
124 :
125 36696 : END SUBROUTINE create_clv_info_section
126 :
127 : ! **************************************************************************************************
128 : !> \brief creates the collective variables for the colvar section
129 : !> \param subsection ...
130 : !> \param section the section to be created
131 : !> \param skip_recursive_colvar ...
132 : !> \author teo
133 : ! **************************************************************************************************
134 36696 : RECURSIVE SUBROUTINE create_colvar_var_section(subsection, section, skip_recursive_colvar)
135 : TYPE(section_type), POINTER :: subsection, section
136 : LOGICAL, INTENT(IN) :: skip_recursive_colvar
137 :
138 36696 : CPASSERT(.NOT. ASSOCIATED(subsection))
139 36696 : CPASSERT(ASSOCIATED(section))
140 :
141 36696 : CALL create_colvar_dist_section(subsection)
142 36696 : CALL section_add_subsection(section, subsection)
143 36696 : CALL section_release(subsection)
144 :
145 36696 : CALL create_colvar_angle_section(subsection)
146 36696 : CALL section_add_subsection(section, subsection)
147 36696 : CALL section_release(subsection)
148 :
149 36696 : CALL create_colvar_torsion_section(subsection)
150 36696 : CALL section_add_subsection(section, subsection)
151 36696 : CALL section_release(subsection)
152 :
153 36696 : CALL create_colvar_coord_section(subsection)
154 36696 : CALL section_add_subsection(section, subsection)
155 36696 : CALL section_release(subsection)
156 :
157 36696 : CALL create_colvar_pop_section(subsection)
158 36696 : CALL section_add_subsection(section, subsection)
159 36696 : CALL section_release(subsection)
160 :
161 36696 : CALL create_colvar_gyr_section(subsection)
162 36696 : CALL section_add_subsection(section, subsection)
163 36696 : CALL section_release(subsection)
164 :
165 36696 : CALL create_colvar_d_pl_section(subsection)
166 36696 : CALL section_add_subsection(section, subsection)
167 36696 : CALL section_release(subsection)
168 :
169 36696 : CALL create_colvar_a_pl_section(subsection)
170 36696 : CALL section_add_subsection(section, subsection)
171 36696 : CALL section_release(subsection)
172 :
173 36696 : CALL create_colvar_rot_section(subsection)
174 36696 : CALL section_add_subsection(section, subsection)
175 36696 : CALL section_release(subsection)
176 :
177 36696 : CALL create_colvar_dfunct_section(subsection)
178 36696 : CALL section_add_subsection(section, subsection)
179 36696 : CALL section_release(subsection)
180 :
181 36696 : CALL create_colvar_qparm_section(subsection)
182 36696 : CALL section_add_subsection(section, subsection)
183 36696 : CALL section_release(subsection)
184 :
185 36696 : CALL create_colvar_hydronium_shell_section(subsection)
186 36696 : CALL section_add_subsection(section, subsection)
187 36696 : CALL section_release(subsection)
188 :
189 36696 : CALL create_colvar_hydronium_dist_section(subsection)
190 36696 : CALL section_add_subsection(section, subsection)
191 36696 : CALL section_release(subsection)
192 :
193 36696 : CALL create_colvar_acid_hyd_dist_section(subsection)
194 36696 : CALL section_add_subsection(section, subsection)
195 36696 : CALL section_release(subsection)
196 :
197 36696 : CALL create_colvar_acid_hyd_shell_section(subsection)
198 36696 : CALL section_add_subsection(section, subsection)
199 36696 : CALL section_release(subsection)
200 :
201 36696 : CALL create_colvar_rmsd_section(subsection)
202 36696 : CALL section_add_subsection(section, subsection)
203 36696 : CALL section_release(subsection)
204 :
205 36696 : CALL create_colvar_xyz_d_section(subsection)
206 36696 : CALL section_add_subsection(section, subsection)
207 36696 : CALL section_release(subsection)
208 :
209 36696 : CALL create_colvar_xyz_od_section(subsection)
210 36696 : CALL section_add_subsection(section, subsection)
211 36696 : CALL section_release(subsection)
212 :
213 36696 : CALL create_colvar_u_section(subsection)
214 36696 : CALL section_add_subsection(section, subsection)
215 36696 : CALL section_release(subsection)
216 :
217 36696 : CALL create_colvar_wc_section(subsection)
218 36696 : CALL section_add_subsection(section, subsection)
219 36696 : CALL section_release(subsection)
220 :
221 36696 : CALL create_colvar_hbp_section(subsection)
222 36696 : CALL section_add_subsection(section, subsection)
223 36696 : CALL section_release(subsection)
224 :
225 36696 : CALL create_colvar_ring_puckering_section(subsection)
226 36696 : CALL section_add_subsection(section, subsection)
227 36696 : CALL section_release(subsection)
228 :
229 36696 : CALL create_colvar_cond_dist_section(subsection)
230 36696 : CALL section_add_subsection(section, subsection)
231 36696 : CALL section_release(subsection)
232 :
233 36696 : IF (.NOT. skip_recursive_colvar) THEN
234 9174 : CALL create_colvar_rpath_section(subsection)
235 9174 : CALL section_add_subsection(section, subsection)
236 9174 : CALL section_release(subsection)
237 :
238 9174 : CALL create_colvar_dpath_section(subsection)
239 9174 : CALL section_add_subsection(section, subsection)
240 9174 : CALL section_release(subsection)
241 :
242 9174 : CALL create_colvar_comb_section(subsection)
243 9174 : CALL section_add_subsection(section, subsection)
244 9174 : CALL section_release(subsection)
245 : END IF
246 :
247 36696 : END SUBROUTINE create_colvar_var_section
248 :
249 : ! **************************************************************************************************
250 : !> \brief collective variables specifying coordination
251 : !> \param section the section to be created
252 : !> \author teo
253 : ! **************************************************************************************************
254 36696 : SUBROUTINE create_colvar_coord_section(section)
255 : TYPE(section_type), POINTER :: section
256 :
257 : TYPE(keyword_type), POINTER :: keyword
258 : TYPE(section_type), POINTER :: subsection
259 :
260 36696 : CPASSERT(.NOT. ASSOCIATED(section))
261 : CALL section_create(section, __LOCATION__, name="coordination", &
262 : description="Section to define the coordination number as a collective variable.", &
263 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
264 :
265 36696 : NULLIFY (subsection, keyword)
266 :
267 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_FROM", &
268 : variants=(/"POINTS_FROM"/), &
269 : description="Specify indexes of atoms/points building the coordination variable. ", &
270 : usage="ATOMS_FROM {integer} {integer} ..", repeats=.TRUE., &
271 73392 : n_var=-1, type_of_var=integer_t)
272 36696 : CALL section_add_keyword(section, keyword)
273 36696 : CALL keyword_release(keyword)
274 :
275 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_TO", &
276 : variants=(/"POINTS_TO"/), &
277 : description="Specify indexes of atoms/points building the coordination variable. ", &
278 : usage="ATOMS_TO {integer} {integer} ..", repeats=.TRUE., &
279 73392 : n_var=-1, type_of_var=integer_t)
280 36696 : CALL section_add_keyword(section, keyword)
281 36696 : CALL keyword_release(keyword)
282 :
283 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_TO_B", &
284 : variants=(/"POINTS_TO_B"/), &
285 : description="For the CV given by the multiplication of two coorination numbers,"// &
286 : " here specify indexes of the third set of atoms/points. ", &
287 : usage="ATOMS_TO_B {integer} {integer} ..", repeats=.TRUE., &
288 73392 : n_var=-1, type_of_var=integer_t)
289 36696 : CALL section_add_keyword(section, keyword)
290 36696 : CALL keyword_release(keyword)
291 :
292 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_FROM", &
293 : description="Specify alternatively kinds of atoms building the coordination variable.", &
294 : usage="KINDS_FROM {CHAR} {CHAR} ..", repeats=.TRUE., &
295 36696 : n_var=-1, type_of_var=char_t)
296 36696 : CALL section_add_keyword(section, keyword)
297 36696 : CALL keyword_release(keyword)
298 :
299 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_TO", &
300 : description="Specify alternatively kinds of atoms building the coordination variable.", &
301 : usage="KINDS_TO {CHAR} {CHAR} ..", repeats=.TRUE., &
302 36696 : n_var=-1, type_of_var=char_t)
303 36696 : CALL section_add_keyword(section, keyword)
304 36696 : CALL keyword_release(keyword)
305 :
306 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_TO_B", &
307 : description="For the CV given by the multiplication of two coorination numbers,"// &
308 : " here specify alternatively kinds of atoms building the coordination variable.", &
309 : usage="KINDS_TO_B {CHAR} {CHAR} ..", repeats=.TRUE., &
310 36696 : n_var=-1, type_of_var=char_t)
311 36696 : CALL section_add_keyword(section, keyword)
312 36696 : CALL keyword_release(keyword)
313 :
314 : ! Must be present in each colvar and handled properly
315 36696 : CALL create_point_section(subsection)
316 36696 : CALL section_add_subsection(section, subsection)
317 36696 : CALL section_release(subsection)
318 :
319 : CALL keyword_create(keyword, __LOCATION__, name="R0", &
320 : variants=(/"R_0"/), &
321 : description="Specify the R0 parameter in the coordination function.", &
322 : usage="R0 {real}", default_r_val=3.0_dp, &
323 73392 : unit_str="bohr", n_var=1)
324 36696 : CALL section_add_keyword(section, keyword)
325 36696 : CALL keyword_release(keyword)
326 :
327 : CALL keyword_create(keyword, __LOCATION__, name="NN", &
328 : variants=(/"EXPON_NUMERATOR"/), &
329 : description="Sets the value of the numerator of the exponential factor"// &
330 : " in the coordination FUNCTION.", &
331 : usage="NN {integer}", default_i_val=6, &
332 73392 : n_var=1)
333 36696 : CALL section_add_keyword(section, keyword)
334 36696 : CALL keyword_release(keyword)
335 :
336 : CALL keyword_create(keyword, __LOCATION__, name="ND", &
337 : variants=(/"EXPON_DENOMINATOR"/), &
338 : description="Sets the value of the denominator of the exponential factor"// &
339 : " in the coordination FUNCTION.", &
340 : usage="ND {integer}", default_i_val=12, &
341 73392 : n_var=1)
342 36696 : CALL section_add_keyword(section, keyword)
343 36696 : CALL keyword_release(keyword)
344 :
345 : CALL keyword_create(keyword, __LOCATION__, name="R0_B", &
346 : variants=(/"R_0_B"/), &
347 : description="For the CV given by the multiplication of two coorination numbers,"// &
348 : " specify the R0 parameter in the second coordination function.", &
349 : usage="R0_B {real}", default_r_val=3.0_dp, &
350 73392 : unit_str="bohr", n_var=1)
351 36696 : CALL section_add_keyword(section, keyword)
352 36696 : CALL keyword_release(keyword)
353 :
354 : CALL keyword_create(keyword, __LOCATION__, name="NN_B", &
355 : variants=(/"EXPON_NUMERATOR_B"/), &
356 : description="For the CV given by the multiplication of two coorination numbers,"// &
357 : " Sets the value of the numerator of the exponential factor"// &
358 : " in the coordination FUNCTION.", &
359 : usage="NN_B {integer}", default_i_val=6, &
360 73392 : n_var=1)
361 36696 : CALL section_add_keyword(section, keyword)
362 36696 : CALL keyword_release(keyword)
363 :
364 : CALL keyword_create(keyword, __LOCATION__, name="ND_B", &
365 : variants=(/"EXPON_DENOMINATOR_B"/), &
366 : description="For the CV given by the multiplication of two coorination numbers,"// &
367 : " Sets the value of the denominator of the exponential factor"// &
368 : " in the coordination FUNCTION.", &
369 : usage="ND_B {integer}", default_i_val=12, &
370 73392 : n_var=1)
371 36696 : CALL section_add_keyword(section, keyword)
372 36696 : CALL keyword_release(keyword)
373 :
374 36696 : END SUBROUTINE create_colvar_coord_section
375 :
376 : ! **************************************************************************************************
377 : !> \brief ...
378 : !> \param section ...
379 : ! **************************************************************************************************
380 36696 : SUBROUTINE create_colvar_cond_dist_section(section)
381 : TYPE(section_type), POINTER :: section
382 :
383 : TYPE(keyword_type), POINTER :: keyword
384 : TYPE(section_type), POINTER :: subsection
385 :
386 36696 : CPASSERT(.NOT. ASSOCIATED(section))
387 : CALL section_create(section, __LOCATION__, name="CONDITIONED_DISTANCE", &
388 : description="Section to define the conditioned distance as a collective variable.", &
389 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
390 :
391 36696 : NULLIFY (subsection, keyword)
392 :
393 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_DISTANCE", &
394 : description="Specify indexes of atoms/points from which the distance is computed. ", &
395 : usage="ATOMS_DISTANCE {integer} {integer} ..", repeats=.TRUE., &
396 36696 : n_var=-1, type_of_var=integer_t)
397 36696 : CALL section_add_keyword(section, keyword)
398 36696 : CALL keyword_release(keyword)
399 :
400 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_FROM", &
401 : variants=(/"POINTS_FROM"/), &
402 : description="Specify indexes of atoms/points building the coordination variable. ", &
403 : usage="ATOMS_FROM {integer} {integer} ..", repeats=.TRUE., &
404 73392 : n_var=-1, type_of_var=integer_t)
405 36696 : CALL section_add_keyword(section, keyword)
406 36696 : CALL keyword_release(keyword)
407 :
408 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_TO", &
409 : variants=(/"POINTS_TO"/), &
410 : description="Specify indexes of atoms/points building the coordination variable. ", &
411 : usage="ATOMS_TO {integer} {integer} ..", repeats=.TRUE., &
412 73392 : n_var=-1, type_of_var=integer_t)
413 36696 : CALL section_add_keyword(section, keyword)
414 36696 : CALL keyword_release(keyword)
415 :
416 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_FROM", &
417 : description="Specify alternatively kinds of atoms building the coordination variable.", &
418 : usage="KINDS_FROM {CHAR} {CHAR} ..", repeats=.TRUE., &
419 36696 : n_var=-1, type_of_var=char_t)
420 36696 : CALL section_add_keyword(section, keyword)
421 36696 : CALL keyword_release(keyword)
422 :
423 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_TO", &
424 : description="Specify alternatively kinds of atoms building the coordination variable.", &
425 : usage="KINDS_TO {CHAR} {CHAR} ..", repeats=.TRUE., &
426 36696 : n_var=-1, type_of_var=char_t)
427 36696 : CALL section_add_keyword(section, keyword)
428 36696 : CALL keyword_release(keyword)
429 :
430 : ! Must be present in each colvar and handled properly
431 36696 : CALL create_point_section(subsection)
432 36696 : CALL section_add_subsection(section, subsection)
433 36696 : CALL section_release(subsection)
434 :
435 : CALL keyword_create(keyword, __LOCATION__, name="R0", &
436 : variants=(/"R_0"/), &
437 : description="Specify the R0 parameter in the coordination function.", &
438 : usage="R0 {real}", default_r_val=3.0_dp, &
439 73392 : unit_str="bohr", n_var=1)
440 36696 : CALL section_add_keyword(section, keyword)
441 36696 : CALL keyword_release(keyword)
442 :
443 : CALL keyword_create(keyword, __LOCATION__, name="NN", &
444 : variants=(/"EXPON_NUMERATOR"/), &
445 : description="Sets the value of the numerator of the exponential factor"// &
446 : " in the coordination FUNCTION.", &
447 : usage="NN {integer}", default_i_val=6, &
448 73392 : n_var=1)
449 36696 : CALL section_add_keyword(section, keyword)
450 36696 : CALL keyword_release(keyword)
451 :
452 : CALL keyword_create(keyword, __LOCATION__, name="ND", &
453 : variants=(/"EXPON_DENOMINATOR"/), &
454 : description="Sets the value of the denominator of the exponential factor"// &
455 : " in the coordination FUNCTION.", &
456 : usage="ND {integer}", default_i_val=12, &
457 73392 : n_var=1)
458 36696 : CALL section_add_keyword(section, keyword)
459 36696 : CALL keyword_release(keyword)
460 :
461 : CALL keyword_create(keyword, __LOCATION__, name="LAMBDA", &
462 : description="Specify the lambda parameter at the exponent of the conditioned distance function.", &
463 : usage="R0 {real}", default_r_val=3.0_dp, &
464 36696 : unit_str="bohr", n_var=1)
465 36696 : CALL section_add_keyword(section, keyword)
466 36696 : CALL keyword_release(keyword)
467 :
468 36696 : END SUBROUTINE create_colvar_cond_dist_section
469 :
470 : ! **************************************************************************************************
471 : !> \brief collective variables specifying population of a specie based on coordination
472 : !> \param section the section to be created
473 : !> \date 01.2009
474 : !> \author Fabio Sterpone
475 : ! **************************************************************************************************
476 36696 : SUBROUTINE create_colvar_pop_section(section)
477 : TYPE(section_type), POINTER :: section
478 :
479 : TYPE(keyword_type), POINTER :: keyword
480 : TYPE(section_type), POINTER :: subsection
481 :
482 36696 : CPASSERT(.NOT. ASSOCIATED(section))
483 : CALL section_create(section, __LOCATION__, name="population", &
484 : description="Section to define the population of specie as a collective variable. "// &
485 : "See also <https://doi.org/10.1021/jp3019588>.", &
486 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
487 :
488 36696 : NULLIFY (subsection, keyword)
489 :
490 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_FROM", &
491 : variants=(/"POINTS_FROM"/), &
492 : description="Specify indexes of atoms/points building the coordination variable. ", &
493 : usage="ATOMS_FROM {integer} {integer} ..", repeats=.TRUE., &
494 73392 : n_var=-1, type_of_var=integer_t)
495 36696 : CALL section_add_keyword(section, keyword)
496 36696 : CALL keyword_release(keyword)
497 :
498 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_TO", &
499 : variants=(/"POINTS_TO"/), &
500 : description="Specify indexes of atoms/points building the coordination variable. ", &
501 : usage="ATOMS_TO {integer} {integer} ..", repeats=.TRUE., &
502 73392 : n_var=-1, type_of_var=integer_t)
503 36696 : CALL section_add_keyword(section, keyword)
504 36696 : CALL keyword_release(keyword)
505 :
506 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_FROM", &
507 : description="Specify alternatively kinds of atoms building the coordination variable.", &
508 : usage="KINDS_FROM {CHAR} {CHAR} ..", repeats=.TRUE., &
509 36696 : n_var=-1, type_of_var=char_t)
510 36696 : CALL section_add_keyword(section, keyword)
511 36696 : CALL keyword_release(keyword)
512 :
513 : CALL keyword_create(keyword, __LOCATION__, name="KINDS_TO", &
514 : description="Specify alternatively kinds of atoms building the coordination variable.", &
515 : usage="KINDS_TO {CHAR} {CHAR} ..", repeats=.TRUE., &
516 36696 : n_var=-1, type_of_var=char_t)
517 36696 : CALL section_add_keyword(section, keyword)
518 36696 : CALL keyword_release(keyword)
519 :
520 : ! Must be present in each colvar and handled properly
521 36696 : CALL create_point_section(subsection)
522 36696 : CALL section_add_subsection(section, subsection)
523 36696 : CALL section_release(subsection)
524 :
525 : CALL keyword_create(keyword, __LOCATION__, name="R0", &
526 : variants=(/"R_0"/), &
527 : description="Specify the R0 parameter in the coordination function.", &
528 : usage="R0 {real}", default_r_val=3.0_dp, &
529 73392 : n_var=1)
530 36696 : CALL section_add_keyword(section, keyword)
531 36696 : CALL keyword_release(keyword)
532 :
533 : CALL keyword_create(keyword, __LOCATION__, name="NN", &
534 : variants=(/"EXPON_NUMERATOR"/), &
535 : description="Sets the value of the numerator of the exponential factor"// &
536 : " in the coordination FUNCTION.", &
537 : usage="NN {integer}", default_i_val=6, &
538 73392 : n_var=1)
539 36696 : CALL section_add_keyword(section, keyword)
540 36696 : CALL keyword_release(keyword)
541 :
542 : CALL keyword_create(keyword, __LOCATION__, name="ND", &
543 : variants=(/"EXPON_DENOMINATOR"/), &
544 : description="Sets the value of the denominator of the exponential factor"// &
545 : " in the coordination FUNCTION.", &
546 : usage="ND {integer}", default_i_val=12, &
547 73392 : n_var=1)
548 36696 : CALL section_add_keyword(section, keyword)
549 36696 : CALL keyword_release(keyword)
550 :
551 : CALL keyword_create(keyword, __LOCATION__, name="n0", &
552 : variants=(/"n_0"/), &
553 : description="Specify the n0 parameter that sets the coordination of the species.", &
554 : usage="n0 {integer}", default_i_val=4, &
555 73392 : n_var=1)
556 36696 : CALL section_add_keyword(section, keyword)
557 36696 : CALL keyword_release(keyword)
558 :
559 : CALL keyword_create(keyword, __LOCATION__, name="SIGMA", &
560 : description="Specify the gaussian width of used to build the population istogram.", &
561 : usage="SIGMA {real}", default_r_val=0.5_dp, &
562 36696 : n_var=1)
563 36696 : CALL section_add_keyword(section, keyword)
564 36696 : CALL keyword_release(keyword)
565 :
566 36696 : END SUBROUTINE create_colvar_pop_section
567 :
568 : ! **************************************************************************************************
569 : !> \brief ...
570 : !> \param section ...
571 : ! **************************************************************************************************
572 36696 : SUBROUTINE create_colvar_gyr_section(section)
573 : TYPE(section_type), POINTER :: section
574 :
575 : TYPE(keyword_type), POINTER :: keyword
576 : TYPE(section_type), POINTER :: subsection
577 :
578 36696 : CPASSERT(.NOT. ASSOCIATED(section))
579 : CALL section_create(section, __LOCATION__, name="GYRATION_RADIUS", &
580 : description="Section to define the gyration radius as a collective variable.", &
581 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
582 :
583 36696 : NULLIFY (subsection, keyword)
584 :
585 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
586 : variants=(/"POINTS"/), &
587 : description="Specify indexes of atoms/points defyining the gyration radius variable. ", &
588 : usage="ATOMS {integer} {integer} ..", repeats=.TRUE., &
589 73392 : n_var=-1, type_of_var=integer_t)
590 36696 : CALL section_add_keyword(section, keyword)
591 36696 : CALL keyword_release(keyword)
592 :
593 : CALL keyword_create(keyword, __LOCATION__, name="KINDS", &
594 : description="Specify alternatively kinds of atoms defining the gyration radius.", &
595 : usage="KINDS {CHAR} {CHAR} ..", repeats=.TRUE., &
596 36696 : n_var=-1, type_of_var=char_t)
597 36696 : CALL section_add_keyword(section, keyword)
598 36696 : CALL keyword_release(keyword)
599 :
600 : ! Must be present in each colvar and handled properly
601 36696 : CALL create_point_section(subsection)
602 36696 : CALL section_add_subsection(section, subsection)
603 36696 : CALL section_release(subsection)
604 :
605 36696 : END SUBROUTINE create_colvar_gyr_section
606 :
607 : ! **************************************************************************************************
608 : !> \brief collective variables specifying torsion
609 : !> \param section the section to be created
610 : !> \author teo
611 : ! **************************************************************************************************
612 36696 : SUBROUTINE create_colvar_dfunct_section(section)
613 : TYPE(section_type), POINTER :: section
614 :
615 : TYPE(keyword_type), POINTER :: keyword
616 : TYPE(section_type), POINTER :: subsection
617 :
618 36696 : CPASSERT(.NOT. ASSOCIATED(section))
619 : CALL section_create(section, __LOCATION__, name="DISTANCE_FUNCTION", &
620 : description="Section to define functions between two distances as collective variables."// &
621 : " The function is defined as d1+coeff*d2", &
622 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
623 :
624 36696 : NULLIFY (keyword, subsection)
625 :
626 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
627 : variants=(/"POINTS"/), &
628 : description="Specifies the indexes of atoms/points for the two bonds d1=(1-2) d2=(3-4).", &
629 : usage="ATOMS {integer} {integer} {integer} {integer}", &
630 73392 : n_var=4, type_of_var=integer_t)
631 36696 : CALL section_add_keyword(section, keyword)
632 36696 : CALL keyword_release(keyword)
633 :
634 : CALL keyword_create(keyword, __LOCATION__, name="COEFFICIENT", &
635 : description="Specifies the coefficient in the function for the constraint."// &
636 : " -1.0 has to be used for distance difference, 1.0 for distance addition", &
637 : usage="COEFFICIENT {real}", &
638 36696 : type_of_var=real_t)
639 36696 : CALL section_add_keyword(section, keyword)
640 36696 : CALL keyword_release(keyword)
641 :
642 : CALL keyword_create(keyword, __LOCATION__, name="PBC", &
643 : description="Whether periodic boundary conditions should be applied on the "// &
644 : "atomic position before computing the colvar or not.", &
645 : usage="PBC", &
646 36696 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
647 36696 : CALL section_add_keyword(section, keyword)
648 36696 : CALL keyword_release(keyword)
649 :
650 : ! Must be present in each colvar and handled properly
651 36696 : CALL create_point_section(subsection)
652 36696 : CALL section_add_subsection(section, subsection)
653 36696 : CALL section_release(subsection)
654 :
655 36696 : END SUBROUTINE create_colvar_dfunct_section
656 :
657 : ! **************************************************************************************************
658 : !> \brief collective variables specifying torsion
659 : !> \param section the section to be created
660 : !> \author teo
661 : ! **************************************************************************************************
662 36696 : SUBROUTINE create_colvar_torsion_section(section)
663 : TYPE(section_type), POINTER :: section
664 :
665 : TYPE(keyword_type), POINTER :: keyword
666 : TYPE(section_type), POINTER :: subsection
667 :
668 36696 : CPASSERT(.NOT. ASSOCIATED(section))
669 : CALL section_create(section, __LOCATION__, name="torsion", &
670 : description="Section to define the torsion as a collective variables.", &
671 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
672 :
673 36696 : NULLIFY (keyword, subsection)
674 :
675 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
676 : variants=(/"POINTS"/), &
677 : description="Specifies the indexes of atoms/points defining the torsion.", &
678 : usage="ATOMS {integer} {integer} {integer} {integer}", &
679 73392 : n_var=4, type_of_var=integer_t)
680 36696 : CALL section_add_keyword(section, keyword)
681 36696 : CALL keyword_release(keyword)
682 :
683 : ! Must be present in each colvar and handled properly
684 36696 : CALL create_point_section(subsection)
685 36696 : CALL section_add_subsection(section, subsection)
686 36696 : CALL section_release(subsection)
687 :
688 36696 : END SUBROUTINE create_colvar_torsion_section
689 :
690 : ! **************************************************************************************************
691 : !> \brief collective variables specifying torsion
692 : !> \param section the section to be created
693 : !> \author teo
694 : ! **************************************************************************************************
695 36696 : SUBROUTINE create_colvar_rot_section(section)
696 : TYPE(section_type), POINTER :: section
697 :
698 : TYPE(keyword_type), POINTER :: keyword
699 : TYPE(section_type), POINTER :: subsection
700 :
701 36696 : CPASSERT(.NOT. ASSOCIATED(section))
702 : CALL section_create(section, __LOCATION__, name="bond_rotation", &
703 : description="Section to define the rotation of a bond/line with respect to"// &
704 : " another bond/line", &
705 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
706 :
707 36696 : NULLIFY (keyword, subsection)
708 :
709 : CALL keyword_create(keyword, __LOCATION__, name="P1_BOND1", &
710 : description="Specifies the index of atom/point defining the first point"// &
711 : " of the first bond/line.", &
712 : usage="P1_BOND1 {integer}", &
713 36696 : n_var=1, type_of_var=integer_t)
714 36696 : CALL section_add_keyword(section, keyword)
715 36696 : CALL keyword_release(keyword)
716 :
717 : CALL keyword_create(keyword, __LOCATION__, name="P2_BOND1", &
718 : description="Specifies the index of atom/point defining the second point"// &
719 : " of the first bond/line.", &
720 : usage="P2_BOND1 {integer}", &
721 36696 : n_var=1, type_of_var=integer_t)
722 36696 : CALL section_add_keyword(section, keyword)
723 36696 : CALL keyword_release(keyword)
724 :
725 : CALL keyword_create(keyword, __LOCATION__, name="P1_BOND2", &
726 : description="Specifies the index of atom/point defining the first point"// &
727 : " of the second bond/line.", &
728 : usage="P1_BOND2 {integer}", &
729 36696 : n_var=1, type_of_var=integer_t)
730 36696 : CALL section_add_keyword(section, keyword)
731 36696 : CALL keyword_release(keyword)
732 :
733 : CALL keyword_create(keyword, __LOCATION__, name="P2_BOND2", &
734 : description="Specifies the index of atom/point defining the second point"// &
735 : " of the second bond/line.", &
736 : usage="P2_BOND2 {integer}", &
737 36696 : n_var=1, type_of_var=integer_t)
738 36696 : CALL section_add_keyword(section, keyword)
739 36696 : CALL keyword_release(keyword)
740 :
741 : ! Must be present in each colvar and handled properly
742 36696 : CALL create_point_section(subsection)
743 36696 : CALL section_add_subsection(section, subsection)
744 36696 : CALL section_release(subsection)
745 :
746 36696 : END SUBROUTINE create_colvar_rot_section
747 :
748 : ! **************************************************************************************************
749 : !> \brief collective variables specifying angles
750 : !> \param section the section to be created
751 : !> \author teo
752 : ! **************************************************************************************************
753 36696 : SUBROUTINE create_colvar_angle_section(section)
754 : TYPE(section_type), POINTER :: section
755 :
756 : TYPE(keyword_type), POINTER :: keyword
757 : TYPE(section_type), POINTER :: subsection
758 :
759 36696 : CPASSERT(.NOT. ASSOCIATED(section))
760 : CALL section_create(section, __LOCATION__, name="angle", &
761 : description="Section to define the angle as a collective variables.", &
762 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
763 36696 : NULLIFY (keyword, subsection)
764 :
765 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
766 : variants=(/"POINTS"/), &
767 : description="Specifies the indexes of atoms/points defining the angle.", &
768 : usage="ATOMS {integer} {integer} {integer}", &
769 73392 : n_var=3, type_of_var=integer_t)
770 36696 : CALL section_add_keyword(section, keyword)
771 36696 : CALL keyword_release(keyword)
772 :
773 : ! Must be present in each colvar and handled properly
774 36696 : CALL create_point_section(subsection)
775 36696 : CALL section_add_subsection(section, subsection)
776 36696 : CALL section_release(subsection)
777 :
778 36696 : END SUBROUTINE create_colvar_angle_section
779 :
780 : ! **************************************************************************************************
781 : !> \brief creates the colvar section regarded to the collective variables dist
782 : !> \param section the section to be created
783 : !> \author teo
784 : ! **************************************************************************************************
785 36696 : SUBROUTINE create_colvar_dist_section(section)
786 : TYPE(section_type), POINTER :: section
787 :
788 : TYPE(keyword_type), POINTER :: keyword
789 : TYPE(section_type), POINTER :: subsection
790 :
791 36696 : CPASSERT(.NOT. ASSOCIATED(section))
792 : CALL section_create(section, __LOCATION__, name="distance", &
793 : description="Section to define the distance as a collective variables.", &
794 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
795 36696 : NULLIFY (keyword, subsection)
796 :
797 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
798 : variants=(/"POINTS"/), &
799 : description="Specifies the indexes of atoms/points defining the distance.", &
800 : usage="ATOMS {integer} {integer}", &
801 73392 : n_var=2, type_of_var=integer_t)
802 36696 : CALL section_add_keyword(section, keyword)
803 36696 : CALL keyword_release(keyword)
804 : CALL keyword_create(keyword, __LOCATION__, name="AXIS", &
805 : description="Define the axes along which the colvar should be evaluated", &
806 : usage="AXIS (XYZ | X | Y | Z | XY| XZ | YZ)", &
807 : enum_c_vals=s2a("XYZ", "X", "Y", "Z", "XY", "XZ", "YZ"), &
808 : enum_i_vals=(/do_clv_xyz, do_clv_x, do_clv_y, do_clv_z, do_clv_xy, do_clv_xz, do_clv_yz/), &
809 36696 : default_i_val=do_clv_xyz)
810 36696 : CALL section_add_keyword(section, keyword)
811 36696 : CALL keyword_release(keyword)
812 :
813 : CALL keyword_create(keyword, __LOCATION__, name="SIGN", &
814 : description="Whether the distance along one Cartesian axis has to be considered with sign."// &
815 : " This option is valid if only one dimension is selected.", &
816 : usage="SIGN", &
817 36696 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
818 36696 : CALL section_add_keyword(section, keyword)
819 36696 : CALL keyword_release(keyword)
820 :
821 : ! Must be present in each colvar and handled properly
822 36696 : CALL create_point_section(subsection)
823 36696 : CALL section_add_subsection(section, subsection)
824 36696 : CALL section_release(subsection)
825 :
826 36696 : END SUBROUTINE create_colvar_dist_section
827 :
828 : ! **************************************************************************************************
829 : !> \brief creates the colvar section regarded to the collective variables dist
830 : !> \param section the section to be created
831 : !> \author teo
832 : ! **************************************************************************************************
833 36699 : SUBROUTINE create_colvar_xyz_d_section(section)
834 : TYPE(section_type), POINTER :: section
835 :
836 : TYPE(keyword_type), POINTER :: keyword
837 : TYPE(section_type), POINTER :: subsection
838 :
839 36699 : CPASSERT(.NOT. ASSOCIATED(section))
840 : CALL section_create(section, __LOCATION__, name="XYZ_DIAG", &
841 : description="Section to define the distance of an atom from its starting "// &
842 : "position ((X-X(0))^2+(Y-Y(0))^2+(Z-Z(0))^2) or part of its components as a collective variable. "// &
843 : "If absolute_position is specified, instead the CV is represented by the "// &
844 : "instantaneous position of the atom (only available for X, Y or Z components).", &
845 36699 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
846 36699 : NULLIFY (keyword, subsection)
847 :
848 : CALL keyword_create(keyword, __LOCATION__, name="ATOM", &
849 : variants=(/"POINT"/), &
850 : description="Specifies the index of the atom/point.", &
851 : usage="ATOM {integer}", &
852 73398 : n_var=1, type_of_var=integer_t)
853 36699 : CALL section_add_keyword(section, keyword)
854 36699 : CALL keyword_release(keyword)
855 :
856 : CALL keyword_create(keyword, __LOCATION__, name="COMPONENT", &
857 : description="Define the component of the position vector which will be used "// &
858 : "as a colvar.", &
859 : usage="AXIS (XYZ | X | Y | Z | XY| XZ | YZ)", &
860 : enum_c_vals=s2a("XYZ", "X", "Y", "Z", "XY", "XZ", "YZ"), &
861 : enum_i_vals=(/do_clv_xyz, do_clv_x, do_clv_y, do_clv_z, do_clv_xy, do_clv_xz, do_clv_yz/), &
862 36699 : default_i_val=do_clv_xyz)
863 36699 : CALL section_add_keyword(section, keyword)
864 36699 : CALL keyword_release(keyword)
865 :
866 : CALL keyword_create(keyword, __LOCATION__, name="PBC", &
867 : description="Whether periodic boundary conditions should be applied on the "// &
868 : "atomic position before computing the colvar or not.", &
869 : usage="PBC", &
870 36699 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
871 36699 : CALL section_add_keyword(section, keyword)
872 36699 : CALL keyword_release(keyword)
873 :
874 : CALL keyword_create(keyword, __LOCATION__, name="ABSOLUTE_POSITION", &
875 : description="If enabled, the absolute position of the atoms will be used. ", &
876 : usage="ABSOLUTE_POSITION", &
877 36699 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
878 36699 : CALL section_add_keyword(section, keyword)
879 36699 : CALL keyword_release(keyword)
880 :
881 : ! Must be present in each colvar and handled properly
882 36699 : CALL create_point_section(subsection)
883 36699 : CALL section_add_subsection(section, subsection)
884 36699 : CALL section_release(subsection)
885 :
886 36699 : END SUBROUTINE create_colvar_xyz_d_section
887 :
888 : ! **************************************************************************************************
889 : !> \brief creates the colvar section regarded to the collective variables dist
890 : !> \param section the section to be created
891 : !> \author teo
892 : ! **************************************************************************************************
893 36699 : SUBROUTINE create_colvar_xyz_od_section(section)
894 : TYPE(section_type), POINTER :: section
895 :
896 : TYPE(keyword_type), POINTER :: keyword
897 : TYPE(section_type), POINTER :: subsection
898 :
899 36699 : CPASSERT(.NOT. ASSOCIATED(section))
900 : CALL section_create(section, __LOCATION__, name="XYZ_OUTERDIAG", &
901 : description="Section to define the cross term (XA-XA(0))*(XB-XB(0))+(XA-XA(0))*(YB-YB(0))"// &
902 : " or part of its components as a collective variable. The final term is given by the product"// &
903 : " of the components of A with the components of B.", &
904 36699 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
905 36699 : NULLIFY (keyword, subsection)
906 :
907 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
908 : variants=(/"POINTS"/), &
909 : description="Specifies the index of the atoms/points A and B.", &
910 : usage="ATOMS {integer} {integer}", &
911 73398 : n_var=2, type_of_var=integer_t)
912 36699 : CALL section_add_keyword(section, keyword)
913 36699 : CALL keyword_release(keyword)
914 :
915 : CALL keyword_create(keyword, __LOCATION__, name="COMPONENT_A", &
916 : description="Define the component of the position vector which will be used "// &
917 : "as a colvar for atom A.", &
918 : usage="AXIS (XYZ | X | Y | Z | XY| XZ | YZ)", &
919 : enum_c_vals=s2a("XYZ", "X", "Y", "Z", "XY", "XZ", "YZ"), &
920 : enum_i_vals=(/do_clv_xyz, do_clv_x, do_clv_y, do_clv_z, do_clv_xy, do_clv_xz, do_clv_yz/), &
921 36699 : default_i_val=do_clv_xyz)
922 36699 : CALL section_add_keyword(section, keyword)
923 36699 : CALL keyword_release(keyword)
924 :
925 : CALL keyword_create(keyword, __LOCATION__, name="COMPONENT_B", &
926 : description="Define the component of the position vector which will be used "// &
927 : "as a colvar for atom B.", &
928 : usage="AXIS (XYZ | X | Y | Z | XY| XZ | YZ)", &
929 : enum_c_vals=s2a("XYZ", "X", "Y", "Z", "XY", "XZ", "YZ"), &
930 : enum_i_vals=(/do_clv_xyz, do_clv_x, do_clv_y, do_clv_z, do_clv_xy, do_clv_xz, do_clv_yz/), &
931 36699 : default_i_val=do_clv_xyz)
932 36699 : CALL section_add_keyword(section, keyword)
933 36699 : CALL keyword_release(keyword)
934 :
935 : CALL keyword_create(keyword, __LOCATION__, name="PBC", &
936 : description="Whether periodic boundary conditions should be applied on the "// &
937 : "atomic position before computing the colvar or not.", &
938 : usage="PBC", &
939 36699 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
940 36699 : CALL section_add_keyword(section, keyword)
941 36699 : CALL keyword_release(keyword)
942 :
943 : ! Must be present in each colvar and handled properly
944 36699 : CALL create_point_section(subsection)
945 36699 : CALL section_add_subsection(section, subsection)
946 36699 : CALL section_release(subsection)
947 :
948 36699 : END SUBROUTINE create_colvar_xyz_od_section
949 :
950 : ! **************************************************************************************************
951 : !> \brief energy as collective variable
952 : !> \param section the section to be created
953 : !> \author Sebastiano Caravati
954 : ! **************************************************************************************************
955 36696 : SUBROUTINE create_colvar_u_section(section)
956 : TYPE(section_type), POINTER :: section
957 :
958 : TYPE(keyword_type), POINTER :: keyword
959 : TYPE(section_type), POINTER :: subsection
960 :
961 36696 : CPASSERT(.NOT. ASSOCIATED(section))
962 : CALL section_create(section, __LOCATION__, name="u", &
963 : description="Section to define the energy as a generalized collective variable.", &
964 36696 : n_keywords=0, n_subsections=0, repeats=.FALSE.)
965 :
966 36696 : NULLIFY (subsection, keyword)
967 : CALL section_create(subsection, __LOCATION__, name="MIXED", &
968 : description="This section allows to use any function of the energy subsystems"// &
969 : " in a mixed_env calculation as a collective variable.", &
970 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
971 :
972 : CALL keyword_create(keyword, __LOCATION__, name="ENERGY_FUNCTION", &
973 : description="Specifies the functional form of the collective variable in mathematical notation.", &
974 : usage="ENERGY_FUNCTION (E1+E2-LOG(E1/E2))", type_of_var=lchar_t, &
975 36696 : n_var=1)
976 36696 : CALL section_add_keyword(subsection, keyword)
977 36696 : CALL keyword_release(keyword)
978 :
979 : CALL keyword_create(keyword, __LOCATION__, name="VARIABLES", &
980 : description="Defines the variables of the functional form. To allow an efficient"// &
981 : " mapping the order of the energy variables will be considered identical to the"// &
982 : " order of the force_eval in the force_eval_order list.", &
983 : usage="VARIABLES x", type_of_var=char_t, &
984 36696 : n_var=-1)
985 36696 : CALL section_add_keyword(subsection, keyword)
986 36696 : CALL keyword_release(keyword)
987 :
988 : CALL keyword_create(keyword, __LOCATION__, name="PARAMETERS", &
989 : description="Defines the parameters of the functional form", &
990 : usage="PARAMETERS a b D", type_of_var=char_t, &
991 36696 : n_var=-1, repeats=.TRUE.)
992 36696 : CALL section_add_keyword(subsection, keyword)
993 36696 : CALL keyword_release(keyword)
994 :
995 : CALL keyword_create(keyword, __LOCATION__, name="VALUES", &
996 : description="Defines the values of parameter of the functional form", &
997 : usage="VALUES ", type_of_var=real_t, &
998 36696 : n_var=-1, repeats=.TRUE., unit_str="internal_cp2k")
999 36696 : CALL section_add_keyword(subsection, keyword)
1000 36696 : CALL keyword_release(keyword)
1001 :
1002 : CALL keyword_create(keyword, __LOCATION__, name="UNITS", &
1003 : description="Optionally, allows to define valid CP2K unit strings for each parameter value. "// &
1004 : "It is assumed that the corresponding parameter value is specified in this unit.", &
1005 : usage="UNITS angstrom eV*angstrom^-1 angstrom^1 K", type_of_var=char_t, &
1006 36696 : n_var=-1, repeats=.TRUE.)
1007 36696 : CALL section_add_keyword(subsection, keyword)
1008 36696 : CALL keyword_release(keyword)
1009 :
1010 : CALL keyword_create(keyword, __LOCATION__, name="DX", &
1011 : description="Parameter used for computing the derivative with the Ridders' method.", &
1012 36696 : usage="DX <REAL>", default_r_val=0.1_dp, unit_str="bohr")
1013 36696 : CALL section_add_keyword(subsection, keyword)
1014 36696 : CALL keyword_release(keyword)
1015 :
1016 : CALL keyword_create(keyword, __LOCATION__, name="ERROR_LIMIT", &
1017 : description="Checks that the error in computing the derivative is not larger than "// &
1018 : "the value set; in case error is larger a warning message is printed.", &
1019 36696 : usage="ERROR_LIMIT <REAL>", default_r_val=1.0E-12_dp)
1020 36696 : CALL section_add_keyword(subsection, keyword)
1021 36696 : CALL keyword_release(keyword)
1022 :
1023 36696 : CALL section_add_subsection(section, subsection)
1024 36696 : CALL section_release(subsection)
1025 :
1026 36696 : END SUBROUTINE create_colvar_u_section
1027 :
1028 : ! **************************************************************************************************
1029 : !> \brief creates the colvar section regarded to the collective variables distance
1030 : !> of a point from a plane
1031 : !> \param section the section to be created
1032 : !> \author teo
1033 : ! **************************************************************************************************
1034 36696 : SUBROUTINE create_colvar_d_pl_section(section)
1035 : TYPE(section_type), POINTER :: section
1036 :
1037 : TYPE(keyword_type), POINTER :: keyword
1038 : TYPE(section_type), POINTER :: subsection
1039 :
1040 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1041 : CALL section_create(section, __LOCATION__, name="distance_point_plane", &
1042 : description="Section to define the distance of a point from a plane "// &
1043 : "as a collective variables.", &
1044 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1045 36696 : NULLIFY (keyword, subsection)
1046 :
1047 : CALL keyword_create(keyword, __LOCATION__, name="PBC", &
1048 : description="Whether periodic boundary conditions should be applied on the "// &
1049 : "atomic position before computing the colvar or not.", &
1050 : usage="PBC", &
1051 36696 : default_l_val=.TRUE., lone_keyword_l_val=.TRUE.)
1052 36696 : CALL section_add_keyword(section, keyword)
1053 36696 : CALL keyword_release(keyword)
1054 :
1055 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_PLANE", &
1056 : variants=(/"POINTS_PLANE"/), &
1057 : description="Specifies the indexes of atoms/points defining the plane.", &
1058 : usage="ATOMS_PLANE <INTEGER> <INTEGER> <INTEGER>", &
1059 73392 : n_var=3, type_of_var=integer_t)
1060 36696 : CALL section_add_keyword(section, keyword)
1061 36696 : CALL keyword_release(keyword)
1062 :
1063 : CALL keyword_create(keyword, __LOCATION__, name="ATOM_POINT", &
1064 : variants=(/"POINT_POINT"/), &
1065 : description="Specifies the atom/point index defining the point.", &
1066 : usage="ATOM_POINT <INTEGER>", &
1067 73392 : n_var=1, type_of_var=integer_t)
1068 36696 : CALL section_add_keyword(section, keyword)
1069 36696 : CALL keyword_release(keyword)
1070 :
1071 : ! Must be present in each colvar and handled properly
1072 36696 : CALL create_point_section(subsection)
1073 36696 : CALL section_add_subsection(section, subsection)
1074 36696 : CALL section_release(subsection)
1075 :
1076 36696 : END SUBROUTINE create_colvar_d_pl_section
1077 :
1078 : ! **************************************************************************************************
1079 : !> \brief creates the colvar section regarded to the collective variables
1080 : !> angles betweem two planes
1081 : !> \param section the section to be created
1082 : !> \author teo
1083 : ! **************************************************************************************************
1084 36696 : SUBROUTINE create_colvar_a_pl_section(section)
1085 : TYPE(section_type), POINTER :: section
1086 :
1087 : TYPE(keyword_type), POINTER :: keyword
1088 : TYPE(section_type), POINTER :: subsection
1089 :
1090 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1091 : CALL section_create(section, __LOCATION__, name="angle_plane_plane", &
1092 : description="This section defines the angle between two planes "// &
1093 : "as a collective variables.", &
1094 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1095 36696 : NULLIFY (keyword, subsection)
1096 :
1097 : CALL section_create(subsection, __LOCATION__, name="PLANE", &
1098 : description="This section defines the plane. When using this colvar, "// &
1099 : "two plane section must be defined!", &
1100 36696 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
1101 :
1102 : CALL keyword_create(keyword, __LOCATION__, name="DEF_TYPE", &
1103 : description="Specify how the plane is defined: either by 3 atoms or by a fixed normal "// &
1104 : "vector. At least one plane must be defined through atoms.", &
1105 : usage="DEF_TYPE ATOMS", &
1106 : default_i_val=plane_def_atoms, &
1107 : enum_c_vals=s2a("ATOMS", "VECTOR"), &
1108 : enum_desc=s2a("Plane defined by the position of 3 atoms", &
1109 : "Plane defined by a fixed normal vector"), &
1110 36696 : enum_i_vals=(/plane_def_atoms, plane_def_vec/))
1111 36696 : CALL section_add_keyword(subsection, keyword)
1112 36696 : CALL keyword_release(keyword)
1113 :
1114 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
1115 : description="Specifies the indexes of 3 atoms/points defining the plane.", &
1116 : usage="ATOMS <INTEGER> <INTEGER> <INTEGER>", &
1117 36696 : n_var=3, type_of_var=integer_t)
1118 36696 : CALL section_add_keyword(subsection, keyword)
1119 36696 : CALL keyword_release(keyword)
1120 :
1121 : CALL keyword_create(keyword, __LOCATION__, name="NORMAL_VECTOR", &
1122 : description="Alternatively to 3 atoms/points one can define one of the two, "// &
1123 : "planes by defining its NORMAL vector.", &
1124 : usage="NORMAL_VECTOR 0.0 1.0 0.0", &
1125 36696 : n_var=3, type_of_var=real_t)
1126 36696 : CALL section_add_keyword(subsection, keyword)
1127 36696 : CALL keyword_release(keyword)
1128 36696 : CALL section_add_subsection(section, subsection)
1129 36696 : CALL section_release(subsection)
1130 :
1131 : ! Must be present in each colvar and handled properly
1132 36696 : CALL create_point_section(subsection)
1133 36696 : CALL section_add_subsection(section, subsection)
1134 36696 : CALL section_release(subsection)
1135 36696 : END SUBROUTINE create_colvar_a_pl_section
1136 :
1137 : ! **************************************************************************************************
1138 : !> \brief create a geometrical point as a function of several atom coordinates
1139 : !> \param section the section to be created
1140 : !> \author teo
1141 : ! **************************************************************************************************
1142 623838 : SUBROUTINE create_point_section(section)
1143 : TYPE(section_type), POINTER :: section
1144 :
1145 : TYPE(keyword_type), POINTER :: keyword
1146 :
1147 623838 : CPASSERT(.NOT. ASSOCIATED(section))
1148 : CALL section_create(section, __LOCATION__, name="POINT", &
1149 : description="Enables the possibility to use geometrical centers instead of single atoms"// &
1150 : " to define colvars", &
1151 623838 : n_keywords=1, n_subsections=0, repeats=.TRUE.)
1152 :
1153 623838 : NULLIFY (keyword)
1154 :
1155 : CALL keyword_create(keyword, __LOCATION__, name="TYPE", &
1156 : description="Chooses the type of geometrical point", &
1157 : usage="type (GEO_CENTER|FIX_POINT)", &
1158 : enum_c_vals=s2a("GEO_CENTER", "FIX_POINT"), &
1159 : enum_desc=s2a("Computes the geometrical center of the listed atoms", &
1160 : "Defines a fixed point in space"), &
1161 : enum_i_vals=(/do_clv_geo_center, do_clv_fix_point/), &
1162 623838 : default_i_val=do_clv_geo_center)
1163 623838 : CALL section_add_keyword(section, keyword)
1164 623838 : CALL keyword_release(keyword)
1165 :
1166 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
1167 : description="Specifies the indexes of atoms defining the geometrical center", &
1168 : usage="ATOMS {integer} {integer} {integer} {integer}", &
1169 623838 : n_var=-1, type_of_var=integer_t, repeats=.TRUE.)
1170 623838 : CALL section_add_keyword(section, keyword)
1171 623838 : CALL keyword_release(keyword)
1172 :
1173 : CALL keyword_create( &
1174 : keyword, __LOCATION__, name="WEIGHTS", &
1175 : description="Specifies the weights for a weighted geometrical center. Default is 1/natoms for every atom", &
1176 : usage="WEIGHTS {real} {real} {real} {real}", &
1177 623838 : n_var=-1, type_of_var=real_t, repeats=.TRUE.)
1178 623838 : CALL section_add_keyword(section, keyword)
1179 623838 : CALL keyword_release(keyword)
1180 :
1181 : CALL keyword_create(keyword, __LOCATION__, name="XYZ", &
1182 : description="Specifies the xyz of the fixed point (if the case)", &
1183 : usage="XYZ {real} {real} {real}", &
1184 : n_var=3, type_of_var=real_t, unit_str="bohr", &
1185 623838 : repeats=.FALSE.)
1186 623838 : CALL section_add_keyword(section, keyword)
1187 623838 : CALL keyword_release(keyword)
1188 :
1189 623838 : END SUBROUTINE create_point_section
1190 :
1191 : ! **************************************************************************************************
1192 : !> \brief collective variables specifying torsion
1193 : !> \param section the section to be created
1194 : !> \author teo
1195 : ! **************************************************************************************************
1196 36696 : SUBROUTINE create_colvar_qparm_section(section)
1197 : TYPE(section_type), POINTER :: section
1198 :
1199 : TYPE(keyword_type), POINTER :: keyword
1200 : TYPE(section_type), POINTER :: subsection
1201 :
1202 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1203 : CALL section_create(section, __LOCATION__, name="qparm", &
1204 : description="Section to define the Q parameter (crystalline order parameter) as a collective variable.", &
1205 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1206 :
1207 36696 : NULLIFY (keyword, subsection)
1208 :
1209 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_FROM", &
1210 : variants=(/"POINTS_FROM"/), &
1211 : description="Specify indexes of atoms/points building the coordination variable. ", &
1212 : usage="ATOMS_FROM {integer} {integer} ..", repeats=.TRUE., &
1213 73392 : n_var=-1, type_of_var=integer_t)
1214 36696 : CALL section_add_keyword(section, keyword)
1215 36696 : CALL keyword_release(keyword)
1216 :
1217 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS_TO", &
1218 : variants=(/"POINTS_TO"/), &
1219 : description="Specify indexes of atoms/points building the coordination variable. ", &
1220 : usage="ATOMS_TO {integer} {integer} ..", repeats=.TRUE., &
1221 73392 : n_var=-1, type_of_var=integer_t)
1222 36696 : CALL section_add_keyword(section, keyword)
1223 36696 : CALL keyword_release(keyword)
1224 :
1225 : CALL keyword_create(keyword, __LOCATION__, name="RCUT", &
1226 : description="Specifies the distance cutoff for neighbors. "// &
1227 : "Cutoff function is exactly zero for all neighbors beyond RCUT.", &
1228 : usage="RCUT {real}", &
1229 36696 : n_var=1, unit_str="angstrom", type_of_var=real_t)
1230 36696 : CALL section_add_keyword(section, keyword)
1231 36696 : CALL keyword_release(keyword)
1232 :
1233 : CALL keyword_create(keyword, __LOCATION__, name="INCLUDE_IMAGES", &
1234 : description="Whether to include periodic images of ATOMS_TO into the neighbor list.", &
1235 : usage="INCLUDE_IMAGES", &
1236 36696 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1237 36696 : CALL section_add_keyword(section, keyword)
1238 36696 : CALL keyword_release(keyword)
1239 :
1240 : CALL keyword_create(keyword, __LOCATION__, name="RSTART", &
1241 : description="Specifies the distance cutoff for neighbors. "// &
1242 : "Cutoff function is exactly 1 for all neighbors closer than RSTART.", &
1243 : usage="RSTART {real}", &
1244 36696 : n_var=1, unit_str="angstrom", type_of_var=real_t)
1245 36696 : CALL section_add_keyword(section, keyword)
1246 36696 : CALL keyword_release(keyword)
1247 :
1248 : CALL keyword_create(keyword, __LOCATION__, name="L", &
1249 : description="Specifies the L spherical harmonics from Ylm.", &
1250 : usage="L {integer}", &
1251 36696 : n_var=1, type_of_var=integer_t)
1252 36696 : CALL section_add_keyword(section, keyword)
1253 36696 : CALL keyword_release(keyword)
1254 :
1255 : !CALL keyword_create(keyword, __LOCATION__, name="ALPHA", &
1256 : ! description="Specifies the width of the Fermi-Dirac style smearing around RCUT.", &
1257 : ! usage="ALPHA {real}", unit_str="angstrom^-1", default_r_val=0.0_dp)
1258 : !CALL section_add_keyword(section, keyword)
1259 : !CALL keyword_release(keyword)
1260 :
1261 : ! Must be present in each colvar and handled properly
1262 36696 : CALL create_point_section(subsection)
1263 36696 : CALL section_add_subsection(section, subsection)
1264 36696 : CALL section_release(subsection)
1265 :
1266 36696 : END SUBROUTINE create_colvar_qparm_section
1267 :
1268 : ! **************************************************************************************************
1269 : !> \brief collective variables specifying hydronium solvation
1270 : !> \param section the section to be created
1271 : !> \author Marcel Baer
1272 : ! **************************************************************************************************
1273 36696 : SUBROUTINE create_colvar_hydronium_shell_section(section)
1274 : TYPE(section_type), POINTER :: section
1275 :
1276 : TYPE(keyword_type), POINTER :: keyword
1277 :
1278 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1279 : CALL section_create(section, __LOCATION__, name="HYDRONIUM_SHELL", &
1280 : description="Section to define the formation of a hydronium as a"// &
1281 : " collective variable. Number of oxygens in the 1st shell of the"// &
1282 : " hydronium. Adapted from Equation (3) in Supplementary Info of"// &
1283 : " J. Am. Chem. Soc.,128, 2006, 11318, i.e. omitting the cutoff function"// &
1284 : " and summing only over the oxygens of water.", &
1285 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1286 :
1287 36696 : NULLIFY (keyword)
1288 :
1289 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS", &
1290 : description="Specifies indexes of atoms building the coordination variable."// &
1291 : " Oxygens of the water molecules.", &
1292 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1293 36696 : n_var=-1, type_of_var=integer_t)
1294 36696 : CALL section_add_keyword(section, keyword)
1295 36696 : CALL keyword_release(keyword)
1296 :
1297 : CALL keyword_create(keyword, __LOCATION__, name="HYDROGENS", &
1298 : description="Specifies indexes of atoms building the coordination variable."// &
1299 : " Hydrogens of the water molecules.", &
1300 : usage="HYDROGENS {integer} {integer} ..", repeats=.TRUE., &
1301 36696 : n_var=-1, type_of_var=integer_t)
1302 36696 : CALL section_add_keyword(section, keyword)
1303 36696 : CALL keyword_release(keyword)
1304 :
1305 : CALL keyword_create(keyword, __LOCATION__, name="ROO", &
1306 : description="Specifies the rc parameter in the coordination function:"// &
1307 : " number of oxygens per water oxygen.", &
1308 : usage="ROO {real}", default_r_val=cp_unit_to_cp2k(value=3.0_dp, &
1309 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1310 36696 : CALL section_add_keyword(section, keyword)
1311 36696 : CALL keyword_release(keyword)
1312 :
1313 : CALL keyword_create(keyword, __LOCATION__, name="pOO", &
1314 : variants=(/"EXPON_NUMERATORA"/), &
1315 : description="Sets the value of the numerator of the exponential factor"// &
1316 : " in the coordination function: number of oxygens per water oxygen.", &
1317 : usage="pOO {integer}", default_i_val=6, &
1318 73392 : n_var=1)
1319 36696 : CALL section_add_keyword(section, keyword)
1320 36696 : CALL keyword_release(keyword)
1321 :
1322 : CALL keyword_create(keyword, __LOCATION__, name="qOO", &
1323 : variants=(/"EXPON_DENOMINATORA"/), &
1324 : description="Sets the value of the denominator of the exponential factor"// &
1325 : " in the coordination function: number of oxygens per water oxygen.", &
1326 : usage="qOO {integer}", default_i_val=12, &
1327 73392 : n_var=1)
1328 36696 : CALL section_add_keyword(section, keyword)
1329 36696 : CALL keyword_release(keyword)
1330 :
1331 : CALL keyword_create(keyword, __LOCATION__, name="ROH", &
1332 : description="Specifies the rc parameter in the coordination function:"// &
1333 : " number of hydrogens per water molecule.", &
1334 : usage="ROH {real}", default_r_val=cp_unit_to_cp2k(value=3.0_dp, &
1335 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1336 36696 : CALL section_add_keyword(section, keyword)
1337 36696 : CALL keyword_release(keyword)
1338 :
1339 : CALL keyword_create(keyword, __LOCATION__, name="pOH", &
1340 : variants=(/"EXPON_NUMERATORB"/), &
1341 : description="Sets the value of the numerator of the exponential factor"// &
1342 : " in the coordination function: number of hydrogens per water molecule.", &
1343 : usage="pOH {integer}", default_i_val=6, &
1344 73392 : n_var=1)
1345 36696 : CALL section_add_keyword(section, keyword)
1346 36696 : CALL keyword_release(keyword)
1347 :
1348 : CALL keyword_create(keyword, __LOCATION__, name="qOH", &
1349 : variants=(/"EXPON_DENOMINATORB"/), &
1350 : description="Sets the value of the denominator of the exponential factor"// &
1351 : " in the coordination function: number of hydrogens per water molecule.", &
1352 : usage="qOH {integer}", default_i_val=12, &
1353 73392 : n_var=1)
1354 36696 : CALL section_add_keyword(section, keyword)
1355 36696 : CALL keyword_release(keyword)
1356 :
1357 : CALL keyword_create(keyword, __LOCATION__, name="NH", &
1358 : description="Specifies the NH parameter in the M function.", &
1359 : usage="NH {real}", default_r_val=3.0_dp, &
1360 36696 : n_var=1)
1361 36696 : CALL section_add_keyword(section, keyword)
1362 36696 : CALL keyword_release(keyword)
1363 :
1364 : CALL keyword_create(keyword, __LOCATION__, name="pM", &
1365 : variants=(/"EXPON_NUMERATOR"/), &
1366 : description="Sets the value of the numerator of the exponential factor"// &
1367 : " in the M function.", &
1368 : usage="pM {integer}", default_i_val=8, &
1369 73392 : n_var=1)
1370 36696 : CALL section_add_keyword(section, keyword)
1371 36696 : CALL keyword_release(keyword)
1372 :
1373 : CALL keyword_create(keyword, __LOCATION__, name="qM", &
1374 : variants=(/"EXPON_DENOMINATOR"/), &
1375 : description="Sets the value of the denominator of the exponential factor"// &
1376 : " in the M function.", &
1377 : usage="qM {integer}", default_i_val=16, &
1378 73392 : n_var=1)
1379 36696 : CALL section_add_keyword(section, keyword)
1380 36696 : CALL keyword_release(keyword)
1381 :
1382 : CALL keyword_create(keyword, __LOCATION__, name="LAMBDA", &
1383 : description="Specify the LAMBDA parameter in the hydronium function.", &
1384 : usage="LAMBDA {real}", default_r_val=10.0_dp, &
1385 36696 : n_var=1)
1386 36696 : CALL section_add_keyword(section, keyword)
1387 36696 : CALL keyword_release(keyword)
1388 :
1389 36696 : END SUBROUTINE create_colvar_hydronium_shell_section
1390 :
1391 : ! **************************************************************************************************
1392 : !> \brief collective variables specifying the distance between hydronium and hydroxide ion
1393 : !> \param section the section to be created
1394 : !> \author Dorothea Golze
1395 : ! **************************************************************************************************
1396 36696 : SUBROUTINE create_colvar_hydronium_dist_section(section)
1397 : TYPE(section_type), POINTER :: section
1398 :
1399 : TYPE(keyword_type), POINTER :: keyword
1400 :
1401 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1402 : CALL section_create(section, __LOCATION__, name="HYDRONIUM_DISTANCE", &
1403 : description="Section to define the formation of a hydronium as a"// &
1404 : " collective variable. Distance between hydronium and hydroxide ion"// &
1405 : " Experimental at this point, i.e. not proved to be an effective"// &
1406 : " collective variable.", &
1407 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1408 :
1409 36696 : NULLIFY (keyword)
1410 :
1411 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS", &
1412 : description="Specifies indexes of atoms building the coordination variable."// &
1413 : " Oxygens of the water molecules.", &
1414 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1415 36696 : n_var=-1, type_of_var=integer_t)
1416 36696 : CALL section_add_keyword(section, keyword)
1417 36696 : CALL keyword_release(keyword)
1418 :
1419 : CALL keyword_create(keyword, __LOCATION__, name="HYDROGENS", &
1420 : description="Specifies indexes of atoms building the coordination variable."// &
1421 : " Hydrogens of the water molecules.", &
1422 : usage="HYDROGENS {integer} {integer} ..", repeats=.TRUE., &
1423 36696 : n_var=-1, type_of_var=integer_t)
1424 36696 : CALL section_add_keyword(section, keyword)
1425 36696 : CALL keyword_release(keyword)
1426 :
1427 : CALL keyword_create(keyword, __LOCATION__, name="ROH", &
1428 : description="Specifies the rc parameter in the coordination function:"// &
1429 : " number of hydrogens per water molecule.", &
1430 : usage="ROH {real}", default_r_val=cp_unit_to_cp2k(value=2.4_dp, &
1431 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1432 36696 : CALL section_add_keyword(section, keyword)
1433 36696 : CALL keyword_release(keyword)
1434 :
1435 : CALL keyword_create(keyword, __LOCATION__, name="pOH", &
1436 : description="Sets the value of the numerator of the exponential factor"// &
1437 : " in the coordination function: number of hydrogens per water molecule.", &
1438 : usage="pOH {integer}", default_i_val=6, &
1439 36696 : n_var=1)
1440 36696 : CALL section_add_keyword(section, keyword)
1441 36696 : CALL keyword_release(keyword)
1442 :
1443 : CALL keyword_create(keyword, __LOCATION__, name="qOH", &
1444 : description="Sets the value of the denominator of the exponential factor"// &
1445 : " in the coordination function: number of hydrogens per water molecule.", &
1446 : usage="qOH {integer}", default_i_val=12, &
1447 36696 : n_var=1)
1448 36696 : CALL section_add_keyword(section, keyword)
1449 36696 : CALL keyword_release(keyword)
1450 :
1451 : CALL keyword_create(keyword, __LOCATION__, name="NH", &
1452 : description="Specifies the NH parameter in the M function.", &
1453 : usage="NH {real}", default_r_val=2.2_dp, &
1454 36696 : n_var=1)
1455 36696 : CALL section_add_keyword(section, keyword)
1456 36696 : CALL keyword_release(keyword)
1457 :
1458 : CALL keyword_create(keyword, __LOCATION__, name="pM", &
1459 : description="Sets the value of the numerator of the exponential factor"// &
1460 : " in the M function.", &
1461 : usage="pM {integer}", default_i_val=8, &
1462 36696 : n_var=1)
1463 36696 : CALL section_add_keyword(section, keyword)
1464 36696 : CALL keyword_release(keyword)
1465 :
1466 : CALL keyword_create(keyword, __LOCATION__, name="qM", &
1467 : description="Sets the value of the denominator of the exponential factor"// &
1468 : " in the M function.", &
1469 : usage="qM {integer}", default_i_val=16, &
1470 36696 : n_var=1)
1471 36696 : CALL section_add_keyword(section, keyword)
1472 36696 : CALL keyword_release(keyword)
1473 :
1474 : CALL keyword_create(keyword, __LOCATION__, name="NN", &
1475 : description="Specifies the NN parameter in the F function.", &
1476 : usage="NN {real}", default_r_val=1.5_dp, &
1477 36696 : n_var=1)
1478 36696 : CALL section_add_keyword(section, keyword)
1479 36696 : CALL keyword_release(keyword)
1480 :
1481 : CALL keyword_create(keyword, __LOCATION__, name="pF", &
1482 : description="Sets the value of the numerator of the exponential factor"// &
1483 : " in the F function.", &
1484 : usage="pF {integer}", default_i_val=8, &
1485 36696 : n_var=1)
1486 36696 : CALL section_add_keyword(section, keyword)
1487 36696 : CALL keyword_release(keyword)
1488 :
1489 : CALL keyword_create(keyword, __LOCATION__, name="qF", &
1490 : description="Sets the value of the denominator of the exponential factor"// &
1491 : " in the F function.", &
1492 : usage="qF {integer}", default_i_val=16, &
1493 36696 : n_var=1)
1494 36696 : CALL section_add_keyword(section, keyword)
1495 36696 : CALL keyword_release(keyword)
1496 :
1497 : CALL keyword_create(keyword, __LOCATION__, name="LAMBDA", &
1498 : description="Specify the LAMBDA parameter in the hydronium function.", &
1499 : usage="LAMBDA {real}", default_r_val=20.0_dp, &
1500 36696 : n_var=1)
1501 36696 : CALL section_add_keyword(section, keyword)
1502 36696 : CALL keyword_release(keyword)
1503 :
1504 36696 : END SUBROUTINE create_colvar_hydronium_dist_section
1505 :
1506 : ! **************************************************************************************************
1507 : !> \brief collective variables specifying the solvation of carboxylic acid;
1508 : !> distance between hydronium ion and acetate ion; Equation (2) in
1509 : !> Supplementary Information of J. Am. Chem. Soc.,128, 2006, 11318
1510 : !> \param section the section to be created
1511 : !> \author Dorothea Golze
1512 : ! **************************************************************************************************
1513 36696 : SUBROUTINE create_colvar_acid_hyd_dist_section(section)
1514 : TYPE(section_type), POINTER :: section
1515 :
1516 : TYPE(keyword_type), POINTER :: keyword
1517 :
1518 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1519 : CALL section_create(section, __LOCATION__, name="ACID_HYDRONIUM_DISTANCE", &
1520 : description="Section to define the dissociation of a carboxylic acid in"// &
1521 : " water. Distance between hydronium ion and acetate ion. Equation (2)"// &
1522 : " in Supplementary Info of J. Am. Chem. Soc.,128, 2006, 11318.", &
1523 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1524 :
1525 36696 : NULLIFY (keyword)
1526 :
1527 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS_WATER", &
1528 : description="Specifies indexes of atoms building the coordination variable."// &
1529 : " Oxygens of the water molecules. ", &
1530 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1531 36696 : n_var=-1, type_of_var=integer_t)
1532 36696 : CALL section_add_keyword(section, keyword)
1533 36696 : CALL keyword_release(keyword)
1534 :
1535 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS_ACID", &
1536 : description="Specifies indexes of atoms building the coordination variable."// &
1537 : " Oxygens of the carboxyl groups.", &
1538 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1539 36696 : n_var=-1, type_of_var=integer_t)
1540 36696 : CALL section_add_keyword(section, keyword)
1541 36696 : CALL keyword_release(keyword)
1542 :
1543 : CALL keyword_create(keyword, __LOCATION__, name="HYDROGENS", &
1544 : description="Specifies indexes of atoms building the coordination variable."// &
1545 : " Hydrogens of the water molecules and of the carboxyl groups.", &
1546 : usage="HYDROGENS {integer} {integer} ..", repeats=.TRUE., &
1547 36696 : n_var=-1, type_of_var=integer_t)
1548 36696 : CALL section_add_keyword(section, keyword)
1549 36696 : CALL keyword_release(keyword)
1550 :
1551 : CALL keyword_create(keyword, __LOCATION__, name="pWOH", &
1552 : description="Sets the value of the numerator of the exponential factor"// &
1553 : " in the coordination function: number of hydrogens per water molecule.", &
1554 : usage="pWOH {integer}", default_i_val=8, &
1555 36696 : n_var=1)
1556 36696 : CALL section_add_keyword(section, keyword)
1557 36696 : CALL keyword_release(keyword)
1558 :
1559 : CALL keyword_create(keyword, __LOCATION__, name="qWOH", &
1560 : description="Sets the value of the denominator of the exponential factor"// &
1561 : " in the coordination function: number of hydrogens per water molecule.", &
1562 : usage="qWOH {integer}", default_i_val=16, &
1563 36696 : n_var=1)
1564 36696 : CALL section_add_keyword(section, keyword)
1565 36696 : CALL keyword_release(keyword)
1566 :
1567 : CALL keyword_create(keyword, __LOCATION__, name="RWOH", &
1568 : description="Specify the rc parameter in the coordination function:"// &
1569 : " number of hydrogens per water molecule.", &
1570 : usage="RWOH {real}", default_r_val=cp_unit_to_cp2k(value=2.4_dp, &
1571 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1572 36696 : CALL section_add_keyword(section, keyword)
1573 36696 : CALL keyword_release(keyword)
1574 :
1575 : CALL keyword_create(keyword, __LOCATION__, name="pAOH", &
1576 : description="Sets the value of the numerator of the exponential factor"// &
1577 : " in the coordination function: number of hydrogens per carboxyl group.", &
1578 : usage="pAOH {integer}", default_i_val=6, &
1579 36696 : n_var=1)
1580 36696 : CALL section_add_keyword(section, keyword)
1581 36696 : CALL keyword_release(keyword)
1582 :
1583 : CALL keyword_create(keyword, __LOCATION__, name="qAOH", &
1584 : description="Sets the value of the denominator of the exponential factor"// &
1585 : " in the coordination function: number of hydrogens per carboxyl group.", &
1586 : usage="qAOH {integer}", default_i_val=14, &
1587 36696 : n_var=1)
1588 36696 : CALL section_add_keyword(section, keyword)
1589 36696 : CALL keyword_release(keyword)
1590 :
1591 : CALL keyword_create(keyword, __LOCATION__, name="RAOH", &
1592 : description="Specify the rc parameter in the coordination function:"// &
1593 : " number of hydrogens per carboxyl group.", &
1594 : usage="RAOH {real}", default_r_val=cp_unit_to_cp2k(value=2.4_dp, &
1595 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1596 36696 : CALL section_add_keyword(section, keyword)
1597 36696 : CALL keyword_release(keyword)
1598 :
1599 : CALL keyword_create(keyword, __LOCATION__, name="pCUT", &
1600 : description="Sets the value of the numerator of the exponential factor"// &
1601 : " in the cutoff function.", &
1602 : usage="pCUT {integer}", default_i_val=6, &
1603 36696 : n_var=1)
1604 36696 : CALL section_add_keyword(section, keyword)
1605 36696 : CALL keyword_release(keyword)
1606 :
1607 : CALL keyword_create(keyword, __LOCATION__, name="qCUT", &
1608 : description="Sets the value of the denominator of the exponential factor"// &
1609 : " in the cutoff function.", &
1610 : usage="qCUT {integer}", default_i_val=12, &
1611 36696 : n_var=1)
1612 36696 : CALL section_add_keyword(section, keyword)
1613 36696 : CALL keyword_release(keyword)
1614 :
1615 : CALL keyword_create(keyword, __LOCATION__, name="NC", &
1616 : description="Specifies the NC parameter in the cutoff function.", &
1617 : usage="NC {real}", default_r_val=0.56_dp, &
1618 36696 : n_var=1)
1619 36696 : CALL section_add_keyword(section, keyword)
1620 36696 : CALL keyword_release(keyword)
1621 :
1622 : CALL keyword_create(keyword, __LOCATION__, name="LAMBDA", &
1623 : variants=(/"LAMBDA"/), &
1624 : description="Specifies the LAMBDA parameter carboxylic acid function.", &
1625 : usage="LAMBDA {real}", default_r_val=20.0_dp, &
1626 73392 : n_var=1)
1627 36696 : CALL section_add_keyword(section, keyword)
1628 36696 : CALL keyword_release(keyword)
1629 :
1630 36696 : END SUBROUTINE create_colvar_acid_hyd_dist_section
1631 :
1632 : ! **************************************************************************************************
1633 : !> \brief collective variables specifying the solvation of carboxylic acid;
1634 : !> number of oxygens in the 1st shell of the hydronium; Equation (3) in
1635 : !> Supplementary Information of J. Am. Chem. Soc.,128, 2006, 11318
1636 : !> \param section the section to be created
1637 : !> \author Dorothea Golze
1638 : ! **************************************************************************************************
1639 36696 : SUBROUTINE create_colvar_acid_hyd_shell_section(section)
1640 : TYPE(section_type), POINTER :: section
1641 :
1642 : TYPE(keyword_type), POINTER :: keyword
1643 :
1644 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1645 : CALL section_create(section, __LOCATION__, name="ACID_HYDRONIUM_SHELL", &
1646 : description="Section to define the dissociation of a carboxylic acid in"// &
1647 : " water. Number of oxygens in the 1st shell of the hydronium. Equation (3)"// &
1648 : " in Supplementary Info of J. Am. Chem. Soc.,128, 2006, 11318. Similar to"// &
1649 : " the HYDRONIUM colvar, but with modification for the acid.", &
1650 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1651 :
1652 36696 : NULLIFY (keyword)
1653 :
1654 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS_WATER", &
1655 : description="Specifies indexes of atoms building the coordination variable."// &
1656 : " Oxygens of the water molecules. ", &
1657 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1658 36696 : n_var=-1, type_of_var=integer_t)
1659 36696 : CALL section_add_keyword(section, keyword)
1660 36696 : CALL keyword_release(keyword)
1661 :
1662 : CALL keyword_create(keyword, __LOCATION__, name="OXYGENS_ACID", &
1663 : description="Specifies indexes of atoms building the coordination variable."// &
1664 : " Oxygens of the carboxyl groups.", &
1665 : usage="OXYGENS {integer} {integer} ..", repeats=.TRUE., &
1666 36696 : n_var=-1, type_of_var=integer_t)
1667 36696 : CALL section_add_keyword(section, keyword)
1668 36696 : CALL keyword_release(keyword)
1669 :
1670 : CALL keyword_create(keyword, __LOCATION__, name="HYDROGENS", &
1671 : description="Specifies indexes of atoms building the coordination variable."// &
1672 : " Hydrogens of the water molecules and of the carboxyl groups.", &
1673 : usage="HYDROGENS {integer} {integer} ..", repeats=.TRUE., &
1674 36696 : n_var=-1, type_of_var=integer_t)
1675 36696 : CALL section_add_keyword(section, keyword)
1676 36696 : CALL keyword_release(keyword)
1677 :
1678 : CALL keyword_create(keyword, __LOCATION__, name="pWOH", &
1679 : description="Sets the value of the numerator of the exponential factor"// &
1680 : " in the coordination function: number of hydrogens per water molecule.", &
1681 : usage="pWOH {integer}", default_i_val=8, &
1682 36696 : n_var=1)
1683 36696 : CALL section_add_keyword(section, keyword)
1684 36696 : CALL keyword_release(keyword)
1685 :
1686 : CALL keyword_create(keyword, __LOCATION__, name="qWOH", &
1687 : description="Sets the value of the denominator of the exponential factor"// &
1688 : " in the coordination function: number of hydrogens per water molecule.", &
1689 : usage="qWOH {integer}", default_i_val=16, &
1690 36696 : n_var=1)
1691 36696 : CALL section_add_keyword(section, keyword)
1692 36696 : CALL keyword_release(keyword)
1693 :
1694 : CALL keyword_create(keyword, __LOCATION__, name="RWOH", &
1695 : description="Specifies the rc parameter in the coordination function:"// &
1696 : " number of hydrogens per water molecule.", &
1697 : usage="RWOH {real}", default_r_val=cp_unit_to_cp2k(value=2.4_dp, &
1698 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1699 36696 : CALL section_add_keyword(section, keyword)
1700 36696 : CALL keyword_release(keyword)
1701 :
1702 : CALL keyword_create(keyword, __LOCATION__, name="pAOH", &
1703 : description="Sets the value of the numerator of the exponential factor"// &
1704 : " in the coordination function: number of hydrogens per carboxyl group.", &
1705 : usage="pAOH {integer}", default_i_val=6, &
1706 36696 : n_var=1)
1707 36696 : CALL section_add_keyword(section, keyword)
1708 36696 : CALL keyword_release(keyword)
1709 :
1710 : CALL keyword_create(keyword, __LOCATION__, name="qAOH", &
1711 : description="Sets the value of the denominator of the exponential factor"// &
1712 : " in the coordination function: number of hydrogens per carboxyl group.", &
1713 : usage="qAOH {integer}", default_i_val=14, &
1714 36696 : n_var=1)
1715 36696 : CALL section_add_keyword(section, keyword)
1716 36696 : CALL keyword_release(keyword)
1717 :
1718 : CALL keyword_create(keyword, __LOCATION__, name="RAOH", &
1719 : description="Specifies the rc parameter in the coordination function:"// &
1720 : " number of hydrogens per carboxyl group.", &
1721 : usage="RAOH {real}", default_r_val=cp_unit_to_cp2k(value=2.4_dp, &
1722 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1723 36696 : CALL section_add_keyword(section, keyword)
1724 36696 : CALL keyword_release(keyword)
1725 :
1726 : CALL keyword_create(keyword, __LOCATION__, name="pOO", &
1727 : description="Sets the value of the numerator of the exponential factor"// &
1728 : " in the coordination function: number of oxygens per water oxygen.", &
1729 : usage="pOO {integer}", default_i_val=6, &
1730 36696 : n_var=1)
1731 36696 : CALL section_add_keyword(section, keyword)
1732 36696 : CALL keyword_release(keyword)
1733 :
1734 : CALL keyword_create(keyword, __LOCATION__, name="qOO", &
1735 : description="Sets the value of the denominator of the exponential factor"// &
1736 : " in the coordination function: number of oxygens per water oxygen.", &
1737 : usage="qOO {integer}", default_i_val=12, &
1738 36696 : n_var=1)
1739 36696 : CALL section_add_keyword(section, keyword)
1740 36696 : CALL keyword_release(keyword)
1741 :
1742 : CALL keyword_create(keyword, __LOCATION__, name="ROO", &
1743 : description="Specifies the rc parameter in the coordination function:"// &
1744 : " number of oxygens per water oxygen.", &
1745 : usage="ROO {real}", default_r_val=cp_unit_to_cp2k(value=5.5_dp, &
1746 36696 : unit_str="bohr"), unit_str="bohr", n_var=1)
1747 36696 : CALL section_add_keyword(section, keyword)
1748 36696 : CALL keyword_release(keyword)
1749 :
1750 : CALL keyword_create(keyword, __LOCATION__, name="pM", &
1751 : description="Sets the value of the numerator of the exponential factor"// &
1752 : " in the M function.", &
1753 : usage="pM {integer}", default_i_val=8, &
1754 36696 : n_var=1)
1755 36696 : CALL section_add_keyword(section, keyword)
1756 36696 : CALL keyword_release(keyword)
1757 :
1758 : CALL keyword_create(keyword, __LOCATION__, name="qM", &
1759 : description="Sets the value of the denominator of the exponential factor"// &
1760 : " in the M function.", &
1761 : usage="qM {integer}", default_i_val=16, &
1762 36696 : n_var=1)
1763 36696 : CALL section_add_keyword(section, keyword)
1764 36696 : CALL keyword_release(keyword)
1765 :
1766 : CALL keyword_create(keyword, __LOCATION__, name="NH", &
1767 : description="Specifies the NH parameter in the M function.", &
1768 : usage="NH {real}", default_r_val=2.2_dp, &
1769 36696 : n_var=1)
1770 36696 : CALL section_add_keyword(section, keyword)
1771 36696 : CALL keyword_release(keyword)
1772 :
1773 : CALL keyword_create(keyword, __LOCATION__, name="pCUT", &
1774 : description="Sets the value of the numerator of the exponential factor"// &
1775 : " in the cutoff function.", &
1776 : usage="pCUT {integer}", default_i_val=6, &
1777 36696 : n_var=1)
1778 36696 : CALL section_add_keyword(section, keyword)
1779 36696 : CALL keyword_release(keyword)
1780 :
1781 : CALL keyword_create(keyword, __LOCATION__, name="qCUT", &
1782 : description="Sets the value of the denominator of the exponential factor"// &
1783 : " in the cutoff function.", &
1784 : usage="qCUT {integer}", default_i_val=12, &
1785 36696 : n_var=1)
1786 36696 : CALL section_add_keyword(section, keyword)
1787 36696 : CALL keyword_release(keyword)
1788 :
1789 : CALL keyword_create(keyword, __LOCATION__, name="NC", &
1790 : description="Specifies the NC parameter in the cutoff function.", &
1791 : usage="NC {real}", default_r_val=0.9_dp, &
1792 36696 : n_var=1)
1793 36696 : CALL section_add_keyword(section, keyword)
1794 36696 : CALL keyword_release(keyword)
1795 :
1796 : CALL keyword_create(keyword, __LOCATION__, name="LAMBDA", &
1797 : variants=(/"LAMBDA"/), &
1798 : description="Specifies the LAMBDA parameter carboxylic acid function.", &
1799 : usage="LAMBDA {real}", default_r_val=10.0_dp, &
1800 73392 : n_var=1)
1801 36696 : CALL section_add_keyword(section, keyword)
1802 36696 : CALL keyword_release(keyword)
1803 :
1804 36696 : END SUBROUTINE create_colvar_acid_hyd_shell_section
1805 :
1806 : ! **************************************************************************************************
1807 : !> \brief ...
1808 : !> \param section ...
1809 : ! **************************************************************************************************
1810 36696 : SUBROUTINE create_colvar_rmsd_section(section)
1811 : TYPE(section_type), POINTER :: section
1812 :
1813 : TYPE(keyword_type), POINTER :: keyword
1814 : TYPE(section_type), POINTER :: subsection, subsubsection
1815 :
1816 36696 : CPASSERT(.NOT. ASSOCIATED(section))
1817 : CALL section_create(section, __LOCATION__, name="rmsd", &
1818 : description="Section to define a CV as function of RMSD computed with respect to"// &
1819 : " given reference configurations. For 2 configurations the colvar is equal to:"// &
1820 : " ss = (RMSDA-RMSDB)/(RMSDA+RMSDB), while if only 1 configuration is given, then the"// &
1821 : " colvar is just the RMSD from that frame.", &
1822 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
1823 :
1824 36696 : NULLIFY (keyword, subsection, subsubsection)
1825 : CALL keyword_create(keyword, __LOCATION__, name="SUBSET_TYPE", &
1826 : description="Define the subsytem used to compute the RMSD. With ALL the displacements"// &
1827 : " are mass-weighted, with LIST all weights are set to 1,"// &
1828 : " with WEIGHT_LIST a list of weights is expected from input.", &
1829 : usage="SUBSET_TYPE ALL", &
1830 : enum_c_vals=s2a("ALL", "LIST", "WEIGHT_LIST"), &
1831 : enum_i_vals=(/rmsd_all, rmsd_list, rmsd_weightlist/), &
1832 36696 : default_i_val=rmsd_all)
1833 36696 : CALL section_add_keyword(section, keyword)
1834 36696 : CALL keyword_release(keyword)
1835 :
1836 : CALL keyword_create(keyword, __LOCATION__, name="ALIGN_FRAMES", &
1837 : description="Whether the reference frames should be aligned to minimize the RMSD", &
1838 : usage="ALIGN_FRAME", &
1839 36696 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1840 36696 : CALL section_add_keyword(section, keyword)
1841 36696 : CALL keyword_release(keyword)
1842 :
1843 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
1844 : description="Specify indexes of atoms building the subset. ", &
1845 : usage="ATOMS {integer} {integer} ..", repeats=.TRUE., &
1846 36696 : n_var=-1, type_of_var=integer_t)
1847 36696 : CALL section_add_keyword(section, keyword)
1848 36696 : CALL keyword_release(keyword)
1849 :
1850 : CALL keyword_create(keyword, __LOCATION__, name="WEIGHTS", &
1851 : description="Specify weights of atoms building the subset. It is used only with WEIGHT_LIST ", &
1852 : usage="weightS {real} {real} ..", repeats=.TRUE., &
1853 36696 : n_var=-1, type_of_var=real_t)
1854 36696 : CALL section_add_keyword(section, keyword)
1855 36696 : CALL keyword_release(keyword)
1856 :
1857 : CALL section_create(subsection, __LOCATION__, name="FRAME", &
1858 : description="Specify coordinates of the frame (number of frames can be either 1 or 2)", &
1859 36696 : repeats=.TRUE.)
1860 :
1861 : CALL keyword_create(keyword, __LOCATION__, name="COORD_FILE_NAME", &
1862 : description="Name of the xyz file with coordinates (alternative to &COORD section)", &
1863 : usage="COORD_FILE_NAME <CHAR>", &
1864 36696 : default_lc_val="")
1865 36696 : CALL section_add_keyword(subsection, keyword)
1866 36696 : CALL keyword_release(keyword)
1867 :
1868 36696 : CALL create_coord_section_cv(subsubsection, "RMSD")
1869 36696 : CALL section_add_subsection(subsection, subsubsection)
1870 36696 : CALL section_release(subsubsection)
1871 :
1872 36696 : CALL section_add_subsection(section, subsection)
1873 36696 : CALL section_release(subsection)
1874 :
1875 36696 : END SUBROUTINE create_colvar_rmsd_section
1876 :
1877 : ! **************************************************************************************************
1878 : !> \brief collective variables specifying the space orthogonal to the reaction path
1879 : !> in the space spanned by the involved collective coordinates
1880 : !> \param section the section to be created
1881 : !> \author fschiff
1882 : ! **************************************************************************************************
1883 9174 : SUBROUTINE create_colvar_rpath_section(section)
1884 : TYPE(section_type), POINTER :: section
1885 :
1886 9174 : CPASSERT(.NOT. ASSOCIATED(section))
1887 : CALL section_create(section, __LOCATION__, name="REACTION_PATH", &
1888 : description="Section defining a one dimensional reaction path in an Q-dimensional space of colvars. "// &
1889 : "Constraining this colvar, allows to sample the space orthogonal to the reaction path, "// &
1890 : "both in the Q-dimensional colvar and 3N-Q remaining coordinates. "// &
1891 : "For the details of the function see cited literature.", &
1892 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
1893 18348 : citations=(/Branduardi2007/))
1894 :
1895 9174 : CALL keywords_colvar_path(section)
1896 9174 : END SUBROUTINE create_colvar_rpath_section
1897 :
1898 : ! **************************************************************************************************
1899 : !> \brief Distance from reaction path
1900 : !> \param section the section to be created
1901 : !> \author 01.2010
1902 : ! **************************************************************************************************
1903 9174 : SUBROUTINE create_colvar_dpath_section(section)
1904 : TYPE(section_type), POINTER :: section
1905 :
1906 9174 : CPASSERT(.NOT. ASSOCIATED(section))
1907 : CALL section_create(section, __LOCATION__, name="DISTANCE_FROM_PATH", &
1908 : description="Section defining the distance from a one dimensional reaction "// &
1909 : "path in an Q-dimensional space of colvars. "// &
1910 : "Constraining this colvar, allows to sample the space equidistant to the reaction path, "// &
1911 : "both in the Q-dimensional colvar and 3N-Q remaining coordinates. "// &
1912 : "For the details of the function see cited literature.", &
1913 : n_keywords=1, n_subsections=0, repeats=.FALSE., &
1914 18348 : citations=(/Branduardi2007/))
1915 :
1916 9174 : CALL keywords_colvar_path(section)
1917 9174 : END SUBROUTINE create_colvar_dpath_section
1918 :
1919 : ! **************************************************************************************************
1920 : !> \brief Section describinf keywords for both reaction path and distance from reaction path
1921 : !> \param section the section to be created
1922 : !> \author 01.2010
1923 : ! **************************************************************************************************
1924 18348 : SUBROUTINE keywords_colvar_path(section)
1925 :
1926 : TYPE(section_type), POINTER :: section
1927 :
1928 : TYPE(keyword_type), POINTER :: keyword
1929 : TYPE(section_type), POINTER :: print_key, subsection, subsubsection
1930 :
1931 18348 : NULLIFY (keyword, subsection, subsubsection, print_key)
1932 18348 : CALL create_colvar_section(subsection, skip_recursive_colvar=.TRUE.)
1933 18348 : CALL section_add_subsection(section, subsection)
1934 18348 : CALL section_release(subsection)
1935 :
1936 : CALL keyword_create(keyword, __LOCATION__, name="DISTANCES_RMSD", &
1937 : description=" ", &
1938 : usage="DISTANCES_RMSD T", &
1939 18348 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1940 18348 : CALL section_add_keyword(section, keyword)
1941 18348 : CALL keyword_release(keyword)
1942 :
1943 : CALL keyword_create(keyword, __LOCATION__, name="RMSD", &
1944 : description=" ", &
1945 : usage="RMSD T", &
1946 18348 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1947 18348 : CALL section_add_keyword(section, keyword)
1948 18348 : CALL keyword_release(keyword)
1949 :
1950 : CALL keyword_create(keyword, __LOCATION__, name="SUBSET_TYPE", &
1951 : description="Define the subsytem used to compute the RMSD", &
1952 : usage="SUBSET_TYPE ALL", &
1953 : enum_c_vals=s2a("ALL", "LIST"), &
1954 : enum_i_vals=(/rmsd_all, rmsd_list/), &
1955 18348 : default_i_val=rmsd_all)
1956 18348 : CALL section_add_keyword(section, keyword)
1957 18348 : CALL keyword_release(keyword)
1958 :
1959 : CALL keyword_create(keyword, __LOCATION__, name="ALIGN_FRAMES", &
1960 : description="Whether the reference frames should be aligned to minimize the RMSD", &
1961 : usage="ALIGN_FRAME", &
1962 18348 : default_l_val=.FALSE., lone_keyword_l_val=.TRUE.)
1963 18348 : CALL section_add_keyword(section, keyword)
1964 18348 : CALL keyword_release(keyword)
1965 :
1966 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
1967 : description="Specify indexes of atoms building the subset. ", &
1968 : usage="ATOMS {integer} {integer} ..", repeats=.TRUE., &
1969 18348 : n_var=-1, type_of_var=integer_t)
1970 18348 : CALL section_add_keyword(section, keyword)
1971 18348 : CALL keyword_release(keyword)
1972 :
1973 : CALL section_create(subsection, __LOCATION__, name="FRAME", &
1974 : description="Specify coordinates of the frame", &
1975 18348 : repeats=.TRUE.)
1976 :
1977 : CALL keyword_create(keyword, __LOCATION__, name="COORD_FILE_NAME", &
1978 : description="Name of the xyz file with coordinates (alternative to &COORD section)", &
1979 : usage="COORD_FILE_NAME <CHAR>", &
1980 18348 : default_lc_val="")
1981 18348 : CALL section_add_keyword(subsection, keyword)
1982 18348 : CALL keyword_release(keyword)
1983 :
1984 18348 : CALL create_coord_section_cv(subsubsection, "RMSD")
1985 18348 : CALL section_add_subsection(subsection, subsubsection)
1986 18348 : CALL section_release(subsubsection)
1987 :
1988 18348 : CALL section_add_subsection(section, subsection)
1989 18348 : CALL section_release(subsection)
1990 :
1991 : CALL keyword_create(keyword, __LOCATION__, name="FUNCTION", &
1992 : description="Specifies the ith element of the vector valued function that defines the reaction path. "// &
1993 : "This keyword needs to repeat exactly Q times, and the order must match the order of the colvars. "// &
1994 : "The VARIABLE (e.g. T) which parametrises the curve can be used as the target of a constraint.", &
1995 : usage="FUNCTION (sin(T+2)+2*T)", type_of_var=lchar_t, &
1996 18348 : n_var=1, default_lc_val="0", repeats=.TRUE.)
1997 18348 : CALL section_add_keyword(section, keyword)
1998 18348 : CALL keyword_release(keyword)
1999 :
2000 : CALL keyword_create(keyword, __LOCATION__, name="VARIABLE", &
2001 : description="Specifies the name of the variable that parametrises the FUNCTION "// &
2002 : "defining the reaction path.", &
2003 : usage="VARIABLE T", type_of_var=char_t, &
2004 18348 : n_var=1, repeats=.FALSE.)
2005 18348 : CALL section_add_keyword(section, keyword)
2006 18348 : CALL keyword_release(keyword)
2007 :
2008 : CALL keyword_create( &
2009 : keyword, __LOCATION__, name="LAMBDA", &
2010 : description="Specifies the exponent of the Gaussian used in the integral representation of the colvar. "// &
2011 : "The shape of the space orthogonal to the reaction path is defined by this choice. "// &
2012 : "In the limit of large values, it is given by the plane orthogonal to the path. "// &
2013 : "In practice, modest values are required for stable numerical integration.", &
2014 : usage="LAMBDA {real}", &
2015 18348 : type_of_var=real_t, default_r_val=5.0_dp)
2016 18348 : CALL section_add_keyword(section, keyword)
2017 18348 : CALL keyword_release(keyword)
2018 :
2019 : CALL keyword_create(keyword, __LOCATION__, name="STEP_SIZE", &
2020 : description="Step size in the numerical integration, "// &
2021 : "a few thousand points are common, and the proper number also depends on LAMBDA.", &
2022 : usage="STEP_SIZE {real}", &
2023 18348 : type_of_var=real_t, default_r_val=0.01_dp)
2024 18348 : CALL section_add_keyword(section, keyword)
2025 18348 : CALL keyword_release(keyword)
2026 :
2027 : CALL keyword_create(keyword, __LOCATION__, name="RANGE", &
2028 : description="The range of VARIABLE used for the parametrisation.", &
2029 : usage="RANGE <REAL> <REAL>", &
2030 18348 : n_var=2, type_of_var=real_t)
2031 18348 : CALL section_add_keyword(section, keyword)
2032 18348 : CALL keyword_release(keyword)
2033 :
2034 : CALL cp_print_key_section_create( &
2035 : print_key, __LOCATION__, name="MAP", &
2036 : description="Activating this print key will print once a file with the values of the FUNCTION on a grid "// &
2037 : "of COLVAR values in a specified range. "// &
2038 : "GRID_SPACING and RANGE for every COLVAR has to be specified again in the same order as they are in the input.", &
2039 18348 : print_level=high_print_level, filename="PATH")
2040 :
2041 : CALL keyword_create(keyword, __LOCATION__, name="RANGE", &
2042 : description="The range of of the grid of the COLVAR.", &
2043 : usage="RANGE <REAL> <REAL>", &
2044 18348 : n_var=2, type_of_var=real_t, repeats=.TRUE.)
2045 18348 : CALL section_add_keyword(print_key, keyword)
2046 18348 : CALL keyword_release(keyword)
2047 :
2048 : CALL keyword_create(keyword, __LOCATION__, name="GRID_SPACING", &
2049 : description="Distance between two gridpoints for the grid on the COLVAR", &
2050 : usage="STEP_SIZE {real}", repeats=.TRUE., &
2051 18348 : type_of_var=real_t, default_r_val=0.01_dp)
2052 18348 : CALL section_add_keyword(print_key, keyword)
2053 18348 : CALL keyword_release(keyword)
2054 :
2055 18348 : CALL section_add_subsection(section, print_key)
2056 18348 : CALL section_release(print_key)
2057 :
2058 18348 : END SUBROUTINE keywords_colvar_path
2059 :
2060 : ! **************************************************************************************************
2061 : !> \brief Colvar allowing a combination of COLVARS
2062 : !> \param section the section to be created
2063 : !> \author Teodoro Laino [tlaino] - 12.2008
2064 : ! **************************************************************************************************
2065 9174 : SUBROUTINE create_colvar_comb_section(section)
2066 : TYPE(section_type), POINTER :: section
2067 :
2068 : TYPE(keyword_type), POINTER :: keyword
2069 : TYPE(section_type), POINTER :: subsection
2070 :
2071 9174 : CPASSERT(.NOT. ASSOCIATED(section))
2072 : CALL section_create(section, __LOCATION__, name="COMBINE_COLVAR", &
2073 : description="Allows the possibility to combine several COLVARs into one COLVAR "// &
2074 : "with a generic function.", &
2075 9174 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
2076 :
2077 9174 : NULLIFY (keyword, subsection)
2078 9174 : CALL create_colvar_section(subsection, skip_recursive_colvar=.TRUE.)
2079 9174 : CALL section_add_subsection(section, subsection)
2080 9174 : CALL section_release(subsection)
2081 :
2082 : CALL keyword_create(keyword, __LOCATION__, name="FUNCTION", &
2083 : description="Specifies the function used to combine different COLVARs into one.", &
2084 : ! **************************************************************************************************
2085 : !> \brief ...
2086 : !> \param CV1^2 ...
2087 : !> \param CV2^2 ...
2088 : !> \param " ...
2089 : !> \param type_of_var=lchar_t ...
2090 : !> \param n_var=1 ...
2091 : !> \param error=error ...
2092 : ! **************************************************************************************************
2093 : usage="FUNCTION SQRT(CV1^2+CV2^2)", type_of_var=lchar_t, &
2094 9174 : n_var=1)
2095 9174 : CALL section_add_keyword(section, keyword)
2096 9174 : CALL keyword_release(keyword)
2097 :
2098 : CALL keyword_create(keyword, __LOCATION__, name="VARIABLES", &
2099 : description="Specifies the name of the variable that parametrises the FUNCTION "// &
2100 : "defining how COLVARS should be combined. The matching follows the same order of the "// &
2101 : "COLVARS definition in the input file.", &
2102 9174 : usage="VARIABLE CV1 CV2 CV3", type_of_var=char_t, n_var=-1, repeats=.FALSE.)
2103 9174 : CALL section_add_keyword(section, keyword)
2104 9174 : CALL keyword_release(keyword)
2105 :
2106 : CALL keyword_create(keyword, __LOCATION__, name="PARAMETERS", &
2107 : description="Defines the parameters of the functional form", &
2108 : usage="PARAMETERS a b D", type_of_var=char_t, &
2109 9174 : n_var=-1, repeats=.TRUE.)
2110 9174 : CALL section_add_keyword(section, keyword)
2111 9174 : CALL keyword_release(keyword)
2112 :
2113 : CALL keyword_create(keyword, __LOCATION__, name="VALUES", &
2114 : description="Defines the values of parameter of the functional form", &
2115 : usage="VALUES ", type_of_var=real_t, &
2116 9174 : n_var=-1, repeats=.TRUE., unit_str="internal_cp2k")
2117 9174 : CALL section_add_keyword(section, keyword)
2118 9174 : CALL keyword_release(keyword)
2119 :
2120 : CALL keyword_create(keyword, __LOCATION__, name="DX", &
2121 : description="Parameter used for computing the derivative of the combination "// &
2122 : "of COLVARs with the Ridders' method.", &
2123 9174 : usage="DX <REAL>", default_r_val=0.1_dp, unit_str="bohr")
2124 9174 : CALL section_add_keyword(section, keyword)
2125 9174 : CALL keyword_release(keyword)
2126 :
2127 : CALL keyword_create(keyword, __LOCATION__, name="ERROR_LIMIT", &
2128 : description="Checks that the error in computing the derivative is not larger than "// &
2129 : "the value set; in case error is larger a warning message is printed.", &
2130 9174 : usage="ERROR_LIMIT <REAL>", default_r_val=1.0E-12_dp)
2131 9174 : CALL section_add_keyword(section, keyword)
2132 9174 : CALL keyword_release(keyword)
2133 :
2134 9174 : END SUBROUTINE create_colvar_comb_section
2135 :
2136 : ! **************************************************************************************************
2137 : !> \brief Creates the coord section
2138 : !> \param section the section to create
2139 : !> \param name ...
2140 : !> \author teo
2141 : ! **************************************************************************************************
2142 55044 : SUBROUTINE create_coord_section_cv(section, name)
2143 : TYPE(section_type), POINTER :: section
2144 : CHARACTER(LEN=*), INTENT(IN) :: name
2145 :
2146 : TYPE(keyword_type), POINTER :: keyword
2147 :
2148 55044 : CPASSERT(.NOT. ASSOCIATED(section))
2149 : CALL section_create(section, __LOCATION__, name="coord", &
2150 : description="The positions for "//TRIM(name)//" used for restart", &
2151 55044 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
2152 55044 : NULLIFY (keyword)
2153 :
2154 : CALL keyword_create(keyword, __LOCATION__, name="_DEFAULT_KEYWORD_", &
2155 : description="Specify positions of the system", repeats=.TRUE., &
2156 55044 : usage="{Real} ...", type_of_var=real_t, n_var=-1)
2157 55044 : CALL section_add_keyword(section, keyword)
2158 55044 : CALL keyword_release(keyword)
2159 :
2160 55044 : END SUBROUTINE create_coord_section_cv
2161 :
2162 : ! **************************************************************************************************
2163 : !> \brief collective variables specifying h bonds
2164 : !> \param section the section to be created
2165 : !> \author alin m elena
2166 : ! **************************************************************************************************
2167 36696 : SUBROUTINE create_colvar_wc_section(section)
2168 : TYPE(section_type), POINTER :: section
2169 :
2170 : TYPE(keyword_type), POINTER :: keyword
2171 : TYPE(section_type), POINTER :: subsection
2172 :
2173 36696 : CPASSERT(.NOT. ASSOCIATED(section))
2174 : CALL section_create(section, __LOCATION__, name="wc", &
2175 : description="Section to define the hbond wannier centre as a collective variables.", &
2176 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
2177 36696 : NULLIFY (keyword, subsection)
2178 :
2179 : CALL keyword_create(keyword, __LOCATION__, name="RCUT", &
2180 : description="Parameter used for computing the cutoff radius for searching "// &
2181 : "the wannier centres around an atom", &
2182 : usage="RCUT <REAL>", default_r_val=0.529177208590000_dp, unit_str="angstrom", &
2183 36696 : type_of_var=real_t, repeats=.FALSE.)
2184 36696 : CALL section_add_keyword(section, keyword)
2185 36696 : CALL keyword_release(keyword)
2186 :
2187 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2188 : variants=(/"POINTS"/), &
2189 : description="Specifies the indexes of atoms/points defining the bond (Od, H, Oa).", &
2190 : usage="ATOMS {integer} {integer} {integer}", &
2191 73392 : n_var=3, type_of_var=integer_t, repeats=.TRUE.)
2192 36696 : CALL section_add_keyword(section, keyword)
2193 36696 : CALL keyword_release(keyword)
2194 :
2195 : ! Must be present in each colvar and handled properly
2196 36696 : CALL create_point_section(subsection)
2197 36696 : CALL section_add_subsection(section, subsection)
2198 36696 : CALL section_release(subsection)
2199 :
2200 36696 : END SUBROUTINE create_colvar_wc_section
2201 :
2202 : ! **************************************************************************************************
2203 : !> \brief collective variables specifying h bonds= wire
2204 : !> \param section the section to be created
2205 : !> \author alin m elena
2206 : ! **************************************************************************************************
2207 36696 : SUBROUTINE create_colvar_hbp_section(section)
2208 : TYPE(section_type), POINTER :: section
2209 :
2210 : TYPE(keyword_type), POINTER :: keyword
2211 : TYPE(section_type), POINTER :: subsection
2212 :
2213 36696 : CPASSERT(.NOT. ASSOCIATED(section))
2214 : CALL section_create(section, __LOCATION__, name="hbp", &
2215 : description="Section to define the hbond wannier centre as a collective variables.", &
2216 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
2217 36696 : NULLIFY (keyword, subsection)
2218 :
2219 : CALL keyword_create(keyword, __LOCATION__, name="RCUT", &
2220 : description="Parameter used for computing the cutoff radius for searching "// &
2221 : "the wannier centres around an atom", &
2222 : usage="RCUT <REAL>", default_r_val=0.529177208590000_dp, unit_str="angstrom", &
2223 36696 : type_of_var=real_t, repeats=.FALSE.)
2224 36696 : CALL section_add_keyword(section, keyword)
2225 36696 : CALL keyword_release(keyword)
2226 :
2227 : CALL keyword_create(keyword, __LOCATION__, name="SHIFT", &
2228 : description="Parameter used for shifting each term in the sum ", &
2229 : usage="SHIFT <REAL>", default_r_val=0.5_dp, &
2230 36696 : type_of_var=real_t, repeats=.FALSE.)
2231 36696 : CALL section_add_keyword(section, keyword)
2232 36696 : CALL keyword_release(keyword)
2233 :
2234 : CALL keyword_create(keyword, __LOCATION__, name="NPOINTS", &
2235 : description="The number of points in the path", &
2236 : usage="NPOINTS {integer}", default_i_val=-1, &
2237 36696 : n_var=1, type_of_var=integer_t, repeats=.FALSE.)
2238 36696 : CALL section_add_keyword(section, keyword)
2239 36696 : CALL keyword_release(keyword)
2240 :
2241 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2242 : variants=(/"POINTS"/), &
2243 : description="Specifies the indexes of atoms/points defining the bond (Od, H, Oa).", &
2244 : usage="ATOMS {integer} {integer} {integer}", &
2245 73392 : n_var=3, type_of_var=integer_t, repeats=.TRUE.)
2246 36696 : CALL section_add_keyword(section, keyword)
2247 36696 : CALL keyword_release(keyword)
2248 :
2249 : ! Must be present in each colvar and handled properly
2250 36696 : CALL create_point_section(subsection)
2251 36696 : CALL section_add_subsection(section, subsection)
2252 36696 : CALL section_release(subsection)
2253 :
2254 36696 : END SUBROUTINE create_colvar_hbp_section
2255 :
2256 : ! **************************************************************************************************
2257 : !> \brief collective variables specifying ring puckering
2258 : !> \brief D. Cremer and J.A. Pople, JACS 97 1354 (1975)
2259 : !> \param section the section to be created
2260 : !> \author Marcel Baer
2261 : ! **************************************************************************************************
2262 36696 : SUBROUTINE create_colvar_ring_puckering_section(section)
2263 : TYPE(section_type), POINTER :: section
2264 :
2265 : TYPE(keyword_type), POINTER :: keyword
2266 : TYPE(section_type), POINTER :: subsection
2267 :
2268 36696 : CPASSERT(.NOT. ASSOCIATED(section))
2269 : CALL section_create(section, __LOCATION__, name="RING_PUCKERING", &
2270 : description="Section to define general ring puckering collective variables.", &
2271 36696 : n_keywords=1, n_subsections=0, repeats=.FALSE.)
2272 :
2273 36696 : NULLIFY (keyword, subsection)
2274 :
2275 : CALL keyword_create(keyword, __LOCATION__, name="ATOMS", &
2276 : variants=(/"POINTS"/), &
2277 : description="Specifies the indexes of atoms/points defining the ring. "// &
2278 : "At least 4 Atoms are needed.", &
2279 : usage="ATOMS {integer} {integer} {integer} ..", &
2280 73392 : n_var=-1, type_of_var=integer_t)
2281 36696 : CALL section_add_keyword(section, keyword)
2282 36696 : CALL keyword_release(keyword)
2283 :
2284 : CALL keyword_create(keyword, __LOCATION__, name="COORDINATE", &
2285 : description="Indicate the coordinate to be used. Follow the Cremer-Pople definition for a N ring. "// &
2286 : "0 is the total puckering variable Q, "// &
2287 : "2..[N/2] are puckering coordinates. "// &
2288 : "-2..-[N/2-1] are puckering angles.", &
2289 : usage="COORDINATE {integer}", default_i_val=0, &
2290 36696 : n_var=1)
2291 36696 : CALL section_add_keyword(section, keyword)
2292 36696 : CALL keyword_release(keyword)
2293 :
2294 : ! Must be present in each colvar and handled properly
2295 36696 : CALL create_point_section(subsection)
2296 36696 : CALL section_add_subsection(section, subsection)
2297 36696 : CALL section_release(subsection)
2298 :
2299 36696 : END SUBROUTINE create_colvar_ring_puckering_section
2300 :
2301 : ! **************************************************************************************************
2302 :
2303 : END MODULE input_cp2k_colvar
|