(git:495eafe)
Loading...
Searching...
No Matches
input_cp2k_qmmm.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2026 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief creates the qmmm section of the input
10!> \note
11!> moved out of input_cp2k
12!> \par History
13!> 10.2005 split out of input_cp2k
14!> \author teo & fawzi
15! **************************************************************************************************
17 USE bibliography, ONLY: bernstein2009,&
19 golze2013,&
20 laino2005,&
22 USE cell_types, ONLY: use_perd_none
32 USE cp_units, ONLY: cp_unit_to_cp2k
33 USE input_constants, ONLY: &
59 USE input_val_types, ONLY: char_t,&
60 integer_t,&
61 lchar_t,&
62 logical_t,&
63 real_t
64 USE kinds, ONLY: dp
65 USE pw_spline_utils, ONLY: no_precond,&
71 USE string_utilities, ONLY: s2a
72#include "./base/base_uses.f90"
73
74 IMPLICIT NONE
75 PRIVATE
76
77 LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
78 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'input_cp2k_qmmm'
79
80 PUBLIC :: create_qmmm_section
81
82!***
83CONTAINS
84
85! **************************************************************************************************
86!> \brief Creates the QM/MM section
87!> \param section the section to create
88!> \author teo
89! **************************************************************************************************
90 SUBROUTINE create_qmmm_section(section)
91 TYPE(section_type), POINTER :: section
92
93 TYPE(keyword_type), POINTER :: keyword
94 TYPE(section_type), POINTER :: subsection
95
96 cpassert(.NOT. ASSOCIATED(section))
97 CALL section_create(section, __location__, name="qmmm", &
98 description="Controls QM/MM calculations, including mechanical, electrostatic, "// &
99 "Gaussian-expanded, and periodic embedding options.", &
100 n_keywords=6, n_subsections=3, repeats=.false., &
101 citations=[laino2005, laino2006])
102
103 NULLIFY (keyword, subsection)
104 CALL keyword_create(keyword, __location__, name="E_COUPL", &
105 variants=s2a("QMMM_COUPLING", "ECOUPL"), &
106 description="Selects the QM-MM coupling model used for the electrostatic interaction.", &
107 usage="E_COUPL GAUSS", &
108 enum_c_vals=s2a("NONE", "COULOMB", "GAUSS", "S-WAVE", "POINT_CHARGE"), &
110 enum_desc=s2a("Mechanical coupling (i.e. classical point charge based)", &
111 "Using analytical 1/r potential (Coulomb) - not available for GPW/GAPW", &
112 "Using fast Gaussian expansion of the electrostatic potential (Erf(r/rc)/r) "// &
113 "- not available for DFTB.", &
114 "Using fast Gaussian expansion of the s-wave electrostatic potential", &
115 "Using quantum mechanics derived point charges interacting with MM charges"), &
116 default_i_val=do_qmmm_none)
117 CALL section_add_keyword(section, keyword)
118 CALL keyword_release(keyword)
119
120 CALL keyword_create(keyword, __location__, name="MM_POTENTIAL_FILE_NAME", &
121 description="Name of the file containing the potential expansion in gaussians. See the "// &
122 "USE_GEEP_LIB keyword.", &
123 usage="MM_POTENTIAL_FILE_NAME {filename}", &
124 default_lc_val="MM_POTENTIAL")
125 CALL section_add_keyword(section, keyword)
126 CALL keyword_release(keyword)
127
128 CALL keyword_create(keyword, __location__, name="use_geep_lib", &
129 description=" This keyword enables the use of the internal GEEP library to generate the "// &
130 "gaussian expansion of the MM potential. Using this keyword there's no need to provide "// &
131 "the MM_POTENTIAL_FILENAME. It expects a number from 2 to 15 (the number of gaussian functions"// &
132 " to be used in the expansion.", &
133 usage="use_geep_lib INTEGER", &
134 default_i_val=0)
135 CALL section_add_keyword(section, keyword)
136 CALL keyword_release(keyword)
137
138 CALL keyword_create(keyword, __location__, name="nocompatibility", &
139 description="This keyword disables the compatibility of QM/MM "// &
140 "potential between CPMD and CP2K implementations. The compatibility"// &
141 " is achieved using an MM potential of the form: Erf[x/rc]/x + (1/rc -2/(pi^1/2*rc))*Exp[-(x/rc)^2]. "// &
142 "This keyword has effect only selecting GAUSS E_COUPLING type.", &
143 usage="nocompatibility LOGICAL", &
144 default_l_val=.false., lone_keyword_l_val=.true.)
145 CALL section_add_keyword(section, keyword)
146 CALL keyword_release(keyword)
147
148 CALL keyword_create(keyword, __location__, name="eps_mm_rspace", &
149 description="Set the threshold for the collocation of the GEEP gaussian functions. "// &
150 "this keyword affects only the GAUSS E_COUPLING.", &
151 usage="eps_mm_rspace real", &
152 default_r_val=1.0e-10_dp)
153 CALL section_add_keyword(section, keyword)
154 CALL keyword_release(keyword)
155
156 CALL keyword_create(keyword, __location__, name="SPHERICAL_CUTOFF", &
157 description="Set the spherical cutoff for the QMMM electrostatic interaction. "// &
158 "This acts like a charge multiplicative factor dependent on cutoff. For MM atoms "// &
159 "farther than the SPHERICAL_CUTOFF(1) their charge is zero. The switch is performed "// &
160 "with a smooth function: 0.5*(1-TANH((r-[SPH_CUT(1)-20*SPH_CUT(2)])/(SPH_CUT(2)))). "// &
161 "Two values are required: the first one is the distance cutoff. The second one controls "// &
162 "the stiffness of the smoothing.", &
163 usage="SPHERICAL_CUTOFF <REAL>", default_r_vals=[-1.0_dp, 0.0_dp], n_var=2, &
164 unit_str="angstrom")
165 CALL section_add_keyword(section, keyword)
166 CALL keyword_release(keyword)
167
168 CALL keyword_create(keyword, __location__, name="parallel_scheme", &
169 description="Chooses the parallel_scheme for the long range Potential ", &
170 usage="parallel_scheme (ATOM|GRID)", &
171 enum_c_vals=s2a("ATOM", "GRID"), &
172 enum_desc=s2a("parallelizes on atoms. grids replicated. "// &
173 "Replication of the grids can be quite expensive memory wise if running on a system "// &
174 "with limited memory per core. The grid option may be preferred in this case.", &
175 "parallelizes on grid slices. atoms replicated."), &
176 enum_i_vals=[do_par_atom, do_par_grid], &
177 default_i_val=do_par_atom)
178 CALL section_add_keyword(section, keyword)
179 CALL keyword_release(keyword)
180
181 ! Centering keywords
182 CALL keyword_create(keyword, __location__, name="CENTER", &
183 description="This keyword sets when the QM system is automatically "// &
184 "centered. Default is EVERY_STEP.", &
185 usage="center (EVERY_STEP|SETUP_ONLY|NEVER)", &
186 enum_c_vals=s2a("EVERY_STEP", "SETUP_ONLY", "NEVER"), &
187 enum_desc=s2a("Re-center every step", &
188 "Center at first step only", &
189 "Never center"), &
191 default_i_val=do_qmmm_center_every_step)
192 CALL section_add_keyword(section, keyword)
193 CALL keyword_release(keyword)
194
195 CALL keyword_create(keyword, __location__, name="CENTER_TYPE", &
196 description="This keyword specifies how to do the QM system centering.", &
197 usage="center_type (MAX_MINUS_MIN|PBC_AWARE_MAX_MINUS_MIN)", &
198 enum_c_vals=s2a("MAX_MINUS_MIN", "PBC_AWARE_MAX_MINUS_MIN"), &
199 enum_desc=s2a("Center of box defined by maximum coordinate minus minimum coordinate", &
200 "PBC-aware centering (useful for &QMMM&FORCE_MIXING)"), &
202 default_i_val=do_qmmm_center_max_minus_min)
203 CALL section_add_keyword(section, keyword)
204 CALL keyword_release(keyword)
205
206 CALL keyword_create(keyword, __location__, name="CENTER_GRID", &
207 description="This keyword specifies whether the QM system is centered in units of the grid spacing.", &
208 usage="CENTER_GRID LOGICAL", &
209 default_l_val=.false.)
210 CALL section_add_keyword(section, keyword)
211 CALL keyword_release(keyword)
212
213 CALL keyword_create(keyword, __location__, name="initial_translation_vector", &
214 description="This keyword specify the initial translation vector to be applied to the system.", &
215 usage="initial_translation_vector <REAL> <REAL> <REAL>", &
216 n_var=3, default_r_vals=[0.0_dp, 0.0_dp, 0.0_dp])
217 CALL section_add_keyword(section, keyword)
218 CALL keyword_release(keyword)
219
220 CALL keyword_create( &
221 keyword, __location__, name="DELTA_CHARGE", &
222 description="Additional net charge relative to that specified in DFT section. Used automatically by force mixing", &
223 usage="DELTA_CHARGE q", default_i_val=0, &
224 n_var=1, type_of_var=integer_t, repeats=.false.)
225 CALL section_add_keyword(section, keyword)
226 CALL keyword_release(keyword)
227
228 ! NB: remember to create these
229 CALL create_qmmm_force_mixing_section(subsection)
230 CALL section_add_subsection(section, subsection)
231 CALL section_release(subsection)
232
233 CALL create_qmmm_qm_kinds(subsection)
234 CALL section_add_subsection(section, subsection)
235 CALL section_release(subsection)
236
237 CALL create_qmmm_mm_kinds(subsection)
238 CALL section_add_subsection(section, subsection)
239 CALL section_release(subsection)
240
241 CALL create_cell_section(subsection, periodic=use_perd_none)
242 CALL section_add_subsection(section, subsection)
243 CALL section_release(subsection)
244
245 CALL create_qmmm_periodic_section(subsection)
246 CALL section_add_subsection(section, subsection)
247 CALL section_release(subsection)
248
249 CALL create_qmmm_link_section(subsection)
250 CALL section_add_subsection(section, subsection)
251 CALL section_release(subsection)
252
253 CALL create_qmmm_interp_section(subsection)
254 CALL section_add_subsection(section, subsection)
255 CALL section_release(subsection)
256
257 CALL create_qmmm_forcefield_section(subsection)
258 CALL section_add_subsection(section, subsection)
259 CALL section_release(subsection)
260
261 CALL create_qmmm_walls_section(subsection)
262 CALL section_add_subsection(section, subsection)
263 CALL section_release(subsection)
264
265 CALL create_qmmm_image_charge_section(subsection)
266 CALL section_add_subsection(section, subsection)
267 CALL section_release(subsection)
268
269 CALL create_print_qmmm_section(subsection)
270 CALL section_add_subsection(section, subsection)
271 CALL section_release(subsection)
272
273 END SUBROUTINE create_qmmm_section
274
275! **************************************************************************************************
276!> \brief Input section to create MM kinds sections
277!> \param section ...
278!> \author tlaino
279! **************************************************************************************************
280 SUBROUTINE create_qmmm_mm_kinds(section)
281 TYPE(section_type), POINTER :: section
282
283 TYPE(keyword_type), POINTER :: keyword
284
285 NULLIFY (keyword)
286 cpassert(.NOT. ASSOCIATED(section))
287 CALL section_create(section, __location__, name="MM_KIND", &
288 description="Information about the MM kind in the QM/MM scheme", &
289 n_keywords=2, n_subsections=0, repeats=.true.)
290
291 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
292 description="The MM kind", usage="O", n_var=1, type_of_var=char_t)
293 CALL section_add_keyword(section, keyword)
294 CALL keyword_release(keyword)
295
296 CALL keyword_create(keyword, __location__, name="RADIUS", &
297 description="Specifies the radius of the atomic kinds", &
298 usage="RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom", &
299 default_r_val=cp_unit_to_cp2k(radius_qmmm_default, "angstrom"))
300 CALL section_add_keyword(section, keyword)
301 CALL keyword_release(keyword)
302
303 CALL keyword_create(keyword, __location__, name="CORR_RADIUS", &
304 description="Specifies the correction radius of the atomic kinds"// &
305 " The correction radius is connected to the use of the compatibility keyword.", &
306 usage="CORR_RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom")
307 CALL section_add_keyword(section, keyword)
308 CALL keyword_release(keyword)
309
310 END SUBROUTINE create_qmmm_mm_kinds
311
312! **************************************************************************************************
313!> \brief Input section to create FORCE_MIXING sections
314!> \param section ...
315!> \author noam
316! **************************************************************************************************
317 SUBROUTINE create_qmmm_force_mixing_section(section)
318 TYPE(section_type), POINTER :: section
319
320 TYPE(keyword_type), POINTER :: keyword
321 TYPE(section_type), POINTER :: link_subsection, print_key, &
322 qm_kinds_subsection, subsection
323
324 NULLIFY (keyword)
325 cpassert(.NOT. ASSOCIATED(section))
326 CALL section_create(section, __location__, name="FORCE_MIXING", &
327 description="This section enables and defines parameters for force-mixing based QM/MM,"// &
328 " which actually does two conventional QM/MM calculations, on a small"// &
329 " and a large QM region, and combines the MM forces from one and QM"// &
330 " forces from the other to create a complete set of forces. Energy is"// &
331 " not conserved (although the QM/MM energy from the large QM region calculation is reported)"// &
332 " so a proper thermostat (i.e. massive, and able to handle dissipation, such as"// &
333 " Adaptive Langevin (AD_LANGEVIN)) must be used. For some propagation algorithms"// &
334 " (NVT and REFTRAJ MD ensembles) algorithm is adaptive,"// &
335 " including molecules hysteretically based on their instantaneous distance from the core region."// &
336 " Information on core/QM/buffer labels can be written in PDB file using"// &
337 " MOTION&PRINT&FORCE_MIXING_LABELS. Will fail if calculation requires a"// &
338 " meaningfull stress, or an energy that is consistent with the forces."// &
339 " For GEO_OPT this means"// &
340 " only MOTION&GEO_OPT&TYPE CG, MOTION&GEO_OPT&CG&LINE_SEARCH&TYPE 2PNT, and"// &
341 " MOTION&GEO_OPT&CG&LINE_SEARCH&2PNT&LINMIN_GRAD_ONLY T", &
342 n_keywords=5, n_subsections=3, repeats=.false., &
343 citations=[bernstein2009, bernstein2012])
344
345 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
346 description="Enables force-mixing", &
347 default_l_val=.false., lone_keyword_l_val=.true.)
348 CALL section_add_keyword(section, keyword)
349 CALL keyword_release(keyword)
350
351 CALL keyword_create(keyword, __location__, name="MOMENTUM_CONSERVATION_TYPE", &
352 description="How to apply force to get momentum conservation", &
353 usage="MOMENTUM_CONSERVATION_TYPE <type>", &
354 enum_c_vals=s2a("NONE", "EQUAL_F", "EQUAL_A"), &
356 enum_desc=s2a("No momentum conservation", &
357 "Equal force on each atom", &
358 "Equal acceleration on each atom"), &
359 default_i_val=do_fm_mom_conserv_equal_a)
360 CALL section_add_keyword(section, keyword)
361 CALL keyword_release(keyword)
362
363 CALL keyword_create(keyword, __location__, name="MOMENTUM_CONSERVATION_REGION", &
364 description="Region to apply correction force to for momentum conservation", &
365 usage="MOMENTUM_CONSERVATION_REGION <label>", &
366 enum_c_vals=s2a("CORE", "QM", "BUFFER"), &
368 enum_desc=s2a("Apply to QM core region", &
369 "Apply to full QM (dynamics) region", &
370 "Apply to QM+buffer regions"), &
371 default_i_val=do_fm_mom_conserv_qm)
372 CALL section_add_keyword(section, keyword)
373 CALL keyword_release(keyword)
374
375 CALL keyword_create(keyword, __location__, name="R_CORE", &
376 description="Specify the inner and outer radii of core QM region."// &
377 " All molecules with any atoms within this distance (hysteretically) of any atoms"// &
378 " specified as QM in enclosing QM/MM section will be core QM atoms in the force-mixing calculation.", &
379 usage="R_CORE <real> <real>", n_var=2, type_of_var=real_t, &
380 default_r_vals=[cp_unit_to_cp2k(0.0_dp, "angstrom"), &
381 cp_unit_to_cp2k(0.0_dp, "angstrom")], &
382 unit_str="angstrom")
383 CALL section_add_keyword(section, keyword)
384 CALL keyword_release(keyword)
385
386 CALL keyword_create(keyword, __location__, name="R_QM", &
387 description="Specify the inner and outer radii of QM dynamics region."// &
388 " All molecules with atoms within this distance (hysteretically) of any atoms in"// &
389 " core will follow QM dynamics in the force-mixing calculation.", &
390 usage="R_QM <real> <real>", n_var=2, type_of_var=real_t, &
391 default_r_vals=[cp_unit_to_cp2k(0.5_dp, "angstrom"), &
392 cp_unit_to_cp2k(1.0_dp, "angstrom")], &
393 unit_str="angstrom")
394 CALL section_add_keyword(section, keyword)
395 CALL keyword_release(keyword)
396
397 CALL keyword_create(keyword, __location__, name="QM_EXTENDED_SEED_IS_ONLY_CORE_LIST", &
398 description="Makes the extended QM zone be defined hysterestically"// &
399 " by distance from QM core list (i.e. atoms specified explicitly by"// &
400 " user) instead of from full QM core region (specified by user + hysteretic"// &
401 " selection + unbreakable bonds)", &
402 usage="QM_EXTENDED_SEED_IS_ONLY_CORE_LIST <logical>", n_var=1, type_of_var=logical_t, &
403 default_l_val=.false., repeats=.false.)
404 CALL section_add_keyword(section, keyword)
405 CALL keyword_release(keyword)
406
407 CALL keyword_create(keyword, __location__, name="R_BUF", &
408 description="Specify the inner and outer radii of buffer region."// &
409 " All atoms within this distance (hysteretically) of any QM atoms"// &
410 " will be buffer atoms in the force-mixing calculation.", &
411 usage="R_BUF <real> <real>", n_var=2, type_of_var=real_t, &
412 default_r_vals=[cp_unit_to_cp2k(0.5_dp, "angstrom"), &
413 cp_unit_to_cp2k(1.0_dp, "angstrom")], &
414 unit_str="angstrom")
415 CALL section_add_keyword(section, keyword)
416 CALL keyword_release(keyword)
417
418 CALL keyword_create(keyword, __location__, name="QM_KIND_ELEMENT_MAPPING", &
419 description="Mapping from elements to QM_KINDs for adaptively included atoms.", &
420 usage="QM_KIND_ELEMENT_MAPPING {El} {QM_KIND}", &
421 n_var=2, type_of_var=char_t, repeats=.true.)
422 CALL section_add_keyword(section, keyword)
423 CALL keyword_release(keyword)
424
425 CALL keyword_create(keyword, __location__, name="MAX_N_QM", &
426 description="Maximum number of QM atoms, for detection of runaway adaptive selection.", &
427 usage="MAX_N_QM int", default_i_val=300, &
428 n_var=1, type_of_var=integer_t, repeats=.false.)
429 CALL section_add_keyword(section, keyword)
430 CALL keyword_release(keyword)
431
432 CALL keyword_create(keyword, __location__, name="ADAPTIVE_EXCLUDE_MOLECULES", &
433 description="List of molecule names to exclude from adaptive regions (e.g. big things like proteins)", &
434 usage="ADAPTIVE_EXCLUDE_MOLECULES molec1 molec2 ...", &
435 n_var=-1, type_of_var=char_t, repeats=.false.)
436 CALL section_add_keyword(section, keyword)
437 CALL keyword_release(keyword)
438
439 CALL keyword_create(keyword, __location__, name="EXTENDED_DELTA_CHARGE", &
440 description="Additional net charge in extended region relative to core (core charge is"// &
441 " specified in DFT section, as usual for a convetional QM/MM calculation)", &
442 usage="EXTENDED_DELTA_CHARGE q", default_i_val=0, &
443 n_var=1, type_of_var=integer_t, repeats=.false.)
444 CALL section_add_keyword(section, keyword)
445 CALL keyword_release(keyword)
446
447 ! QM_NON_ADAPTIVE subsection
448 NULLIFY (subsection)
449 CALL section_create(subsection, __location__, name="QM_NON_ADAPTIVE", &
450 description="List of atoms always in QM region, non-adaptively", &
451 n_keywords=0, n_subsections=1, repeats=.true.)
452
453 NULLIFY (qm_kinds_subsection)
454 CALL create_qmmm_qm_kinds(qm_kinds_subsection)
455 CALL section_add_subsection(subsection, qm_kinds_subsection)
456 CALL section_release(qm_kinds_subsection)
457
458 CALL section_add_subsection(section, subsection)
459 CALL section_release(subsection)
460
461 ! BUFFER_NON_ADAPTIVE subsection
462 NULLIFY (subsection)
463 CALL section_create(subsection, __location__, name="BUFFER_NON_ADAPTIVE", &
464 description="List of atoms always in buffer region, non-adaptively, and any needed LINK sections", &
465 n_keywords=0, n_subsections=1, repeats=.true.)
466
467 NULLIFY (qm_kinds_subsection)
468 CALL create_qmmm_qm_kinds(qm_kinds_subsection)
469 CALL section_add_subsection(subsection, qm_kinds_subsection)
470 CALL section_release(qm_kinds_subsection)
471 NULLIFY (link_subsection)
472 CALL create_qmmm_link_section(link_subsection)
473 CALL section_add_subsection(subsection, link_subsection)
474 CALL section_release(link_subsection)
475
476 CALL section_add_subsection(section, subsection)
477 CALL section_release(subsection)
478
479 ![NB] also need a list?
480 ![NB] maybe not list+links , but some sort of link template
481 ![NB] also, breakable bonds?
482 ! BUFFER_LINKS subsection
483 NULLIFY (subsection)
484 CALL section_create( &
485 subsection, __location__, name="BUFFER_LINKS", &
486 description="Information about possible links for automatic covalent bond breaking for the buffer QM/MM calculation. "// &
487 "Ignored - need to implement buffer selection by atom and walking of connectivity data.", &
488 n_keywords=0, n_subsections=1, repeats=.true.)
489
490 NULLIFY (link_subsection)
491 CALL create_qmmm_link_section(link_subsection)
492 CALL section_add_subsection(subsection, link_subsection)
493 CALL section_release(link_subsection)
494
495 CALL section_add_subsection(section, subsection)
496 CALL section_release(subsection)
497
498 ! RESTART_INFO subsection
499 NULLIFY (subsection)
500 CALL section_create(subsection, __location__, name="RESTART_INFO", &
501 description="This section provides information about old force-mixing indices and labels, "// &
502 "for restarts.", &
503 n_keywords=2, n_subsections=0, repeats=.false.)
504
505 CALL keyword_create(keyword, __location__, name="INDICES", &
506 description="Indices of atoms in previous step QM regions.", &
507 usage="INDICES 1 2 ...", &
508 n_var=-1, type_of_var=integer_t, repeats=.true.)
509 CALL section_add_keyword(subsection, keyword)
510 CALL keyword_release(keyword)
511
512 CALL keyword_create(keyword, __location__, name="LABELS", &
513 description="Labels of atoms in previous step QM regions.", &
514 usage="LABELS 1 1 ...", &
515 n_var=-1, type_of_var=integer_t, repeats=.true.)
516 CALL section_add_keyword(subsection, keyword)
517 CALL keyword_release(keyword)
518
519 CALL section_add_subsection(section, subsection)
520 CALL section_release(subsection)
521
522 ! PRINT subsection, with keys for neighbor list
523 CALL section_create(subsection, __location__, name="print", &
524 description="Section of possible print options in FORCE_MIXING.", &
525 n_keywords=0, n_subsections=2, repeats=.false.)
526 NULLIFY (print_key)
527 CALL cp_print_key_section_create(print_key, __location__, "SUBCELL", &
528 description="Activates the printing of the subcells used for the "// &
529 "generation of neighbor lists.", unit_str="angstrom", &
530 print_level=high_print_level, filename="__STD_OUT__")
531 CALL section_add_subsection(subsection, print_key)
532 CALL section_release(print_key)
533
534 CALL cp_print_key_section_create(print_key, __location__, "NEIGHBOR_LISTS", &
535 description="Activates the printing of the neighbor lists used"// &
536 " for the hysteretic region calculations.", &
537 print_level=high_print_level, filename="", unit_str="angstrom")
538 CALL section_add_subsection(subsection, print_key)
539 CALL section_release(print_key)
540
541 CALL section_add_subsection(section, subsection)
542 CALL section_release(subsection)
543
544 END SUBROUTINE create_qmmm_force_mixing_section
545
546! **************************************************************************************************
547!> \brief Input section to create QM kinds sections
548!> \param section ...
549!> \author tlaino
550! **************************************************************************************************
551 SUBROUTINE create_qmmm_qm_kinds(section)
552 TYPE(section_type), POINTER :: section
553
554 TYPE(keyword_type), POINTER :: keyword
555
556 NULLIFY (keyword)
557 cpassert(.NOT. ASSOCIATED(section))
558 CALL section_create(section, __location__, name="QM_KIND", &
559 description="Information about the QM kind in the QM/MM scheme", &
560 n_keywords=3, n_subsections=0, repeats=.true.)
561
562 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
563 description="The QM kind", usage="O", n_var=1, type_of_var=char_t)
564 CALL section_add_keyword(section, keyword)
565 CALL keyword_release(keyword)
566
567 CALL keyword_create(keyword, __location__, name="MM_INDEX", &
568 description="The indexes of the MM atoms that have this kind. This keyword can be"// &
569 " repeated several times (useful if you have to specify many indexes).", &
570 usage="MM_INDEX 1 2", &
571 n_var=-1, type_of_var=integer_t, repeats=.true.)
572 CALL section_add_keyword(section, keyword)
573 CALL keyword_release(keyword)
574
575 END SUBROUTINE create_qmmm_qm_kinds
576
577! **************************************************************************************************
578!> \brief Input section to set QM/MM periodic boundary conditions
579!> \param section ...
580!> \author tlaino
581! **************************************************************************************************
582 SUBROUTINE create_qmmm_walls_section(section)
583 TYPE(section_type), POINTER :: section
584
585 TYPE(keyword_type), POINTER :: keyword
586
587 NULLIFY (keyword)
588 cpassert(.NOT. ASSOCIATED(section))
589 CALL section_create(section, __location__, name="WALLS", &
590 description="Enables Walls for the QM box. This can be used to avoid that QM"// &
591 " atoms move out of the QM box.", &
592 n_keywords=0, n_subsections=0, repeats=.false.)
593
594 CALL keyword_create(keyword, __location__, name="WALL_SKIN", &
595 description="Specify the value of the skin of the Wall in each dimension. "// &
596 "The wall's effect is felt when atoms fall within the skin of the Wall.", &
597 usage="WALL_SKIN <real> <real> <real>", n_var=3, type_of_var=real_t, &
598 default_r_vals=[cp_unit_to_cp2k(0.5_dp, "angstrom"), &
599 cp_unit_to_cp2k(0.5_dp, "angstrom"), &
600 cp_unit_to_cp2k(0.5_dp, "angstrom")], &
601 unit_str="angstrom")
602 CALL section_add_keyword(section, keyword)
603 CALL keyword_release(keyword)
604
605 CALL keyword_create(keyword, __location__, name="TYPE", &
606 description="Specifies the type of wall", &
607 usage="TYPE REFLECTIVE", &
608 enum_c_vals=s2a("NONE", "REFLECTIVE", "QUADRATIC"), &
610 enum_desc=s2a("No Wall around QM box", &
611 "Reflective Wall around QM box", &
612 "Quadratic Wall around QM box"), &
613 default_i_val=do_qmmm_wall_reflective)
614 CALL section_add_keyword(section, keyword)
615 CALL keyword_release(keyword)
616
617 CALL keyword_create(keyword, __location__, name="K", &
618 description="Specify the value of the the force constant for the quadratic wall", &
619 usage="K <real>", unit_str='internal_cp2k', &
620 type_of_var=real_t)
621 CALL section_add_keyword(section, keyword)
622 CALL keyword_release(keyword)
623
624 END SUBROUTINE create_qmmm_walls_section
625
626! ****************************************************************************
627!> \brief Input section for QM/MM image charge calculations
628!> \param section ...
629!> \author Dorothea Golze
630! **************************************************************************************************
631 SUBROUTINE create_qmmm_image_charge_section(section)
632 TYPE(section_type), POINTER :: section
633
634 TYPE(keyword_type), POINTER :: keyword
635 TYPE(section_type), POINTER :: subsection
636
637 NULLIFY (keyword, subsection)
638 cpassert(.NOT. ASSOCIATED(section))
639 CALL section_create(section, __location__, name="IMAGE_CHARGE", &
640 description="Inclusion of polarization effects within the image charge "// &
641 "approach for systems where QM molecules are physisorbed on e.g. metal "// &
642 "surfaces described by MM. This correction introduces only a very small overhead. "// &
643 "QM box size has to be equal to MM box size.", &
644 n_keywords=3, n_subsections=1, repeats=.false., &
645 citations=[golze2013])
646
647 CALL keyword_create(keyword, __location__, name="MM_ATOM_LIST", &
648 description="List of MM atoms carrying an induced Gaussian charge. "// &
649 "If this keyword is not given, all MM atoms will carry an image charge.", &
650 usage="MM_ATOM_LIST 1 2 3 or 1..3 ", n_var=-1, type_of_var=integer_t, &
651 repeats=.true.)
652 CALL section_add_keyword(section, keyword)
653 CALL keyword_release(keyword)
654
655 CALL keyword_create(keyword, __location__, name="WIDTH", &
656 description="Specifies the value of the width of the (induced) Gaussian "// &
657 "charge distribution carried by each MM atom.", &
658 usage="WIDTH <real> ", n_var=1, type_of_var=real_t, &
659 default_r_val=cp_unit_to_cp2k(value=3.0_dp, unit_str="angstrom^-2"), &
660 unit_str="angstrom^-2")
661 CALL section_add_keyword(section, keyword)
662 CALL keyword_release(keyword)
663
664 CALL keyword_create(keyword, __location__, name="EXT_POTENTIAL", &
665 description="External potential applied to the metal electrode ", &
666 usage="EXT_POTENTIAL <real> ", n_var=1, type_of_var=real_t, &
667 default_r_val=0.0_dp, &
668 unit_str="volt")
669 CALL section_add_keyword(section, keyword)
670 CALL keyword_release(keyword)
671
672 CALL keyword_create(keyword, __location__, name="DETERM_COEFF", &
673 description="Specifies how the coefficients are determined.", &
674 usage="DETERM_COEFF ITERATIVE", &
675 enum_c_vals=s2a("CALC_MATRIX", "ITERATIVE"), &
677 enum_desc=s2a("Calculates image matrix and solves linear set of equations", &
678 "Uses an iterative scheme to calculate the coefficients"), &
679 default_i_val=do_qmmm_image_calcmatrix)
680 CALL section_add_keyword(section, keyword)
681 CALL keyword_release(keyword)
682
683 CALL keyword_create(keyword, __location__, name="RESTART_IMAGE_MATRIX", &
684 description="Restart the image matrix. Useful when "// &
685 "calculating coefficients iteratively (the image matrix "// &
686 "is used as preconditioner in that case)", &
687 usage="RESTART_IMAGE_MATRIX", default_l_val=.false., &
688 lone_keyword_l_val=.true.)
689 CALL section_add_keyword(section, keyword)
690 CALL keyword_release(keyword)
691
692 CALL keyword_create(keyword, __location__, name="IMAGE_RESTART_FILE_NAME", &
693 description="File name where to read the image matrix used "// &
694 "as preconditioner in the iterative scheme", &
695 usage="IMAGE_RESTART_FILE_NAME <FILENAME>", &
696 type_of_var=lchar_t)
697 CALL section_add_keyword(section, keyword)
698 CALL keyword_release(keyword)
699
700 CALL keyword_create(keyword, __location__, name="IMAGE_MATRIX_METHOD", &
701 description="Method for calculating the image matrix.", &
702 usage="IMAGE_MATRIX_METHOD MME", &
703 enum_c_vals=s2a("GPW", "MME"), &
704 enum_i_vals=[do_eri_gpw, do_eri_mme], &
705 enum_desc=s2a("Uses Gaussian Plane Wave method [Golze2013]", &
706 "Uses MiniMax-Ewald method (ERI_MME subsection)"), &
707 default_i_val=do_eri_mme)
708 CALL section_add_keyword(section, keyword)
709 CALL keyword_release(keyword)
710
711 ! for qmmm image charges we can afford taking the most accurate minimax approximation
712 CALL create_eri_mme_section(subsection, default_n_minimax=53)
713 CALL section_add_subsection(section, subsection)
714 CALL section_release(subsection)
715
716 END SUBROUTINE create_qmmm_image_charge_section
717
718! **************************************************************************************************
719!> \brief Input section to set QM/MM periodic boundary conditions
720!> \param section ...
721!> \author tlaino
722! **************************************************************************************************
723 SUBROUTINE create_qmmm_periodic_section(section)
724 TYPE(section_type), POINTER :: section
725
726 TYPE(keyword_type), POINTER :: keyword
727 TYPE(section_type), POINTER :: subsection
728
729 NULLIFY (keyword, subsection)
730 cpassert(.NOT. ASSOCIATED(section))
731 CALL section_create(section, __location__, name="PERIODIC", &
732 description="Specify parameters for QM/MM periodic boundary conditions calculations", &
733 n_keywords=0, n_subsections=0, repeats=.false., &
734 citations=[laino2006])
735
736 CALL keyword_create( &
737 keyword, __location__, name="GMAX", &
738 description="Specifies the maximum value of G in the reciprocal space over which perform the Ewald sum.", &
739 usage="GMAX <real>", n_var=1, default_r_val=1.0_dp)
740 CALL section_add_keyword(section, keyword)
741 CALL keyword_release(keyword)
742
743 CALL keyword_create(keyword, __location__, name="REPLICA", &
744 description="Specifies the number of replica to take into consideration for the real part of the "// &
745 "calculation. Default is letting the qmmm module decide how many replica you really need.", &
746 usage="REPLICA <integer>", n_var=1, default_i_val=-1)
747 CALL section_add_keyword(section, keyword)
748 CALL keyword_release(keyword)
749
750 CALL keyword_create(keyword, __location__, name="NGRIDS", &
751 description="Specifies the number of grid points used for the Interpolation of the G-space term", &
752 usage="NGRIDS <integer> <integer> <integer> ", n_var=3, default_i_vals=[50, 50, 50])
753 CALL section_add_keyword(section, keyword)
754 CALL keyword_release(keyword)
755
756 CALL create_multipole_qmmm_section(subsection)
757 CALL section_add_subsection(section, subsection)
758 CALL section_release(subsection)
759
760 CALL create_gspace_interp_section(subsection)
761 CALL section_add_subsection(section, subsection)
762 CALL section_release(subsection)
763
764 CALL create_poisson_section(subsection)
765 CALL section_add_subsection(section, subsection)
766 CALL section_release(subsection)
767
768 CALL cp_print_key_section_create(subsection, __location__, "check_spline", &
769 description="Controls the checking of the G-space term Spline Interpolation.", &
770 print_level=medium_print_level, filename="GSpace-SplInterp")
771 CALL section_add_subsection(section, subsection)
772 CALL section_release(subsection)
773
774 END SUBROUTINE create_qmmm_periodic_section
775
776! **************************************************************************************************
777!> \brief Section to set-up parameters for decoupling using the Bloechl scheme
778!> \param section the section to create
779!> \par History
780!> Dorothea Golze [04.2014] copied from input_cp2k_poisson.F and
781!> enabled switch-on/off
782!> \author teo
783! **************************************************************************************************
784 SUBROUTINE create_multipole_qmmm_section(section)
785 TYPE(section_type), POINTER :: section
786
787 TYPE(keyword_type), POINTER :: keyword
788 TYPE(section_type), POINTER :: subsection
789
790 cpassert(.NOT. ASSOCIATED(section))
791
792 CALL section_create(section, __location__, name="MULTIPOLE", &
793 description="This section is used to set up the decoupling of QM periodic images with "// &
794 "the use of density derived atomic point charges. Switched on by default even if not "// &
795 "explicitly given. Can be switched off if e.g. QM and MM box are of the same size.", &
796 n_keywords=1, n_subsections=0, repeats=.false.)
797
798 NULLIFY (keyword, subsection)
799 CALL keyword_create(keyword, __location__, name="_SECTION_PARAMETERS_", &
800 description="Defines the usage of the multipole section", &
801 usage="ON", &
802 enum_c_vals=s2a("ON", "OFF"), &
804 enum_desc=s2a("switch on MULTIPOLE section", &
805 "switch off MULTIPOLE section"), &
806 default_i_val=do_multipole_section_on, lone_keyword_i_val=do_multipole_section_on)
807 CALL section_add_keyword(section, keyword)
808 CALL keyword_release(keyword)
809
810 CALL keyword_create(keyword, __location__, name="RCUT", &
811 description="Real space cutoff for the Ewald sum.", &
812 usage="RCUT {real}", n_var=1, type_of_var=real_t, &
813 unit_str="angstrom")
814 CALL section_add_keyword(section, keyword)
815 CALL keyword_release(keyword)
816
817 CALL keyword_create(keyword, __location__, name="EWALD_PRECISION", &
818 description="Precision achieved in the Ewald sum.", &
819 usage="EWALD_PRECISION {real}", n_var=1, type_of_var=real_t, &
820 unit_str="hartree", default_r_val=1.0e-6_dp)
821 CALL section_add_keyword(section, keyword)
822 CALL keyword_release(keyword)
823
824 CALL keyword_create(keyword, __location__, name="ANALYTICAL_GTERM", &
825 description="Evaluates the Gterm in the Ewald Scheme analytically instead of using Splines.", &
826 usage="ANALYTICAL_GTERM <LOGICAL>", &
827 default_l_val=.false., lone_keyword_l_val=.true.)
828 CALL section_add_keyword(section, keyword)
829 CALL keyword_release(keyword)
830
831 CALL keyword_create(keyword, __location__, name="NGRIDS", &
832 description="Specifies the number of grid points used for the Interpolation of the G-space term", &
833 usage="NGRIDS <integer> <integer> <integer> ", n_var=3, default_i_vals=[50, 50, 50])
834 CALL section_add_keyword(section, keyword)
835 CALL keyword_release(keyword)
836
837 CALL create_gspace_interp_section(subsection)
838 CALL section_add_subsection(section, subsection)
839 CALL section_release(subsection)
840
841 CALL cp_print_key_section_create(subsection, __location__, "check_spline", &
842 description="Controls the checking of the G-space term Spline Interpolation.", &
843 print_level=medium_print_level, filename="GSpace-SplInterp")
844 CALL section_add_subsection(section, subsection)
845 CALL section_release(subsection)
846
847 CALL cp_print_key_section_create(subsection, __location__, "program_run_info", &
848 description="Controls the printing of basic information during the run", &
849 print_level=low_print_level, add_last=add_last_numeric, filename="__STD_OUT__")
850 CALL section_add_subsection(section, subsection)
851 CALL section_release(subsection)
852
853 END SUBROUTINE create_multipole_qmmm_section
854
855! **************************************************************************************************
856!> \brief creates the qm/mm forcefield section to override to the FF specification
857!> given in the FIST input
858!> \param section ...
859!> \author tlaino
860! **************************************************************************************************
861 SUBROUTINE create_qmmm_forcefield_section(section)
862 TYPE(section_type), POINTER :: section
863
864 TYPE(keyword_type), POINTER :: keyword
865 TYPE(section_type), POINTER :: subsection
866
867 NULLIFY (subsection, keyword)
868 cpassert(.NOT. ASSOCIATED(section))
869 CALL section_create(section, __location__, name="FORCEFIELD", &
870 description="Specify information on the QM/MM forcefield", &
871 n_keywords=0, n_subsections=2, repeats=.true.)
872
873 CALL keyword_create(keyword, __location__, name="MULTIPLE_POTENTIAL", &
874 description="Enables the possibility to define NONBONDED and NONBONDED14 as a"// &
875 " sum of different kinds of potential. Useful for piecewise defined potentials.", &
876 usage="MULTIPLE_POTENTIAL T", default_l_val=.false., lone_keyword_l_val=.true.)
877 CALL section_add_keyword(section, keyword)
878 CALL keyword_release(keyword)
879
880 CALL create_qmmm_ff_nb_section(subsection)
881 CALL section_add_subsection(section, subsection)
882 CALL section_release(subsection)
883
884 CALL create_nonbonded14_section(subsection)
885 CALL section_add_subsection(section, subsection)
886 CALL section_release(subsection)
887
888 END SUBROUTINE create_qmmm_forcefield_section
889
890! **************************************************************************************************
891!> \brief creates the qm/mm forcefield section to override to the FF specification
892!> given in the FIST input - NONBONDED PART
893!> \param section ...
894!> \author tlaino
895! **************************************************************************************************
896 SUBROUTINE create_qmmm_ff_nb_section(section)
897 TYPE(section_type), POINTER :: section
898
899 TYPE(section_type), POINTER :: subsection
900
901 NULLIFY (subsection)
902 cpassert(.NOT. ASSOCIATED(section))
903 CALL section_create(section, __location__, name="NONBONDED", &
904 description="Specify information on the QM/MM non-bonded forcefield", &
905 n_keywords=0, n_subsections=2, repeats=.true.)
906
907 CALL create_lj_section(subsection)
908 CALL section_add_subsection(section, subsection)
909 CALL section_release(subsection)
910
911 CALL create_williams_section(subsection)
912 CALL section_add_subsection(section, subsection)
913 CALL section_release(subsection)
914
915 CALL create_goodwin_section(subsection)
916 CALL section_add_subsection(section, subsection)
917 CALL section_release(subsection)
918
919 CALL create_genpot_section(subsection)
920 CALL section_add_subsection(section, subsection)
921 CALL section_release(subsection)
922
923 END SUBROUTINE create_qmmm_ff_nb_section
924
925! **************************************************************************************************
926!> \brief creates the qm/mm link section
927!> \param section ...
928!> \author tlaino
929! **************************************************************************************************
930 SUBROUTINE create_qmmm_link_section(section)
931 TYPE(section_type), POINTER :: section
932
933 TYPE(keyword_type), POINTER :: keyword
934 TYPE(section_type), POINTER :: subsection
935
936 NULLIFY (keyword, subsection)
937 cpassert(.NOT. ASSOCIATED(section))
938 CALL section_create(section, __location__, name="LINK", &
939 description="Specify information on the QM/MM link treatment", &
940 n_keywords=7, n_subsections=2, repeats=.true.)
941
942 CALL keyword_create(keyword, __location__, name="QM_INDEX", &
943 variants=["QM"], &
944 description="Specifies the index of the QM atom involved in the QM/MM link", &
945 usage="QM_INDEX integer", n_var=1, type_of_var=integer_t)
946 CALL section_add_keyword(section, keyword)
947 CALL keyword_release(keyword)
948
949 CALL keyword_create(keyword, __location__, name="QM_KIND", &
950 description="Specifies the element of the QM capping atom involved in the QM/MM link", &
951 usage="QM_KIND char", n_var=1, type_of_var=char_t, &
952 default_c_val="H")
953 CALL section_add_keyword(section, keyword)
954 CALL keyword_release(keyword)
955
956 CALL keyword_create(keyword, __location__, name="MM_INDEX", &
957 variants=["MM"], &
958 description="Specifies the index of the MM atom involved in the QM/MM link, Default hydrogen.", &
959 usage="MM_INDEX integer", n_var=1, type_of_var=integer_t)
960 CALL section_add_keyword(section, keyword)
961 CALL keyword_release(keyword)
962
963 CALL keyword_create(keyword, __location__, name="RADIUS", &
964 description="Overwrite the specification of the radius only for the MM atom involved in the link. "// &
965 "Default is to use the same radius as for the specified type.", &
966 usage="RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom")
967 CALL section_add_keyword(section, keyword)
968 CALL keyword_release(keyword)
969
970 CALL keyword_create( &
971 keyword, __location__, name="CORR_RADIUS", &
972 description="Overwrite the specification of the correction radius only for the MM atom involved in the link. "// &
973 "Default is to use the same correction radius as for the specified type.", &
974 usage="CORR_RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom")
975 CALL section_add_keyword(section, keyword)
976 CALL keyword_release(keyword)
977
978 CALL keyword_create(keyword, __location__, name="LINK_TYPE", &
979 variants=["LINK ", "TYPE ", "LTYPE"], &
980 description="Specifies the method to use to treat the defined QM/MM link", &
981 usage="LINK_TYPE char", &
982 enum_c_vals=s2a("IMOMM", "GHO", "PSEUDO"), &
984 enum_desc=s2a("Use Integrated Molecular Orbital Molecular Mechanics method", &
985 "Use Generalized Hybrid Orbital method", &
986 "Use a monovalent pseudo-potential"), &
987 default_i_val=do_qmmm_link_imomm)
988 CALL section_add_keyword(section, keyword)
989 CALL keyword_release(keyword)
990
991 CALL keyword_create(keyword, __location__, name="ALPHA_IMOMM", &
992 variants=s2a("ALPHA"), &
993 description="Specifies the scaling factor to be used for projecting the forces "// &
994 "on the capping hydrogen in the IMOMM QM/MM link scheme to the MM atom of the link. "// &
995 "A good guess can be derived from the bond distances of the forcefield: "// &
996 "alpha = r_eq(QM-MM) / r_eq(QM-H).", &
997 usage="ALPHA_IMOMM real", n_var=1, type_of_var=real_t, &
998 default_r_val=alpha_imomm_default)
999 CALL section_add_keyword(section, keyword)
1000 CALL keyword_release(keyword)
1001
1002 CALL keyword_create(keyword, __location__, name="QMMM_SCALE_FACTOR", &
1003 variants=["QMMM_CHARGE_SCALE ", &
1004 "QMMM_CHARGE_FACTOR", &
1005 "QMMM_SCALE_CHARGE "], &
1006 description="Specifies the scaling factor for the MM charge involved in the link QM/MM."// &
1007 " This keyword affects only the QM/MM potential, it doesn't affect the electrostatic in"// &
1008 " the classical part of the code."// &
1009 " Default 1.0 i.e. no charge rescaling of the MM atom of the QM/MM link bond.", &
1010 usage="QMMM_SCALE_FACTOR real", n_var=1, type_of_var=real_t, &
1011 default_r_val=charge_scale_factor)
1012 CALL section_add_keyword(section, keyword)
1013 CALL keyword_release(keyword)
1014
1015 CALL keyword_create(keyword, __location__, name="FIST_SCALE_FACTOR", &
1016 variants=["FIST_CHARGE_SCALE ", &
1017 "FIST_CHARGE_FACTOR", &
1018 "FIST_SCALE_CHARGE "], &
1019 description="Specifies the scaling factor for the MM charge involved in the link QM/MM."// &
1020 " This keyword modifies the MM charge in FIST. The modified charge will be used then also"// &
1021 " for the generation of the QM/MM potential. "// &
1022 "Default 1.0 i.e. no charge rescaling of the MM atom of the QM/MM link bond.", &
1023 usage="FIST_SCALE_FACTOR real", n_var=1, type_of_var=real_t, &
1024 default_r_val=charge_scale_factor)
1025 CALL section_add_keyword(section, keyword)
1026 CALL keyword_release(keyword)
1027
1028 CALL section_create(subsection, __location__, name="MOVE_MM_CHARGE", &
1029 description="Specify information to move a classical charge before the"// &
1030 " QM/MM energies and forces evaluation", &
1031 n_keywords=4, n_subsections=0, repeats=.true.)
1032
1033 CALL keyword_create(keyword, __location__, name="ATOM_INDEX_1", &
1034 variants=["MM1"], &
1035 description="Specifies the index of the MM atom involved in the QM/MM link to be moved", &
1036 usage="ATOM_INDEX_1 integer", n_var=1, type_of_var=integer_t)
1037 CALL section_add_keyword(subsection, keyword)
1038 CALL keyword_release(keyword)
1039
1040 CALL keyword_create(keyword, __location__, name="ATOM_INDEX_2", &
1041 variants=["MM2"], &
1042 description="Specifies the index of the second atom defining the direction along which"// &
1043 " the atom will be moved", &
1044 usage="ATOM_INDEX_2 integer", n_var=1, type_of_var=integer_t)
1045 CALL section_add_keyword(subsection, keyword)
1046 CALL keyword_release(keyword)
1047
1048 CALL keyword_create(keyword, __location__, name="ALPHA", &
1049 description="Specifies the scaling factor that defines the movement along the defined direction", &
1050 usage="ALPHA real", n_var=1, type_of_var=real_t)
1051 CALL section_add_keyword(subsection, keyword)
1052 CALL keyword_release(keyword)
1053
1054 CALL keyword_create(keyword, __location__, name="RADIUS", &
1055 description="Specifies the radius used for the QM/MM electrostatic coupling after movement", &
1056 usage="RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom", default_r_val=0.0_dp)
1057 CALL section_add_keyword(subsection, keyword)
1058 CALL keyword_release(keyword)
1059
1060 CALL keyword_create(keyword, __location__, name="CORR_RADIUS", &
1061 description="Specifies the correction radius used for the QM/MM electrostatic coupling after movement", &
1062 usage="CORR_RADIUS real", n_var=1, type_of_var=real_t, unit_str="angstrom", default_r_val=0.0_dp)
1063 CALL section_add_keyword(subsection, keyword)
1064 CALL keyword_release(keyword)
1065
1066 CALL section_add_subsection(section, subsection)
1067 CALL section_release(subsection)
1068
1069 CALL section_create(subsection, __location__, name="ADD_MM_CHARGE", &
1070 description="Specify information to add a classical charge before the"// &
1071 " QM/MM energies and forces evaluation", &
1072 n_keywords=5, n_subsections=0, repeats=.true.)
1073
1074 CALL keyword_create(keyword, __location__, name="ATOM_INDEX_1", &
1075 variants=["MM1"], &
1076 description="Specifies the index of the first atom defining the direction along which"// &
1077 " the atom will be added", &
1078 usage="ATOM_INDEX_1 integer", n_var=1, type_of_var=integer_t)
1079 CALL section_add_keyword(subsection, keyword)
1080 CALL keyword_release(keyword)
1081
1082 CALL keyword_create(keyword, __location__, name="ATOM_INDEX_2", &
1083 variants=["MM2"], &
1084 description="Specifies the index of the second atom defining the direction along which"// &
1085 " the atom will be added", &
1086 usage="ATOM_INDEX_2 integer", n_var=1, type_of_var=integer_t)
1087 CALL section_add_keyword(subsection, keyword)
1088 CALL keyword_release(keyword)
1089
1090 CALL keyword_create(keyword, __location__, name="ALPHA", &
1091 description="Specifies the scaling factor that defines the movement along the defined direction", &
1092 usage="ALPHA real", n_var=1, type_of_var=real_t)
1093 CALL section_add_keyword(subsection, keyword)
1094 CALL keyword_release(keyword)
1095
1096 CALL keyword_create(keyword, __location__, name="RADIUS", &
1097 description="Specifies the radius used for the QM/MM electrostatic coupling for the added source", &
1098 usage="RADIUS real", n_var=1, unit_str="angstrom", &
1099 default_r_val=cp_unit_to_cp2k(radius_qmmm_default, "angstrom"))
1100 CALL section_add_keyword(subsection, keyword)
1101 CALL keyword_release(keyword)
1102
1103 CALL keyword_create( &
1104 keyword, __location__, name="CORR_RADIUS", &
1105 description="Specifies the correction radius used for the QM/MM electrostatic coupling for the added source", &
1106 usage="CORR_RADIUS real", n_var=1, unit_str="angstrom", &
1107 default_r_val=cp_unit_to_cp2k(radius_qmmm_default, "angstrom"))
1108 CALL section_add_keyword(subsection, keyword)
1109 CALL keyword_release(keyword)
1110
1111 CALL keyword_create(keyword, __location__, name="CHARGE", &
1112 description="Specifies the charge for the added source of QM/MM potential", &
1113 usage="CHARGE real", default_r_val=0.0_dp, n_var=1, type_of_var=real_t)
1114 CALL section_add_keyword(subsection, keyword)
1115 CALL keyword_release(keyword)
1116
1117 CALL section_add_subsection(section, subsection)
1118 CALL section_release(subsection)
1119 END SUBROUTINE create_qmmm_link_section
1120
1121! **************************************************************************************************
1122!> \brief creates the interpolation section
1123!> \param section ...
1124!> \author tlaino
1125! **************************************************************************************************
1126 SUBROUTINE create_qmmm_interp_section(section)
1127 TYPE(section_type), POINTER :: section
1128
1129 TYPE(keyword_type), POINTER :: keyword
1130 TYPE(section_type), POINTER :: print_key
1131
1132 cpassert(.NOT. ASSOCIATED(section))
1133 CALL section_create(section, __location__, name="interpolator", &
1134 description="kind of interpolation used between the multigrids", &
1135 n_keywords=5, n_subsections=0, repeats=.false.)
1136
1137 NULLIFY (keyword, print_key)
1138
1139 CALL keyword_create(keyword, __location__, name="kind", &
1140 description="the interpolator to use", &
1141 usage="kind spline3", &
1142 default_i_val=spline3_nopbc_interp, &
1143 enum_c_vals=s2a("spline3_nopbc"), &
1144 enum_i_vals=[spline3_nopbc_interp])
1145 CALL section_add_keyword(section, keyword)
1146 CALL keyword_release(keyword)
1147
1148 CALL keyword_create(keyword, __location__, name="safe_computation", &
1149 description="if a non unrolled calculation is to be performed in parallel", &
1150 usage="safe_computation OFF", &
1151 default_l_val=.false., &
1152 lone_keyword_l_val=.true.)
1153 CALL section_add_keyword(section, keyword)
1154 CALL keyword_release(keyword)
1155
1156 CALL keyword_create(keyword, __location__, name="aint_precond", &
1157 description="the approximate inverse to use to get the starting point"// &
1158 " for the linear solver of the spline3 methods", &
1159 usage="aint_precond copy", &
1160 default_i_val=precond_spl3_aint, &
1161 enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_precond1", &
1162 "spl3_nopbc_aint2", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1165 CALL section_add_keyword(section, keyword)
1166 CALL keyword_release(keyword)
1167
1168 CALL keyword_create(keyword, __location__, name="precond", &
1169 description="The preconditioner used"// &
1170 " for the linear solver of the spline3 methods", &
1171 usage="precond copy", &
1172 default_i_val=precond_spl3_3, &
1173 enum_c_vals=s2a("copy", "spl3_nopbc_aint1", "spl3_nopbc_precond1", &
1174 "spl3_nopbc_aint2", "spl3_nopbc_precond2", "spl3_nopbc_precond3"), &
1177 CALL section_add_keyword(section, keyword)
1178 CALL keyword_release(keyword)
1179
1180 CALL keyword_create(keyword, __location__, name="eps_x", &
1181 description="accuracy on the solution for spline3 the interpolators", &
1182 usage="eps_x 1.e-15", default_r_val=1.e-10_dp)
1183 CALL section_add_keyword(section, keyword)
1184 CALL keyword_release(keyword)
1185
1186 CALL keyword_create(keyword, __location__, name="eps_r", &
1187 description="accuracy on the residual for spline3 the interpolators", &
1188 usage="eps_r 1.e-15", default_r_val=1.e-10_dp)
1189 CALL section_add_keyword(section, keyword)
1190 CALL keyword_release(keyword)
1191
1192 CALL keyword_create(keyword, __location__, name="max_iter", &
1193 variants=['maxiter'], &
1194 description="the maximum number of iterations", &
1195 usage="max_iter 200", default_i_val=100)
1196 CALL section_add_keyword(section, keyword)
1197 CALL keyword_release(keyword)
1198
1199 NULLIFY (print_key)
1200 CALL cp_print_key_section_create(print_key, __location__, "conv_info", &
1201 description="if convergence information about the linear solver"// &
1202 " of the spline methods should be printed", &
1203 print_level=medium_print_level, each_iter_names=s2a("SPLINE_FIND_COEFFS"), &
1204 each_iter_values=[10], filename="__STD_OUT__", &
1205 add_last=add_last_numeric)
1206 CALL section_add_subsection(section, print_key)
1207 CALL section_release(print_key)
1208
1209 CALL cp_print_key_section_create(print_key, __location__, "spl_coeffs", &
1210 description="outputs a cube with the coefficients calculated for "// &
1211 "the spline interpolation", &
1212 print_level=debug_print_level)
1213 CALL section_add_subsection(section, print_key)
1214 CALL section_release(print_key)
1215 END SUBROUTINE create_qmmm_interp_section
1216
1217! **************************************************************************************************
1218!> \brief Create the print qmmm section
1219!> \param section the section to create
1220!> \author teo
1221! **************************************************************************************************
1222 SUBROUTINE create_print_qmmm_section(section)
1223 TYPE(section_type), POINTER :: section
1224
1225 TYPE(keyword_type), POINTER :: keyword
1226 TYPE(section_type), POINTER :: print_key
1227
1228 cpassert(.NOT. ASSOCIATED(section))
1229 NULLIFY (keyword, print_key)
1230 CALL section_create(section, __location__, name="print", &
1231 description="Section of possible print options specific of the QMMM code.", &
1232 n_keywords=0, n_subsections=1, repeats=.false.)
1233
1234 NULLIFY (print_key)
1235
1236 CALL cp_print_key_section_create(print_key, __location__, "DIPOLE", &
1237 description="Controls the printing of the DIPOLE in a QM/MM calculations."// &
1238 " It requires that the DIPOLE calculations is"// &
1239 " requested both for the QS and for the MM part.", &
1240 print_level=high_print_level, filename="__STD_OUT__")
1241 CALL section_add_subsection(section, print_key)
1242 CALL section_release(print_key)
1243
1244 CALL cp_print_key_section_create(print_key, __location__, "PGF", &
1245 description="Controls the printing of the gaussian expansion basis set of the"// &
1246 " electrostatic potential", &
1247 print_level=high_print_level, filename="__STD_OUT__")
1248 CALL section_add_subsection(section, print_key)
1249 CALL section_release(print_key)
1250
1251 CALL cp_print_key_section_create(print_key, __location__, "POTENTIAL", &
1252 description="Controls the printing of the QMMM potential", &
1253 print_level=high_print_level, filename="MM_ELPOT_QMMM", &
1254 common_iter_levels=1)
1255
1256 CALL keyword_create(keyword, __location__, name="stride", &
1257 description="The stride (X,Y,Z) used to write the cube file "// &
1258 "(larger values result in smaller cube files). You can provide 3 numbers (for X,Y,Z) or"// &
1259 " 1 number valid for all components.", &
1260 usage="STRIDE 2 2 2", n_var=-1, default_i_vals=[2, 2, 2], type_of_var=integer_t)
1261 CALL section_add_keyword(print_key, keyword)
1262 CALL keyword_release(keyword)
1263
1264 CALL section_add_subsection(section, print_key)
1265 CALL section_release(print_key)
1266
1267 CALL cp_print_key_section_create(print_key, __location__, "MM_POTENTIAL", &
1268 description="Controls the printing of the MM unidimensional potential on file", &
1269 print_level=high_print_level, filename="MM_ELPOT", &
1270 common_iter_levels=1)
1271 CALL section_add_subsection(section, print_key)
1272 CALL section_release(print_key)
1273
1274 CALL cp_print_key_section_create(print_key, __location__, "QMMM_MATRIX", &
1275 description="Controls the printing of the QMMM 1 electron Hamiltonian Matrix"// &
1276 " for methods like semiempirical and DFTB", &
1277 print_level=high_print_level, filename="__STD_OUT__", &
1278 common_iter_levels=1)
1279 CALL section_add_subsection(section, print_key)
1280 CALL section_release(print_key)
1281
1282 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_BANNER", &
1283 description="Controls the printing of the banner of the MM program", &
1284 print_level=silent_print_level, filename="__STD_OUT__")
1285 CALL section_add_subsection(section, print_key)
1286 CALL section_release(print_key)
1287
1288 CALL cp_print_key_section_create(print_key, __location__, "PROGRAM_RUN_INFO", &
1289 description="Controls the printing of information regarding the run.", &
1290 print_level=medium_print_level, filename="__STD_OUT__")
1291 CALL section_add_subsection(section, print_key)
1292 CALL section_release(print_key)
1293
1295 print_key, __location__, "PERIODIC_INFO", &
1296 description="Controls the printing of information regarding the periodic boundary condition.", &
1297 print_level=medium_print_level, filename="__STD_OUT__")
1298 CALL section_add_subsection(section, print_key)
1299 CALL section_release(print_key)
1300
1301 CALL cp_print_key_section_create(print_key, __location__, "GRID_INFORMATION", &
1302 description="Controls the printing of information regarding the PW grid structures"// &
1303 " for PERIODIC QM/MM calculations.", &
1304 print_level=medium_print_level, filename="__STD_OUT__")
1305 CALL section_add_subsection(section, print_key)
1306 CALL section_release(print_key)
1307
1308 CALL cp_print_key_section_create(print_key, __location__, "derivatives", &
1309 description="Print all derivatives after QM/MM calculation", &
1310 print_level=high_print_level, filename="__STD_OUT__")
1311 CALL section_add_subsection(section, print_key)
1312 CALL section_release(print_key)
1313
1314 CALL cp_print_key_section_create(print_key, __location__, "qmmm_charges", &
1315 description="Print all charges generating the QM/MM potential", &
1316 print_level=medium_print_level, filename="__STD_OUT__")
1317 CALL section_add_subsection(section, print_key)
1318 CALL section_release(print_key)
1319
1320 CALL cp_print_key_section_create(print_key, __location__, "qmmm_link_info", &
1321 description="Print all information on QM/MM links", &
1322 print_level=medium_print_level, filename="__STD_OUT__")
1323 CALL section_add_subsection(section, print_key)
1324 CALL section_release(print_key)
1325
1326 CALL cp_print_key_section_create(print_key, __location__, "qs_derivatives", &
1327 description="Print QM derivatives after QS calculation", &
1328 print_level=medium_print_level, filename="__STD_OUT__")
1329 CALL section_add_subsection(section, print_key)
1330 CALL section_release(print_key)
1331
1332 CALL cp_print_key_section_create(print_key, __location__, "image_charge_info", &
1333 description="Prints image charge coefficients and detailed energy info", &
1334 print_level=high_print_level, filename="__STD_OUT__")
1335 CALL section_add_subsection(section, print_key)
1336 CALL section_release(print_key)
1337
1338 CALL cp_print_key_section_create(print_key, __location__, "image_charge_restart", &
1339 description="Controls the printing of the restart file for "// &
1340 "the image matrix when using the iterative scheme", &
1341 print_level=low_print_level, add_last=add_last_numeric, filename="RESTART", &
1342 common_iter_levels=3)
1343 CALL section_add_subsection(section, print_key)
1344 CALL section_release(print_key)
1345
1346 END SUBROUTINE create_print_qmmm_section
1347
1348END MODULE input_cp2k_qmmm
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public bernstein2012
integer, save, public bernstein2009
integer, save, public laino2006
integer, save, public golze2013
integer, save, public laino2005
Handles all functions related to the CELL.
Definition cell_types.F:15
integer, parameter, public use_perd_none
Definition cell_types.F:42
Interface to Minimax-Ewald method for periodic ERI's to be used in CP2K.
subroutine, public create_eri_mme_section(section, default_n_minimax)
Create main input section.
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer, parameter, public debug_print_level
integer, parameter, public low_print_level
integer, parameter, public medium_print_level
integer, parameter, public high_print_level
integer, parameter, public add_last_numeric
integer, parameter, public silent_print_level
subroutine, public cp_print_key_section_create(print_key_section, location, name, description, print_level, each_iter_names, each_iter_values, add_last, filename, common_iter_levels, citations, unit_str)
creates a print_key section
utils to manipulate splines on the regular grid of a pw
integer, parameter, public spline3_nopbc_interp
unit conversion facility
Definition cp_units.F:30
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Definition cp_units.F:1149
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_par_atom
integer, parameter, public do_qmmm_image_calcmatrix
integer, parameter, public do_qmmm_none
integer, parameter, public do_qmmm_center_every_step
integer, parameter, public do_qmmm_pcharge
real(kind=dp), parameter, public alpha_imomm_default
integer, parameter, public do_fm_mom_conserv_none
integer, parameter, public do_qmmm_image_iter
integer, parameter, public do_fm_mom_conserv_buffer
integer, parameter, public do_qmmm_link_pseudo
integer, parameter, public do_qmmm_center_pbc_aware
integer, parameter, public do_eri_mme
integer, parameter, public do_fm_mom_conserv_equal_f
real(kind=dp), parameter, public radius_qmmm_default
integer, parameter, public do_qmmm_coulomb
integer, parameter, public do_qmmm_center_never
integer, parameter, public do_multipole_section_off
integer, parameter, public do_qmmm_link_gho
integer, parameter, public do_qmmm_wall_quadratic
integer, parameter, public do_qmmm_wall_reflective
integer, parameter, public do_qmmm_swave
real(kind=dp), parameter, public charge_scale_factor
integer, parameter, public do_fm_mom_conserv_qm
integer, parameter, public do_qmmm_center_setup_only
integer, parameter, public do_multipole_section_on
integer, parameter, public do_fm_mom_conserv_equal_a
integer, parameter, public do_fm_mom_conserv_core
integer, parameter, public do_qmmm_gauss
integer, parameter, public gaussian
integer, parameter, public do_par_grid
integer, parameter, public do_qmmm_wall_none
integer, parameter, public do_qmmm_link_imomm
integer, parameter, public do_qmmm_center_max_minus_min
integer, parameter, public do_eri_gpw
creates the mm section of the input
subroutine, public create_genpot_section(section)
This section specifies the input parameters for a generic potential form.
subroutine, public create_williams_section(section)
This section specifies the input parameters for Williams potential type.
subroutine, public create_goodwin_section(section)
This section specifies the input parameters for Goodwin potential type.
subroutine, public create_nonbonded14_section(section)
This section specifies the input parameters for 1-4 NON-BONDED Interactions.
subroutine, public create_lj_section(section)
This section specifies the input parameters for Lennard-Jones potential type.
function that build the poisson section of the input
subroutine, public create_poisson_section(section)
Creates the Poisson section.
subroutine, public create_gspace_interp_section(section)
creates the interpolation section for the periodic QM/MM
creates the qmmm section of the input
subroutine, public create_qmmm_section(section)
Creates the QM/MM section.
builds the subsystem section of the input
subroutine, public create_cell_section(section, periodic)
creates the cell section
represents keywords in an input
subroutine, public keyword_release(keyword)
releases the given keyword (see doc/ReferenceCounting.html)
subroutine, public keyword_create(keyword, location, name, description, usage, type_of_var, n_var, repeats, variants, default_val, default_l_val, default_r_val, default_lc_val, default_c_val, default_i_val, default_l_vals, default_r_vals, default_c_vals, default_i_vals, lone_keyword_val, lone_keyword_l_val, lone_keyword_r_val, lone_keyword_c_val, lone_keyword_i_val, lone_keyword_l_vals, lone_keyword_r_vals, lone_keyword_c_vals, lone_keyword_i_vals, enum_c_vals, enum_i_vals, enum, enum_strict, enum_desc, unit_str, citations, deprecation_notice, removed)
creates a keyword object
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_create(section, location, name, description, n_keywords, n_subsections, repeats, citations, deprecation_notice)
creates a list of keywords
subroutine, public section_add_keyword(section, keyword)
adds a keyword to the given section
subroutine, public section_add_subsection(section, subsection)
adds a subsection to the given section
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
a wrapper for basic fortran types.
integer, parameter, public real_t
integer, parameter, public lchar_t
integer, parameter, public logical_t
integer, parameter, public char_t
integer, parameter, public integer_t
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
different utils that are useful to manipulate splines on the regular grid of a pw
integer, parameter, public precond_spl3_3
integer, parameter, public precond_spl3_aint
integer, parameter, public no_precond
integer, parameter, public precond_spl3_2
integer, parameter, public precond_spl3_aint2
integer, parameter, public precond_spl3_1
Utilities for string manipulations.
represent a keyword in the input
represent a section of the input file