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