(git:ca6b8f8)
Loading...
Searching...
No Matches
force_fields_input.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!> \par History
10!> Subroutine input_torsions changed (DG) 05-Dec-2000
11!> Output formats changed (DG) 05-Dec-2000
12!> JGH (26-01-2002) : force field parameters stored in tables, not in
13!> matrices. Input changed to have parameters labeled by the position
14!> and not atom pairs (triples etc)
15!> Teo (11.2005) : Moved all information on force field pair_potential to
16!> a much lighter memory structure
17!> Teo 09.2006 : Split all routines force_field I/O in a separate file
18!> \author CJM
19! **************************************************************************************************
23 USE bibliography, ONLY: clabaut2020,&
27 tosi1964a,&
28 tosi1964b,&
30 cite_reference
31 USE cp_files, ONLY: discover_file
43 USE cp_units, ONLY: cp_unit_to_cp2k
47 do_ff_g87,&
48 do_ff_g96,&
60 USE input_val_types, ONLY: val_get,&
62 USE kinds, ONLY: default_path_length,&
64 dp
65 USE mathconstants, ONLY: pi
66 USE mathlib, ONLY: invert_matrix
69 USE pair_potential_types, ONLY: &
79 USE torch_api, ONLY: torch_allow_tf32,&
81#include "./base/base_uses.f90"
82
83 IMPLICIT NONE
84
85 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'force_fields_input'
86
87 PRIVATE
88 PUBLIC :: read_force_field_section, &
94
95CONTAINS
96
97! **************************************************************************************************
98!> \brief Reads the force_field input section
99!> \param ff_section ...
100!> \param mm_section ...
101!> \param ff_type ...
102!> \param para_env ...
103!> \author teo
104! **************************************************************************************************
105 SUBROUTINE read_force_field_section1(ff_section, mm_section, ff_type, para_env)
106 TYPE(section_vals_type), POINTER :: ff_section, mm_section
107 TYPE(force_field_type), INTENT(INOUT) :: ff_type
108 TYPE(mp_para_env_type), POINTER :: para_env
109
110 CHARACTER(LEN=default_string_length), &
111 DIMENSION(:), POINTER :: atm_names
112 INTEGER :: nace, nallegro, nb4, nbends, nbm, nbmhft, nbmhftd, nbonds, nchg, ndeepmd, neam, &
113 ngal, ngal21, ngd, ngp, nimpr, nipbv, nlj, nnequip, nopbend, nshell, nsiepmann, ntab, &
114 ntersoff, ntors, ntot, nubs, nwl
115 LOGICAL :: explicit, unique_spline
116 REAL(KIND=dp) :: min_eps_spline_allowed
117 TYPE(input_info_type), POINTER :: inp_info
118 TYPE(section_vals_type), POINTER :: tmp_section, tmp_section2
119
120 INTEGER::i
121
122 NULLIFY (tmp_section, tmp_section2)
123 inp_info => ff_type%inp_info
124 CALL section_vals_val_get(ff_section, "PARMTYPE", i_val=ff_type%ff_type)
125 CALL section_vals_val_get(ff_section, "EI_SCALE14", r_val=ff_type%ei_scale14)
126 CALL section_vals_val_get(ff_section, "VDW_SCALE14", r_val=ff_type%vdw_scale14)
127 CALL section_vals_val_get(ff_section, "SPLINE%RCUT_NB", r_val=ff_type%rcut_nb)
128 CALL section_vals_val_get(ff_section, "SPLINE%R0_NB", r_val=ff_type%rlow_nb)
129 CALL section_vals_val_get(ff_section, "SPLINE%EPS_SPLINE", r_val=ff_type%eps_spline)
130 CALL section_vals_val_get(ff_section, "SPLINE%EMAX_SPLINE", r_val=ff_type%emax_spline)
131 CALL section_vals_val_get(ff_section, "SPLINE%EMAX_ACCURACY", r_val=ff_type%max_energy)
132 CALL section_vals_val_get(ff_section, "SPLINE%NPOINTS", i_val=ff_type%npoints)
133 CALL section_vals_val_get(ff_section, "IGNORE_MISSING_CRITICAL_PARAMS", l_val=ff_type%ignore_missing_critical)
134 cpassert(ff_type%max_energy <= ff_type%emax_spline)
135 ! Read the parameter file name only if the force field type requires it..
136 SELECT CASE (ff_type%ff_type)
138 CALL section_vals_val_get(ff_section, "PARM_FILE_NAME", c_val=ff_type%ff_file_name)
139 IF (trim(ff_type%ff_file_name) == "") &
140 cpabort("Force Field Parameter's filename is empty! Please check your input file.")
141 CASE (do_ff_undef)
142 ! Do Nothing
143 CASE DEFAULT
144 cpabort("Force field type not implemented")
145 END SELECT
146 ! Numerical Accuracy:
147 ! the factors here should depend on the energy and on the shape of each potential mapped
148 ! with splines. this would make everything un-necessarily complicated. Let's just be safe
149 ! and assume that we cannot achieve an accuracy on the spline 2 orders of magnitude more
150 ! than the smallest representable number (taking into account also the max_energy for the
151 ! spline generation
152 min_eps_spline_allowed = 20.0_dp*max(ff_type%max_energy, 10.0_dp)*epsilon(0.0_dp)
153 IF (ff_type%eps_spline < min_eps_spline_allowed) THEN
154 CALL cp_warn(__location__, &
155 "Requested spline accuracy ("//trim(cp_to_string(ff_type%eps_spline))//" ) "// &
156 "is smaller than the minimum value allowed ("//trim(cp_to_string(min_eps_spline_allowed))// &
157 " ) with the present machine precision ("//trim(cp_to_string(epsilon(0.0_dp)))//" ). "// &
158 "New EPS_SPLINE value ("//trim(cp_to_string(min_eps_spline_allowed))//" ). ")
159 ff_type%eps_spline = min_eps_spline_allowed
160 END IF
161 CALL section_vals_val_get(ff_section, "SHIFT_CUTOFF", l_val=ff_type%shift_cutoff)
162 CALL section_vals_val_get(ff_section, "SPLINE%UNIQUE_SPLINE", l_val=unique_spline)
163 ! Single spline
166
167 CALL section_vals_val_get(ff_section, "MULTIPLE_POTENTIAL", l_val=ff_type%multiple_potential)
168 CALL section_vals_val_get(ff_section, "DO_NONBONDED", l_val=ff_type%do_nonbonded)
169 CALL section_vals_val_get(ff_section, "DO_ELECTROSTATICS", l_val=ff_type%do_electrostatics)
170 tmp_section => section_vals_get_subs_vals(ff_section, "NONBONDED")
171 CALL section_vals_get(tmp_section, explicit=explicit)
172 IF (explicit .AND. ff_type%do_nonbonded) THEN
173 tmp_section2 => section_vals_get_subs_vals(tmp_section, "LENNARD-JONES")
174 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nlj)
175 ntot = 0
176 IF (explicit) THEN
177 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nlj, lj_charmm=.true.)
178 CALL read_lj_section(inp_info%nonbonded, tmp_section2, ntot)
179 END IF
180
181 tmp_section2 => section_vals_get_subs_vals(tmp_section, "WILLIAMS")
182 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nwl)
183 ntot = nlj
184 IF (explicit) THEN
185 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nwl, williams=.true.)
186 CALL read_wl_section(inp_info%nonbonded, tmp_section2, ntot)
187 END IF
188
189 tmp_section2 => section_vals_get_subs_vals(tmp_section, "EAM")
190 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=neam)
191 ntot = nlj + nwl
192 IF (explicit) THEN
193 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + neam, eam=.true.)
194 CALL read_eam_section(inp_info%nonbonded, tmp_section2, ntot, para_env, mm_section)
195 END IF
196
197 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GOODWIN")
198 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngd)
199 ntot = nlj + nwl + neam
200 IF (explicit) THEN
201 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ngd, goodwin=.true.)
202 CALL read_gd_section(inp_info%nonbonded, tmp_section2, ntot)
203 END IF
204
205 tmp_section2 => section_vals_get_subs_vals(tmp_section, "IPBV")
206 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nipbv)
207 ntot = nlj + nwl + neam + ngd
208 IF (explicit) THEN
209 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nipbv, ipbv=.true.)
210 CALL read_ipbv_section(inp_info%nonbonded, tmp_section2, ntot)
211 END IF
212
213 tmp_section2 => section_vals_get_subs_vals(tmp_section, "BMHFT")
214 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nbmhft)
215 ntot = nlj + nwl + neam + ngd + nipbv
216 IF (explicit) THEN
217 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nbmhft, bmhft=.true.)
218 CALL read_bmhft_section(inp_info%nonbonded, tmp_section2, ntot)
219 END IF
220
221 tmp_section2 => section_vals_get_subs_vals(tmp_section, "BMHFTD")
222 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nbmhftd)
223 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft
224 IF (explicit) THEN
225 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nbmhftd, bmhftd=.true.)
226 CALL read_bmhftd_section(inp_info%nonbonded, tmp_section2, ntot)
227 END IF
228
229 tmp_section2 => section_vals_get_subs_vals(tmp_section, "BUCK4RANGES")
230 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nb4)
231 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd
232 IF (explicit) THEN
233 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nb4, buck4r=.true.)
234 CALL read_b4_section(inp_info%nonbonded, tmp_section2, ntot)
235 END IF
236
237 tmp_section2 => section_vals_get_subs_vals(tmp_section, "BUCKMORSE")
238 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nbm)
239 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4
240 IF (explicit) THEN
241 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nbm, buckmo=.true.)
242 CALL read_bm_section(inp_info%nonbonded, tmp_section2, ntot)
243 END IF
244
245 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GENPOT")
246 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngp)
247 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm
248 IF (explicit) THEN
249 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ngp, gp=.true.)
250 CALL read_gp_section(inp_info%nonbonded, tmp_section2, ntot)
251 END IF
252 tmp_section2 => section_vals_get_subs_vals(tmp_section, "TERSOFF")
253 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ntersoff)
254 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp
255 IF (explicit) THEN
256 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ntersoff, tersoff=.true.)
257 CALL read_tersoff_section(inp_info%nonbonded, tmp_section2, ntot, tmp_section2)
258 END IF
259
260 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GAL19")
261 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngal)
262 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff
263 IF (explicit) THEN
264 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ngal, gal=.true.)
265 CALL read_gal_section(inp_info%nonbonded, tmp_section2, ntot, tmp_section2)
266 END IF
267
268 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GAL21")
269 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngal21)
270 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + ngal
271 IF (explicit) THEN
272 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ngal21, gal21=.true.)
273 CALL read_gal21_section(inp_info%nonbonded, tmp_section2, ntot, tmp_section2)
274 END IF
275
276 tmp_section2 => section_vals_get_subs_vals(tmp_section, "SIEPMANN")
277 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nsiepmann)
278 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + ngal + ngal21
279 IF (explicit) THEN
280 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nsiepmann, siepmann=.true.)
281 CALL read_siepmann_section(inp_info%nonbonded, tmp_section2, ntot, tmp_section2)
282 END IF
283
284 tmp_section2 => section_vals_get_subs_vals(tmp_section, "nequip")
285 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nnequip)
286 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + &
287 ngal + ngal21 + nsiepmann
288 IF (explicit) THEN
289 ! avoid repeating the nequip section for each pair
290 CALL section_vals_val_get(tmp_section2, "ATOMS", c_vals=atm_names)
291 nnequip = nnequip - 1 + SIZE(atm_names) + (SIZE(atm_names)*SIZE(atm_names) - SIZE(atm_names))/2
292 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nnequip, nequip=.true.)
293 CALL read_nequip_section(inp_info%nonbonded, tmp_section2, ntot)
294 END IF
295
296 tmp_section2 => section_vals_get_subs_vals(tmp_section, "allegro")
297 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nallegro)
298 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + &
299 ngal + ngal21 + nsiepmann + nnequip
300 IF (explicit) THEN
301 ! avoid repeating the allegro section for each pair
302 CALL section_vals_val_get(tmp_section2, "ATOMS", c_vals=atm_names)
303 nallegro = nallegro - 1 + SIZE(atm_names) + (SIZE(atm_names)*SIZE(atm_names) - SIZE(atm_names))/2
304 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nallegro, allegro=.true.)
305 CALL read_allegro_section(inp_info%nonbonded, tmp_section2, ntot)
306 END IF
307
308 tmp_section2 => section_vals_get_subs_vals(tmp_section, "TABPOT")
309 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ntab)
310 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + &
311 ngal + ngal21 + nsiepmann + nnequip + nallegro
312 IF (explicit) THEN
313 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ntab, tab=.true.)
314 CALL read_tabpot_section(inp_info%nonbonded, tmp_section2, ntot, para_env, mm_section)
315 END IF
316
317 tmp_section2 => section_vals_get_subs_vals(tmp_section, "DEEPMD")
318 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ndeepmd)
319 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + &
320 ngal + ngal21 + nsiepmann + nnequip + nallegro + ntab
321 IF (explicit) THEN
322 ! avoid repeating the deepmd section for each pair
323 CALL section_vals_val_get(tmp_section2, "ATOMS", c_vals=atm_names)
324 ndeepmd = ndeepmd - 1 + SIZE(atm_names) + (SIZE(atm_names)*SIZE(atm_names) - SIZE(atm_names))/2
325 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + ndeepmd, deepmd=.true.)
326 CALL read_deepmd_section(inp_info%nonbonded, tmp_section2, ntot)
327 END IF
328
329 tmp_section2 => section_vals_get_subs_vals(tmp_section, "ACE")
330 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nace)
331 ntot = nlj + nwl + neam + ngd + nipbv + nbmhft + nbmhftd + nb4 + nbm + ngp + ntersoff + &
332 ngal + ngal21 + nsiepmann + nnequip + nallegro + ntab + ndeepmd
333 IF (explicit) THEN
334 ! avoid repeating the ace section for each pair
335 CALL section_vals_val_get(tmp_section2, "ATOMS", c_vals=atm_names)
336 nace = nace - 1 + SIZE(atm_names) + (SIZE(atm_names)*SIZE(atm_names) - SIZE(atm_names))/2
337 CALL pair_potential_reallocate(inp_info%nonbonded, 1, ntot + nace, ace=.true.)
338 CALL read_ace_section(inp_info%nonbonded, tmp_section2, ntot)
339 END IF
340
341 END IF
342
343 tmp_section => section_vals_get_subs_vals(ff_section, "NONBONDED14")
344 CALL section_vals_get(tmp_section, explicit=explicit)
345 IF (explicit .AND. ff_type%do_nonbonded) THEN
346 tmp_section2 => section_vals_get_subs_vals(tmp_section, "LENNARD-JONES")
347 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nlj)
348 ntot = 0
349 IF (explicit) THEN
350 CALL pair_potential_reallocate(inp_info%nonbonded14, 1, ntot + nlj, lj_charmm=.true.)
351 CALL read_lj_section(inp_info%nonbonded14, tmp_section2, ntot)
352 END IF
353 tmp_section2 => section_vals_get_subs_vals(tmp_section, "WILLIAMS")
354 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=nwl)
355 ntot = nlj
356 IF (explicit) THEN
357 CALL pair_potential_reallocate(inp_info%nonbonded14, 1, ntot + nwl, williams=.true.)
358 CALL read_wl_section(inp_info%nonbonded14, tmp_section2, ntot)
359 END IF
360 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GOODWIN")
361 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngd)
362 ntot = nlj + nwl
363 IF (explicit) THEN
364 CALL pair_potential_reallocate(inp_info%nonbonded14, 1, ntot + ngd, goodwin=.true.)
365 CALL read_gd_section(inp_info%nonbonded14, tmp_section2, ntot)
366 END IF
367 tmp_section2 => section_vals_get_subs_vals(tmp_section, "GENPOT")
368 CALL section_vals_get(tmp_section2, explicit=explicit, n_repetition=ngp)
369 ntot = nlj + nwl + ngd
370 IF (explicit) THEN
371 CALL pair_potential_reallocate(inp_info%nonbonded14, 1, ntot + ngp, gp=.true.)
372 CALL read_gp_section(inp_info%nonbonded14, tmp_section2, ntot)
373 END IF
374 END IF
375
376 tmp_section => section_vals_get_subs_vals(ff_section, "CHARGE")
377 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nchg)
378 IF (explicit) THEN
379 ntot = 0
380 CALL reallocate(inp_info%charge_atm, 1, nchg)
381 CALL reallocate(inp_info%charge, 1, nchg)
382 CALL read_chrg_section(inp_info%charge_atm, inp_info%charge, tmp_section, ntot)
383 END IF
384 tmp_section => section_vals_get_subs_vals(ff_section, "DIPOLE")
385 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nchg)
386 IF (explicit) THEN
387 ntot = 0
388 CALL reallocate(inp_info%apol_atm, 1, nchg)
389 CALL reallocate(inp_info%apol, 1, nchg)
390 CALL read_apol_section(inp_info%apol_atm, inp_info%apol, inp_info%damping_list, &
391 tmp_section, ntot)
392 END IF
393 tmp_section => section_vals_get_subs_vals(ff_section, "QUADRUPOLE")
394 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nchg)
395 IF (explicit) THEN
396 ntot = 0
397 CALL reallocate(inp_info%cpol_atm, 1, nchg)
398 CALL reallocate(inp_info%cpol, 1, nchg)
399 CALL read_cpol_section(inp_info%cpol_atm, inp_info%cpol, tmp_section, ntot)
400 END IF
401 tmp_section => section_vals_get_subs_vals(ff_section, "SHELL")
402 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nshell)
403 IF (explicit) THEN
404 ntot = 0
405 CALL shell_p_create(inp_info%shell_list, nshell)
406 CALL read_shell_section(inp_info%shell_list, tmp_section, ntot)
407 END IF
408
409 tmp_section => section_vals_get_subs_vals(ff_section, "BOND")
410 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nbonds)
411 IF (explicit) THEN
412 ntot = 0
413 CALL reallocate(inp_info%bond_kind, 1, nbonds)
414 CALL reallocate(inp_info%bond_a, 1, nbonds)
415 CALL reallocate(inp_info%bond_b, 1, nbonds)
416 CALL reallocate(inp_info%bond_k, 1, 3, 1, nbonds)
417 CALL reallocate(inp_info%bond_r0, 1, nbonds)
418 CALL reallocate(inp_info%bond_cs, 1, nbonds)
419 CALL read_bonds_section(inp_info%bond_kind, inp_info%bond_a, inp_info%bond_b, inp_info%bond_k, &
420 inp_info%bond_r0, inp_info%bond_cs, tmp_section, ntot)
421 END IF
422 tmp_section => section_vals_get_subs_vals(ff_section, "BEND")
423 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nbends)
424 IF (explicit) THEN
425 ntot = 0
426 CALL reallocate(inp_info%bend_kind, 1, nbends)
427 CALL reallocate(inp_info%bend_a, 1, nbends)
428 CALL reallocate(inp_info%bend_b, 1, nbends)
429 CALL reallocate(inp_info%bend_c, 1, nbends)
430 CALL reallocate(inp_info%bend_k, 1, nbends)
431 CALL reallocate(inp_info%bend_theta0, 1, nbends)
432 CALL reallocate(inp_info%bend_cb, 1, nbends)
433 CALL reallocate(inp_info%bend_r012, 1, nbends)
434 CALL reallocate(inp_info%bend_r032, 1, nbends)
435 CALL reallocate(inp_info%bend_kbs12, 1, nbends)
436 CALL reallocate(inp_info%bend_kbs32, 1, nbends)
437 CALL reallocate(inp_info%bend_kss, 1, nbends)
438 IF (ASSOCIATED(inp_info%bend_legendre)) THEN
439 DO i = 1, SIZE(inp_info%bend_legendre)
440 IF (ASSOCIATED(inp_info%bend_legendre(i)%coeffs)) THEN
441 DEALLOCATE (inp_info%bend_legendre(i)%coeffs)
442 NULLIFY (inp_info%bend_legendre(i)%coeffs)
443 END IF
444 END DO
445 DEALLOCATE (inp_info%bend_legendre)
446 NULLIFY (inp_info%bend_legendre)
447 END IF
448 ALLOCATE (inp_info%bend_legendre(1:nbends))
449 DO i = 1, SIZE(inp_info%bend_legendre(1:nbends))
450 NULLIFY (inp_info%bend_legendre(i)%coeffs)
451 inp_info%bend_legendre(i)%order = 0
452 END DO
453 CALL read_bends_section(inp_info%bend_kind, inp_info%bend_a, inp_info%bend_b, inp_info%bend_c, &
454 inp_info%bend_k, inp_info%bend_theta0, inp_info%bend_cb, &
455 inp_info%bend_r012, inp_info%bend_r032, inp_info%bend_kbs12, &
456 inp_info%bend_kbs32, inp_info%bend_kss, &
457 inp_info%bend_legendre, tmp_section, ntot)
458 END IF
459 tmp_section => section_vals_get_subs_vals(ff_section, "BEND")
460 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nubs)
461 IF (explicit) THEN
462 ntot = 0
463 CALL reallocate(inp_info%ub_kind, 1, nubs)
464 CALL reallocate(inp_info%ub_a, 1, nubs)
465 CALL reallocate(inp_info%ub_b, 1, nubs)
466 CALL reallocate(inp_info%ub_c, 1, nubs)
467 CALL reallocate(inp_info%ub_k, 1, 3, 1, nubs)
468 CALL reallocate(inp_info%ub_r0, 1, nubs)
469 CALL read_ubs_section(inp_info%ub_kind, inp_info%ub_a, inp_info%ub_b, inp_info%ub_c, &
470 inp_info%ub_k, inp_info%ub_r0, tmp_section, ntot)
471 END IF
472 tmp_section => section_vals_get_subs_vals(ff_section, "TORSION")
473 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=ntors)
474 IF (explicit) THEN
475 ntot = 0
476 CALL reallocate(inp_info%torsion_kind, 1, ntors)
477 CALL reallocate(inp_info%torsion_a, 1, ntors)
478 CALL reallocate(inp_info%torsion_b, 1, ntors)
479 CALL reallocate(inp_info%torsion_c, 1, ntors)
480 CALL reallocate(inp_info%torsion_d, 1, ntors)
481 CALL reallocate(inp_info%torsion_k, 1, ntors)
482 CALL reallocate(inp_info%torsion_m, 1, ntors)
483 CALL reallocate(inp_info%torsion_phi0, 1, ntors)
484 CALL read_torsions_section(inp_info%torsion_kind, inp_info%torsion_a, inp_info%torsion_b, &
485 inp_info%torsion_c, inp_info%torsion_d, inp_info%torsion_k, inp_info%torsion_phi0, &
486 inp_info%torsion_m, tmp_section, ntot)
487 END IF
488
489 tmp_section => section_vals_get_subs_vals(ff_section, "IMPROPER")
490 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nimpr)
491 IF (explicit) THEN
492 ntot = 0
493 CALL reallocate(inp_info%impr_kind, 1, nimpr)
494 CALL reallocate(inp_info%impr_a, 1, nimpr)
495 CALL reallocate(inp_info%impr_b, 1, nimpr)
496 CALL reallocate(inp_info%impr_c, 1, nimpr)
497 CALL reallocate(inp_info%impr_d, 1, nimpr)
498 CALL reallocate(inp_info%impr_k, 1, nimpr)
499 CALL reallocate(inp_info%impr_phi0, 1, nimpr)
500 CALL read_improper_section(inp_info%impr_kind, inp_info%impr_a, inp_info%impr_b, &
501 inp_info%impr_c, inp_info%impr_d, inp_info%impr_k, inp_info%impr_phi0, &
502 tmp_section, ntot)
503 END IF
504
505 tmp_section => section_vals_get_subs_vals(ff_section, "OPBEND")
506 CALL section_vals_get(tmp_section, explicit=explicit, n_repetition=nopbend)
507 IF (explicit) THEN
508 ntot = 0
509 CALL reallocate(inp_info%opbend_kind, 1, nopbend)
510 CALL reallocate(inp_info%opbend_a, 1, nopbend)
511 CALL reallocate(inp_info%opbend_b, 1, nopbend)
512 CALL reallocate(inp_info%opbend_c, 1, nopbend)
513 CALL reallocate(inp_info%opbend_d, 1, nopbend)
514 CALL reallocate(inp_info%opbend_k, 1, nopbend)
515 CALL reallocate(inp_info%opbend_phi0, 1, nopbend)
516 CALL read_opbend_section(inp_info%opbend_kind, inp_info%opbend_a, inp_info%opbend_b, &
517 inp_info%opbend_c, inp_info%opbend_d, inp_info%opbend_k, inp_info%opbend_phi0, &
518 tmp_section, ntot)
519 END IF
520
521 END SUBROUTINE read_force_field_section1
522
523! **************************************************************************************************
524!> \brief Set up of the IPBV force fields
525!> \param at1 ...
526!> \param at2 ...
527!> \param ipbv ...
528!> \author teo
529! **************************************************************************************************
530 SUBROUTINE set_ipbv_ff(at1, at2, ipbv)
531 CHARACTER(LEN=*), INTENT(IN) :: at1, at2
532 TYPE(ipbv_pot_type), POINTER :: ipbv
533
534 IF ((at1(1:1) == 'O') .AND. (at2(1:1) == 'O')) THEN
535 ipbv%rcore = 0.9_dp ! a.u.
536 ipbv%m = -1.2226442563398141e+11_dp ! Kelvin/a.u.
537 ipbv%b = 1.1791292385486696e+11_dp ! Hartree
538
539 ! Hartree*a.u.^2
540 ipbv%a(2) = 4.786380682394_dp
541 ipbv%a(3) = -1543.407053545_dp
542 ipbv%a(4) = 88783.31188529_dp
543 ipbv%a(5) = -2361200.155376_dp
544 ipbv%a(6) = 35940504.84679_dp
545 ipbv%a(7) = -339762743.6358_dp
546 ipbv%a(8) = 2043874926.466_dp
547 ipbv%a(9) = -7654856796.383_dp
548 ipbv%a(10) = 16195251405.65_dp
549 ipbv%a(11) = -13140392992.18_dp
550 ipbv%a(12) = -9285572894.245_dp
551 ipbv%a(13) = 8756947519.029_dp
552 ipbv%a(14) = 15793297761.67_dp
553 ipbv%a(15) = 12917180227.21_dp
554 ELSEIF (((at1(1:1) == 'O') .AND. (at2(1:1) == 'H')) .OR. &
555 ((at1(1:1) == 'H') .AND. (at2(1:1) == 'O'))) THEN
556 ipbv%rcore = 2.95_dp ! a.u.
557
558 ipbv%m = -0.004025691139759147_dp ! Hartree/a.u.
559 ipbv%b = -2.193731138097428_dp ! Hartree
560 ! Hartree*a.u.^2
561 ipbv%a(2) = -195.7716013277_dp
562 ipbv%a(3) = 15343.78613395_dp
563 ipbv%a(4) = -530864.4586516_dp
564 ipbv%a(5) = 10707934.39058_dp
565 ipbv%a(6) = -140099704.7890_dp
566 ipbv%a(7) = 1250943273.785_dp
567 ipbv%a(8) = -7795458330.676_dp
568 ipbv%a(9) = 33955897217.31_dp
569 ipbv%a(10) = -101135640744.0_dp
570 ipbv%a(11) = 193107995718.7_dp
571 ipbv%a(12) = -193440560940.0_dp
572 ipbv%a(13) = -4224406093.918e0_dp
573 ipbv%a(14) = 217192386506.5e0_dp
574 ipbv%a(15) = -157581228915.5_dp
575 ELSEIF ((at1(1:1) == 'H') .AND. (at2(1:1) == 'H')) THEN
576 ipbv%rcore = 3.165_dp ! a.u.
577 ipbv%m = 0.002639704108787555_dp ! Hartree/a.u.
578 ipbv%b = -0.2735482611857583_dp ! Hartree
579 ! Hartree*a.u.^2
580 ipbv%a(2) = -26.29456010782_dp
581 ipbv%a(3) = 2373.352548248_dp
582 ipbv%a(4) = -93880.43551360_dp
583 ipbv%a(5) = 2154624.884809_dp
584 ipbv%a(6) = -31965151.34955_dp
585 ipbv%a(7) = 322781785.3278_dp
586 ipbv%a(8) = -2271097368.668_dp
587 ipbv%a(9) = 11169163192.90_dp
588 ipbv%a(10) = -37684457778.47_dp
589 ipbv%a(11) = 82562104256.03_dp
590 ipbv%a(12) = -100510435213.4_dp
591 ipbv%a(13) = 24570342714.65e0_dp
592 ipbv%a(14) = 88766181532.94e0_dp
593 ipbv%a(15) = -79705131323.98_dp
594 ELSE
595 cpabort("IPBV only for WATER")
596 END IF
597 END SUBROUTINE set_ipbv_ff
598
599! **************************************************************************************************
600!> \brief Set up of the BMHFT force fields
601!> \param at1 ...
602!> \param at2 ...
603!> \param ft ...
604!> \author teo
605! **************************************************************************************************
606 SUBROUTINE set_bmhft_ff(at1, at2, ft)
607 CHARACTER(LEN=*), INTENT(IN) :: at1, at2
608 TYPE(ft_pot_type), POINTER :: ft
609
610 ft%b = cp_unit_to_cp2k(3.1545_dp, "angstrom^-1")
611 IF ((at1(1:2) == 'NA') .AND. (at2(1:2) == 'NA')) THEN
612 ft%a = cp_unit_to_cp2k(424.097_dp, "eV")
613 ft%c = cp_unit_to_cp2k(1.05_dp, "eV*angstrom^6")
614 ft%d = cp_unit_to_cp2k(0.499_dp, "eV*angstrom^8")
615 ELSEIF (((at1(1:2) == 'NA') .AND. (at2(1:2) == 'CL')) .OR. &
616 ((at1(1:2) == 'CL') .AND. (at2(1:2) == 'NA'))) THEN
617 ft%a = cp_unit_to_cp2k(1256.31_dp, "eV")
618 ft%c = cp_unit_to_cp2k(7.00_dp, "eV*angstrom^6")
619 ft%d = cp_unit_to_cp2k(8.676_dp, "eV*angstrom^8")
620 ELSEIF ((at1(1:2) == 'CL') .AND. (at2(1:2) == 'CL')) THEN
621 ft%a = cp_unit_to_cp2k(3488.998_dp, "eV")
622 ft%c = cp_unit_to_cp2k(72.50_dp, "eV*angstrom^6")
623 ft%d = cp_unit_to_cp2k(145.427_dp, "eV*angstrom^8")
624 ELSE
625 cpabort("BMHFT only for NaCl")
626 END IF
627
628 END SUBROUTINE set_bmhft_ff
629
630! **************************************************************************************************
631!> \brief Set up of the BMHFTD force fields
632!> \author Mathieu Salanne 05.2010
633! **************************************************************************************************
634 SUBROUTINE set_bmhftd_ff()
635
636 cpabort("No default parameters present for BMHFTD")
637
638 END SUBROUTINE set_bmhftd_ff
639
640! **************************************************************************************************
641!> \brief Reads the EAM section
642!> \param nonbonded ...
643!> \param section ...
644!> \param start ...
645!> \param para_env ...
646!> \param mm_section ...
647!> \author teo
648! **************************************************************************************************
649 SUBROUTINE read_eam_section(nonbonded, section, start, para_env, mm_section)
650 TYPE(pair_potential_p_type), POINTER :: nonbonded
651 TYPE(section_vals_type), POINTER :: section
652 INTEGER, INTENT(IN) :: start
653 TYPE(mp_para_env_type), POINTER :: para_env
654 TYPE(section_vals_type), POINTER :: mm_section
655
656 CHARACTER(LEN=default_string_length), &
657 DIMENSION(:), POINTER :: atm_names
658 INTEGER :: isec, n_items
659
660 CALL section_vals_get(section, n_repetition=n_items)
661 DO isec = 1, n_items
662 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
663
664 nonbonded%pot(start + isec)%pot%type = ea_type
665 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
666 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
667 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
668 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
669 CALL section_vals_val_get(section, "PARM_FILE_NAME", i_rep_section=isec, &
670 c_val=nonbonded%pot(start + isec)%pot%set(1)%eam%eam_file_name)
671 CALL read_eam_data(nonbonded%pot(start + isec)%pot%set(1)%eam, para_env, mm_section)
672 nonbonded%pot(start + isec)%pot%rcutsq = nonbonded%pot(start + isec)%pot%set(1)%eam%acutal**2
673 END DO
674 END SUBROUTINE read_eam_section
675
676! **************************************************************************************
677!> \brief Reads the ACE section
678!> \param nonbonded ...
679!> \param section ...
680!> \param start ...
681! **************************************************************************************************
682 SUBROUTINE read_ace_section(nonbonded, section, start)
683 TYPE(pair_potential_p_type), POINTER :: nonbonded
684 TYPE(section_vals_type), POINTER :: section
685 INTEGER, INTENT(IN) :: start
686
687 CHARACTER(LEN=2), ALLOCATABLE, DIMENSION(:) :: ace_atype_symbol
688 CHARACTER(LEN=default_path_length) :: ace_filename
689 CHARACTER(LEN=default_string_length) :: ace_file_name
690 CHARACTER(LEN=default_string_length), &
691 DIMENSION(:), POINTER :: atm_names
692 INTEGER :: ace_ntype, isec, jsec, n_items
693 REAL(KIND=dp), ALLOCATABLE, DIMENSION(:, :) :: rcutall
694 TYPE(ace_model_type) :: model
695
696 n_items = 1
697 isec = 1
698 n_items = isec*n_items
699 CALL section_vals_val_get(section, "ATOMS", c_vals=atm_names)
700
701 ace_ntype = SIZE(atm_names)
702 ALLOCATE (ace_atype_symbol(ace_ntype), rcutall(ace_ntype, ace_ntype))
703 DO isec = 1, ace_ntype
704 ace_atype_symbol(isec) = atm_names(isec) (1:2)
705 END DO
706 CALL section_vals_val_get(section, "POT_FILE_NAME", c_val=ace_file_name)
707
708 ace_filename = discover_file(ace_file_name)
709
710#if defined(__ACE)
711 ! need ace_model_initialize() here somewhere to get rcut
712 CALL ace_model_initialize(ntypec=ace_ntype, symbolc=ace_atype_symbol, &
713 fname=trim(ace_filename), rcutc=rcutall, model=model)
714#else
715 cpabort("CP2K was compiled without ACE library.")
716#endif
717
718 DO isec = 1, SIZE(atm_names)
719 DO jsec = isec, SIZE(atm_names)
720 nonbonded%pot(start + n_items)%pot%type = ace_type
721 nonbonded%pot(start + n_items)%pot%at1 = atm_names(isec)
722 nonbonded%pot(start + n_items)%pot%at2 = atm_names(jsec)
723 CALL uppercase(nonbonded%pot(start + n_items)%pot%at1)
724 CALL uppercase(nonbonded%pot(start + n_items)%pot%at2)
725
726 nonbonded%pot(start + n_items)%pot%set(1)%ace%ace_file_name = ace_filename
727 nonbonded%pot(start + n_items)%pot%set(1)%ace%atom_ace_type = isec
728 nonbonded%pot(start + n_items)%pot%set(1)%ace%model = model
729
730 !using rcutall(isec,jsec) instead of maxval(rcutall) TODO check that
731 !it shouldn't be jsec,isec?
732 nonbonded%pot(start + n_items)%pot%rcutsq = cp_unit_to_cp2k(rcutall(isec, jsec), "angstrom")**2
733
734 n_items = n_items + 1
735 END DO
736 END DO
737 END SUBROUTINE read_ace_section
738
739! **************************************************************************************
740!> \brief Reads the DEEPMD section
741!> \param nonbonded ...
742!> \param section ...
743!> \param start ...
744!> \author teo
745! **************************************************************************************************
746 SUBROUTINE read_deepmd_section(nonbonded, section, start)
747 TYPE(pair_potential_p_type), POINTER :: nonbonded
748 TYPE(section_vals_type), POINTER :: section
749 INTEGER, INTENT(IN) :: start
750
751 CHARACTER(LEN=default_string_length) :: deepmd_file_name
752 CHARACTER(LEN=default_string_length), &
753 DIMENSION(:), POINTER :: atm_names
754 INTEGER :: isec, jsec, n_items
755 INTEGER, DIMENSION(:), POINTER :: atm_deepmd_types
756
757 n_items = 1
758 isec = 1
759 n_items = isec*n_items
760 CALL section_vals_val_get(section, "ATOMS", c_vals=atm_names)
761 CALL section_vals_val_get(section, "ATOMS_DEEPMD_TYPE", i_vals=atm_deepmd_types)
762 CALL section_vals_val_get(section, "POT_FILE_NAME", c_val=deepmd_file_name)
763
764 DO isec = 1, SIZE(atm_names)
765 DO jsec = isec, SIZE(atm_names)
766 nonbonded%pot(start + n_items)%pot%type = deepmd_type
767 nonbonded%pot(start + n_items)%pot%at1 = atm_names(isec)
768 nonbonded%pot(start + n_items)%pot%at2 = atm_names(jsec)
769 CALL uppercase(nonbonded%pot(start + n_items)%pot%at1)
770 CALL uppercase(nonbonded%pot(start + n_items)%pot%at2)
771
772 nonbonded%pot(start + n_items)%pot%set(1)%deepmd%deepmd_file_name = discover_file(deepmd_file_name)
773 nonbonded%pot(start + n_items)%pot%set(1)%deepmd%atom_deepmd_type = atm_deepmd_types(isec)
774 nonbonded%pot(start + n_items)%pot%rcutsq = 0.0_dp
775 n_items = n_items + 1
776 END DO
777 END DO
778 END SUBROUTINE read_deepmd_section
779
780! **************************************************************************************************
781!> \brief Reads the NEQUIP section
782!> \param nonbonded ...
783!> \param section ...
784!> \param start ...
785!> \author Gabriele Tocci
786! **************************************************************************************************
787 SUBROUTINE read_nequip_section(nonbonded, section, start)
788 TYPE(pair_potential_p_type), POINTER :: nonbonded
789 TYPE(section_vals_type), POINTER :: section
790 INTEGER, INTENT(IN) :: start
791
792 CHARACTER(LEN=default_string_length) :: nequip_file_name, unit_cell, &
793 unit_coords, unit_energy, unit_forces
794 CHARACTER(LEN=default_string_length), &
795 DIMENSION(:), POINTER :: atm_names
796 INTEGER :: isec, jsec, n_items
797 TYPE(nequip_pot_type) :: nequip
798
799 n_items = 1
800 isec = 1
801 n_items = isec*n_items
802 CALL section_vals_val_get(section, "ATOMS", c_vals=atm_names)
803 CALL section_vals_val_get(section, "PARM_FILE_NAME", c_val=nequip_file_name)
804 CALL section_vals_val_get(section, "UNIT_COORDS", c_val=unit_coords)
805 CALL section_vals_val_get(section, "UNIT_ENERGY", c_val=unit_energy)
806 CALL section_vals_val_get(section, "UNIT_FORCES", c_val=unit_forces)
807 CALL section_vals_val_get(section, "UNIT_CELL", c_val=unit_cell)
808
809 ! Since reading nequip metadata is expensive we do it only once.
810 nequip%nequip_file_name = discover_file(nequip_file_name)
811 nequip%unit_coords = unit_coords
812 nequip%unit_forces = unit_forces
813 nequip%unit_energy = unit_energy
814 nequip%unit_cell = unit_cell
815 CALL read_nequip_data(nequip)
816 CALL check_cp2k_atom_names_in_torch(atm_names, nequip%type_names_torch)
817
818 DO isec = 1, SIZE(atm_names)
819 DO jsec = isec, SIZE(atm_names)
820 nonbonded%pot(start + n_items)%pot%type = nequip_type
821 nonbonded%pot(start + n_items)%pot%at1 = atm_names(isec)
822 nonbonded%pot(start + n_items)%pot%at2 = atm_names(jsec)
823 CALL uppercase(nonbonded%pot(start + n_items)%pot%at1)
824 CALL uppercase(nonbonded%pot(start + n_items)%pot%at2)
825 nonbonded%pot(start + n_items)%pot%set(1)%nequip = nequip
826 nonbonded%pot(start + n_items)%pot%rcutsq = nequip%rcutsq
827 n_items = n_items + 1
828 END DO
829 END DO
830
831 END SUBROUTINE read_nequip_section
832
833! **************************************************************************************************
834!> \brief Reads the ALLEGRO section
835!> \param nonbonded ...
836!> \param section ...
837!> \param start ...
838!> \author Gabriele Tocci
839! **************************************************************************************************
840 SUBROUTINE read_allegro_section(nonbonded, section, start)
841 TYPE(pair_potential_p_type), POINTER :: nonbonded
842 TYPE(section_vals_type), POINTER :: section
843 INTEGER, INTENT(IN) :: start
844
845 CHARACTER(LEN=default_string_length) :: allegro_file_name, unit_cell, &
846 unit_coords, unit_energy, unit_forces
847 CHARACTER(LEN=default_string_length), &
848 DIMENSION(:), POINTER :: atm_names
849 INTEGER :: isec, jsec, n_items
850 TYPE(allegro_pot_type) :: allegro
851
852 n_items = 1
853 isec = 1
854 n_items = isec*n_items
855 CALL section_vals_val_get(section, "ATOMS", c_vals=atm_names)
856 CALL section_vals_val_get(section, "PARM_FILE_NAME", c_val=allegro_file_name)
857 CALL section_vals_val_get(section, "UNIT_COORDS", c_val=unit_coords)
858 CALL section_vals_val_get(section, "UNIT_ENERGY", c_val=unit_energy)
859 CALL section_vals_val_get(section, "UNIT_FORCES", c_val=unit_forces)
860 CALL section_vals_val_get(section, "UNIT_CELL", c_val=unit_cell)
861
862 ! Since reading allegro metadata is expensive we do it only once.
863 allegro%allegro_file_name = discover_file(allegro_file_name)
864 allegro%unit_coords = unit_coords
865 allegro%unit_forces = unit_forces
866 allegro%unit_energy = unit_energy
867 allegro%unit_cell = unit_cell
868 CALL read_allegro_data(allegro)
869 CALL check_cp2k_atom_names_in_torch(atm_names, allegro%type_names_torch)
870
871 DO isec = 1, SIZE(atm_names)
872 DO jsec = isec, SIZE(atm_names)
873 nonbonded%pot(start + n_items)%pot%type = allegro_type
874 nonbonded%pot(start + n_items)%pot%at1 = atm_names(isec)
875 nonbonded%pot(start + n_items)%pot%at2 = atm_names(jsec)
876 CALL uppercase(nonbonded%pot(start + n_items)%pot%at1)
877 CALL uppercase(nonbonded%pot(start + n_items)%pot%at2)
878 nonbonded%pot(start + n_items)%pot%set(1)%allegro = allegro
879 nonbonded%pot(start + n_items)%pot%rcutsq = allegro%rcutsq
880 n_items = n_items + 1
881 END DO
882 END DO
883 END SUBROUTINE read_allegro_section
884
885! **************************************************************************************************
886!> \brief Reads the LJ section
887!> \param nonbonded ...
888!> \param section ...
889!> \param start ...
890!> \author teo
891! **************************************************************************************************
892 SUBROUTINE read_lj_section(nonbonded, section, start)
893 TYPE(pair_potential_p_type), POINTER :: nonbonded
894 TYPE(section_vals_type), POINTER :: section
895 INTEGER, INTENT(IN) :: start
896
897 CHARACTER(LEN=default_string_length), &
898 DIMENSION(:), POINTER :: atm_names
899 INTEGER :: isec, n_items, n_rep
900 REAL(kind=dp) :: epsilon, rcut, sigma
901
902 CALL section_vals_get(section, n_repetition=n_items)
903 DO isec = 1, n_items
904 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
905 CALL section_vals_val_get(section, "EPSILON", i_rep_section=isec, r_val=epsilon)
906 CALL section_vals_val_get(section, "SIGMA", i_rep_section=isec, r_val=sigma)
907 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
908
909 nonbonded%pot(start + isec)%pot%type = lj_charmm_type
910 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
911 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
912 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
913 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
914 nonbonded%pot(start + isec)%pot%set(1)%lj%epsilon = epsilon
915 nonbonded%pot(start + isec)%pot%set(1)%lj%sigma6 = sigma**6
916 nonbonded%pot(start + isec)%pot%set(1)%lj%sigma12 = sigma**12
917 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
918 !
919 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
920 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
921 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
922 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
923 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
924 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
925 END DO
926 END SUBROUTINE read_lj_section
927
928! **************************************************************************************************
929!> \brief Reads the WILLIAMS section
930!> \param nonbonded ...
931!> \param section ...
932!> \param start ...
933!> \author teo
934! **************************************************************************************************
935 SUBROUTINE read_wl_section(nonbonded, section, start)
936 TYPE(pair_potential_p_type), POINTER :: nonbonded
937 TYPE(section_vals_type), POINTER :: section
938 INTEGER, INTENT(IN) :: start
939
940 CHARACTER(LEN=default_string_length), &
941 DIMENSION(:), POINTER :: atm_names
942 INTEGER :: isec, n_items, n_rep
943 REAL(kind=dp) :: a, b, c, rcut
944
945 CALL section_vals_get(section, n_repetition=n_items)
946 DO isec = 1, n_items
947 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
948 CALL section_vals_val_get(section, "A", i_rep_section=isec, r_val=a)
949 CALL section_vals_val_get(section, "B", i_rep_section=isec, r_val=b)
950 CALL section_vals_val_get(section, "C", i_rep_section=isec, r_val=c)
951 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
952
953 nonbonded%pot(start + isec)%pot%type = wl_type
954 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
955 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
956 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
957 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
958 nonbonded%pot(start + isec)%pot%set(1)%willis%a = a
959 nonbonded%pot(start + isec)%pot%set(1)%willis%b = b
960 nonbonded%pot(start + isec)%pot%set(1)%willis%c = c
961 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
962 !
963 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
964 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
965 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
966 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
967 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
968 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
969 END DO
970 END SUBROUTINE read_wl_section
971
972! **************************************************************************************************
973!> \brief Reads the GOODWIN section
974!> \param nonbonded ...
975!> \param section ...
976!> \param start ...
977!> \author teo
978! **************************************************************************************************
979 SUBROUTINE read_gd_section(nonbonded, section, start)
980 TYPE(pair_potential_p_type), POINTER :: nonbonded
981 TYPE(section_vals_type), POINTER :: section
982 INTEGER, INTENT(IN) :: start
983
984 CHARACTER(LEN=default_string_length), &
985 DIMENSION(:), POINTER :: atm_names
986 INTEGER :: isec, m, mc, n_items, n_rep
987 REAL(kind=dp) :: d, dc, rcut, vr0
988
989 CALL section_vals_get(section, n_repetition=n_items)
990 DO isec = 1, n_items
991 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
992 CALL section_vals_val_get(section, "VR0", i_rep_section=isec, r_val=vr0)
993 CALL section_vals_val_get(section, "D", i_rep_section=isec, r_val=d)
994 CALL section_vals_val_get(section, "DC", i_rep_section=isec, r_val=dc)
995 CALL section_vals_val_get(section, "M", i_rep_section=isec, i_val=m)
996 CALL section_vals_val_get(section, "MC", i_rep_section=isec, i_val=mc)
997 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
998
999 nonbonded%pot(start + isec)%pot%type = gw_type
1000 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1001 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1002 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1003 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1004 nonbonded%pot(start + isec)%pot%set(1)%goodwin%vr0 = vr0
1005 nonbonded%pot(start + isec)%pot%set(1)%goodwin%d = d
1006 nonbonded%pot(start + isec)%pot%set(1)%goodwin%dc = dc
1007 nonbonded%pot(start + isec)%pot%set(1)%goodwin%m = m
1008 nonbonded%pot(start + isec)%pot%set(1)%goodwin%mc = mc
1009 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
1010 !
1011 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1012 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1013 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1014 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1015 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1016 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1017 END DO
1018 END SUBROUTINE read_gd_section
1019
1020! **************************************************************************************************
1021!> \brief Reads the IPBV section
1022!> \param nonbonded ...
1023!> \param section ...
1024!> \param start ...
1025!> \author teo
1026! **************************************************************************************************
1027 SUBROUTINE read_ipbv_section(nonbonded, section, start)
1028 TYPE(pair_potential_p_type), POINTER :: nonbonded
1029 TYPE(section_vals_type), POINTER :: section
1030 INTEGER, INTENT(IN) :: start
1031
1032 CHARACTER(LEN=default_string_length), &
1033 DIMENSION(:), POINTER :: atm_names
1034 INTEGER :: isec, n_items, n_rep
1035 REAL(kind=dp) :: rcut
1036
1037 CALL section_vals_get(section, n_repetition=n_items)
1038 DO isec = 1, n_items
1039 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1040 nonbonded%pot(start + isec)%pot%type = ip_type
1041 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1042 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1043 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1044 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1045 CALL set_ipbv_ff(nonbonded%pot(start + isec)%pot%at1, nonbonded%pot(start + isec)%pot%at2, &
1046 nonbonded%pot(start + isec)%pot%set(1)%ipbv)
1047 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1048 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1049 !
1050 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1051 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1052 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1053 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1054 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1055 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1056 END DO
1057 END SUBROUTINE read_ipbv_section
1058
1059! **************************************************************************************************
1060!> \brief Reads the BMHFT section
1061!> \param nonbonded ...
1062!> \param section ...
1063!> \param start ...
1064!> \author teo
1065! **************************************************************************************************
1066 SUBROUTINE read_bmhft_section(nonbonded, section, start)
1067 TYPE(pair_potential_p_type), POINTER :: nonbonded
1068 TYPE(section_vals_type), POINTER :: section
1069 INTEGER, INTENT(IN) :: start
1070
1071 CHARACTER(LEN=default_string_length), DIMENSION(2) :: map_atoms
1072 CHARACTER(LEN=default_string_length), &
1073 DIMENSION(:), POINTER :: atm_names
1074 INTEGER :: i, isec, n_items, n_rep
1075 REAL(kind=dp) :: rcut
1076
1077 CALL section_vals_get(section, n_repetition=n_items)
1078 DO isec = 1, n_items
1079 CALL cite_reference(tosi1964a)
1080 CALL cite_reference(tosi1964b)
1081 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1082 nonbonded%pot(start + isec)%pot%type = ft_type
1083 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1084 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1085 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1086 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1087
1088 CALL section_vals_val_get(section, "A", i_rep_section=isec, n_rep_val=i)
1089 IF (i == 1) THEN
1090 CALL section_vals_val_get(section, "A", i_rep_section=isec, &
1091 r_val=nonbonded%pot(start + isec)%pot%set(1)%ft%a)
1092 CALL section_vals_val_get(section, "B", i_rep_section=isec, &
1093 r_val=nonbonded%pot(start + isec)%pot%set(1)%ft%b)
1094 CALL section_vals_val_get(section, "C", i_rep_section=isec, &
1095 r_val=nonbonded%pot(start + isec)%pot%set(1)%ft%c)
1096 CALL section_vals_val_get(section, "D", i_rep_section=isec, &
1097 r_val=nonbonded%pot(start + isec)%pot%set(1)%ft%d)
1098 ELSE
1099 CALL section_vals_val_get(section, "MAP_ATOMS", i_rep_section=isec, c_vals=atm_names)
1100 map_atoms = atm_names
1101 CALL uppercase(map_atoms(1))
1102 CALL uppercase(map_atoms(2))
1103 CALL set_bmhft_ff(map_atoms(1), map_atoms(2), nonbonded%pot(start + isec)%pot%set(1)%ft)
1104 END IF
1105 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1106 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1107 !
1108 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1109 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1110 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1111 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1112 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1113 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1114 END DO
1115 END SUBROUTINE read_bmhft_section
1116
1117! **************************************************************************************************
1118!> \brief Reads the BMHFTD section
1119!> \param nonbonded ...
1120!> \param section ...
1121!> \param start ...
1122!> \author Mathieu Salanne 05.2010
1123! **************************************************************************************************
1124 SUBROUTINE read_bmhftd_section(nonbonded, section, start)
1125 TYPE(pair_potential_p_type), POINTER :: nonbonded
1126 TYPE(section_vals_type), POINTER :: section
1127 INTEGER, INTENT(IN) :: start
1128
1129 CHARACTER(LEN=default_string_length), DIMENSION(2) :: map_atoms
1130 CHARACTER(LEN=default_string_length), &
1131 DIMENSION(:), POINTER :: atm_names
1132 INTEGER :: i, isec, n_items, n_rep
1133 REAL(kind=dp) :: rcut
1134 REAL(kind=dp), DIMENSION(:), POINTER :: bd_vals
1135
1136 NULLIFY (bd_vals)
1137
1138 CALL section_vals_get(section, n_repetition=n_items)
1139
1140 DO isec = 1, n_items
1141 CALL cite_reference(tosi1964a)
1142 CALL cite_reference(tosi1964b)
1143 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1144 nonbonded%pot(start + isec)%pot%type = ftd_type
1145 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1146 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1147 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1148 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1149
1150 CALL section_vals_val_get(section, "A", i_rep_section=isec, n_rep_val=i)
1151 IF (i == 1) THEN
1152 CALL section_vals_val_get(section, "A", i_rep_section=isec, &
1153 r_val=nonbonded%pot(start + isec)%pot%set(1)%ftd%a)
1154 CALL section_vals_val_get(section, "B", i_rep_section=isec, &
1155 r_val=nonbonded%pot(start + isec)%pot%set(1)%ftd%b)
1156 CALL section_vals_val_get(section, "C", i_rep_section=isec, &
1157 r_val=nonbonded%pot(start + isec)%pot%set(1)%ftd%c)
1158 CALL section_vals_val_get(section, "D", i_rep_section=isec, &
1159 r_val=nonbonded%pot(start + isec)%pot%set(1)%ftd%d)
1160 CALL section_vals_val_get(section, "BD", i_rep_section=isec, r_vals=bd_vals)
1161 IF (ASSOCIATED(bd_vals)) THEN
1162 SELECT CASE (SIZE(bd_vals))
1163 CASE (0)
1164 cpabort("No values specified for parameter BD in section &BMHFTD")
1165 CASE (1)
1166 nonbonded%pot(start + isec)%pot%set(1)%ftd%bd(1:2) = bd_vals(1)
1167 CASE (2)
1168 nonbonded%pot(start + isec)%pot%set(1)%ftd%bd(1:2) = bd_vals(1:2)
1169 CASE DEFAULT
1170 cpabort("Too many values specified for parameter BD in section &BMHFTD")
1171 END SELECT
1172 ELSE
1173 cpabort("Parameter BD in section &BMHFTD was not specified")
1174 END IF
1175 ELSE
1176 CALL section_vals_val_get(section, "MAP_ATOMS", i_rep_section=isec, c_vals=atm_names)
1177 map_atoms = atm_names
1178 CALL uppercase(map_atoms(1))
1179 CALL uppercase(map_atoms(2))
1180 CALL set_bmhftd_ff()
1181 END IF
1182 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1183 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1184 !
1185 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1186 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1187 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1188 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1189 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1190 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1191 END DO
1192 END SUBROUTINE read_bmhftd_section
1193
1194! **************************************************************************************************
1195!> \brief Reads the Buckingham 4 Ranges potential section
1196!> \param nonbonded ...
1197!> \param section ...
1198!> \param start ...
1199!> \par History
1200!> MK (11.11.2010): Automatic fit of the (default) polynomial coefficients
1201!> \author MI,MK
1202! **************************************************************************************************
1203 SUBROUTINE read_b4_section(nonbonded, section, start)
1204
1205 TYPE(pair_potential_p_type), POINTER :: nonbonded
1206 TYPE(section_vals_type), POINTER :: section
1207 INTEGER, INTENT(IN) :: start
1208
1209 CHARACTER(LEN=default_string_length), &
1210 DIMENSION(:), POINTER :: atm_names
1211 INTEGER :: i, ir, isec, n_items, n_rep, np1, np2
1212 LOGICAL :: explicit_poly1, explicit_poly2
1213 REAL(kind=dp) :: a, b, c, eval_error, r1, r2, r3, rcut
1214 REAL(kind=dp), DIMENSION(10) :: v, x
1215 REAL(kind=dp), DIMENSION(10, 10) :: p, p_inv
1216 REAL(kind=dp), DIMENSION(:), POINTER :: coeff1, coeff2, list
1217
1218 NULLIFY (coeff1)
1219 NULLIFY (coeff2)
1220
1221 CALL section_vals_get(section, n_repetition=n_items)
1222
1223 DO isec = 1, n_items
1224 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1225 CALL section_vals_val_get(section, "A", i_rep_section=isec, r_val=a)
1226 CALL section_vals_val_get(section, "B", i_rep_section=isec, r_val=b)
1227 CALL section_vals_val_get(section, "C", i_rep_section=isec, r_val=c)
1228 CALL section_vals_val_get(section, "R1", i_rep_section=isec, r_val=r1)
1229 CALL section_vals_val_get(section, "R2", i_rep_section=isec, r_val=r2)
1230 CALL section_vals_val_get(section, "R3", i_rep_section=isec, r_val=r3)
1231 CALL section_vals_val_get(section, "POLY1", explicit=explicit_poly1, n_rep_val=n_rep)
1232 ! Check if polynomial coefficients were specified for range 2 and 3 explicitly
1233 IF (explicit_poly1) THEN
1234 np1 = 0
1235 DO ir = 1, n_rep
1236 NULLIFY (list)
1237 CALL section_vals_val_get(section, "POLY1", i_rep_val=ir, r_vals=list)
1238 IF (ASSOCIATED(list)) THEN
1239 CALL reallocate(coeff1, 0, np1 + SIZE(list) - 1)
1240 DO i = 1, SIZE(list)
1241 coeff1(i + np1 - 1) = list(i)
1242 END DO
1243 np1 = np1 + SIZE(list)
1244 END IF
1245 END DO
1246 END IF
1247 CALL section_vals_val_get(section, "POLY2", explicit=explicit_poly2, n_rep_val=n_rep)
1248 IF (explicit_poly2) THEN
1249 np2 = 0
1250 DO ir = 1, n_rep
1251 NULLIFY (list)
1252 CALL section_vals_val_get(section, "POLY2", i_rep_val=ir, r_vals=list)
1253 IF (ASSOCIATED(list)) THEN
1254 CALL reallocate(coeff2, 0, np2 + SIZE(list) - 1)
1255 DO i = 1, SIZE(list)
1256 coeff2(i + np2 - 1) = list(i)
1257 END DO
1258 np2 = np2 + SIZE(list)
1259 END IF
1260 END DO
1261 END IF
1262 ! Default is a 5th/3rd-order polynomial fit
1263 IF ((.NOT. explicit_poly1) .OR. (.NOT. explicit_poly2)) THEN
1264 ! Build matrix p and vector v to calculate the polynomial coefficients
1265 ! in the vector x from p*x = v
1266 p(:, :) = 0.0_dp
1267 ! Row 1: Match the 5th-order polynomial and the potential at r1
1268 p(1, 1) = 1.0_dp
1269 DO i = 2, 6
1270 p(1, i) = p(1, i - 1)*r1
1271 END DO
1272 ! Row 2: Match the first derivatives of the 5th-order polynomial and the potential at r1
1273 DO i = 2, 6
1274 p(2, i) = real(i - 1, kind=dp)*p(1, i - 1)
1275 END DO
1276 ! Row 3: Match the second derivatives of the 5th-order polynomial and the potential at r1
1277 DO i = 3, 6
1278 p(3, i) = real(i - 1, kind=dp)*p(2, i - 1)
1279 END DO
1280 ! Row 4: Match the 5th-order and the 3rd-order polynomials at r2
1281 p(4, 1) = 1.0_dp
1282 DO i = 2, 6
1283 p(4, i) = p(4, i - 1)*r2
1284 END DO
1285 p(4, 7) = -1.0_dp
1286 DO i = 8, 10
1287 p(4, i) = p(4, i - 1)*r2
1288 END DO
1289 ! Row 5: Match the first derivatives of the 5th-order and the 3rd-order polynomials at r2
1290 DO i = 2, 6
1291 p(5, i) = real(i - 1, kind=dp)*p(4, i - 1)
1292 END DO
1293 DO i = 8, 10
1294 p(5, i) = real(i - 7, kind=dp)*p(4, i - 1)
1295 END DO
1296 ! Row 6: Match the second derivatives of the 5th-order and the 3rd-order polynomials at r2
1297 DO i = 3, 6
1298 p(6, i) = real(i - 1, kind=dp)*p(5, i - 1)
1299 END DO
1300 DO i = 9, 10
1301 p(6, i) = real(i - 7, kind=dp)*p(5, i - 1)
1302 END DO
1303 ! Row 7: Minimum at r2, ie. the first derivative of the 3rd-order polynomial has to be zero at r2
1304 DO i = 8, 10
1305 p(7, i) = -p(5, i)
1306 END DO
1307 ! Row 8: Match the 3rd-order polynomial and the potential at r3
1308 p(8, 7) = 1.0_dp
1309 DO i = 8, 10
1310 p(8, i) = p(8, i - 1)*r3
1311 END DO
1312 ! Row 9: Match the first derivatives of the 3rd-order polynomial and the potential at r3
1313 DO i = 8, 10
1314 p(9, i) = real(i - 7, kind=dp)*p(8, i - 1)
1315 END DO
1316 ! Row 10: Match the second derivatives of the 3rd-order polynomial and the potential at r3
1317 DO i = 9, 10
1318 p(10, i) = real(i - 7, kind=dp)*p(9, i - 1)
1319 END DO
1320 ! Build the vector v
1321 v(1) = a*exp(-b*r1)
1322 v(2) = -b*v(1)
1323 v(3) = -b*v(2)
1324 v(4:7) = 0.0_dp
1325 v(8) = -c/p(8, 10)**2 ! = -c/r3**6
1326 v(9) = -6.0_dp*v(8)/r3
1327 v(10) = -7.0_dp*v(9)/r3
1328 ! Calculate p_inv the inverse of the matrix p
1329 p_inv(:, :) = 0.0_dp
1330 CALL invert_matrix(p, p_inv, eval_error)
1331 IF (eval_error >= 1.0e-8_dp) &
1332 CALL cp_warn(__location__, &
1333 "The polynomial fit for the BUCK4RANGES potential is only accurate to "// &
1334 trim(cp_to_string(eval_error)))
1335 ! Get the 6 coefficients of the 5th-order polynomial -> x(1:6)
1336 ! and the 4 coefficients of the 3rd-order polynomial -> x(7:10)
1337 x(:) = matmul(p_inv(:, :), v(:))
1338 ELSE
1339 x(:) = 0.0_dp
1340 END IF
1341
1342 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1343
1344 nonbonded%pot(start + isec)%pot%type = b4_type
1345 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1346 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1347 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1348 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1349 nonbonded%pot(start + isec)%pot%set(1)%buck4r%a = a
1350 nonbonded%pot(start + isec)%pot%set(1)%buck4r%b = b
1351 nonbonded%pot(start + isec)%pot%set(1)%buck4r%c = c
1352 nonbonded%pot(start + isec)%pot%set(1)%buck4r%r1 = r1
1353 nonbonded%pot(start + isec)%pot%set(1)%buck4r%r2 = r2
1354 nonbonded%pot(start + isec)%pot%set(1)%buck4r%r3 = r3
1355 IF ((.NOT. explicit_poly1) .OR. (.NOT. explicit_poly2)) THEN
1356 nonbonded%pot(start + isec)%pot%set(1)%buck4r%npoly1 = 5
1357 nonbonded%pot(start + isec)%pot%set(1)%buck4r%poly1(0:5) = x(1:6)
1358 nonbonded%pot(start + isec)%pot%set(1)%buck4r%npoly2 = 3
1359 nonbonded%pot(start + isec)%pot%set(1)%buck4r%poly2(0:3) = x(7:10)
1360 ELSE
1361 nonbonded%pot(start + isec)%pot%set(1)%buck4r%npoly1 = np1 - 1
1362 cpassert(np1 - 1 <= 10)
1363 nonbonded%pot(start + isec)%pot%set(1)%buck4r%poly1(0:np1 - 1) = coeff1(0:np1 - 1)
1364 nonbonded%pot(start + isec)%pot%set(1)%buck4r%npoly2 = np2 - 1
1365 cpassert(np2 - 1 <= 10)
1366 nonbonded%pot(start + isec)%pot%set(1)%buck4r%poly2(0:np2 - 1) = coeff2(0:np2 - 1)
1367 END IF
1368 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
1369
1370 IF (ASSOCIATED(coeff1)) THEN
1371 DEALLOCATE (coeff1)
1372 END IF
1373 IF (ASSOCIATED(coeff2)) THEN
1374 DEALLOCATE (coeff2)
1375 END IF
1376 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1377 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1378 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1379 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1380 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1381 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1382 END DO
1383
1384 END SUBROUTINE read_b4_section
1385
1386! **************************************************************************************************
1387!> \brief Reads the GENPOT - generic potential section
1388!> \param nonbonded ...
1389!> \param section ...
1390!> \param start ...
1391!> \author Teodoro Laino - 10.2006
1392! **************************************************************************************************
1393 SUBROUTINE read_gp_section(nonbonded, section, start)
1394 TYPE(pair_potential_p_type), POINTER :: nonbonded
1395 TYPE(section_vals_type), POINTER :: section
1396 INTEGER, INTENT(IN) :: start
1397
1398 CHARACTER(LEN=default_string_length), &
1399 DIMENSION(:), POINTER :: atm_names
1400 INTEGER :: isec, n_items, n_rep
1401 REAL(kind=dp) :: rcut
1402
1403 CALL section_vals_get(section, n_repetition=n_items)
1404 DO isec = 1, n_items
1405 NULLIFY (atm_names)
1406 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1407 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1408 nonbonded%pot(start + isec)%pot%type = gp_type
1409 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1410 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1411 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
1412 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1413 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1414 ! Parse the genpot info
1415 CALL get_generic_info(section, "FUNCTION", nonbonded%pot(start + isec)%pot%set(1)%gp%potential, &
1416 nonbonded%pot(start + isec)%pot%set(1)%gp%parameters, &
1417 nonbonded%pot(start + isec)%pot%set(1)%gp%values, &
1418 size_variables=1, i_rep_sec=isec)
1419 nonbonded%pot(start + isec)%pot%set(1)%gp%variables = nonbonded%pot(start + isec)%pot%set(1)%gp%parameters(1)
1420 !
1421 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1422 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1423 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1424 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1425 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1426 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1427 END DO
1428 END SUBROUTINE read_gp_section
1429
1430! **************************************************************************************************
1431!> \brief Reads the tersoff section
1432!> \param nonbonded ...
1433!> \param section ...
1434!> \param start ...
1435!> \param tersoff_section ...
1436!> \author ikuo
1437! **************************************************************************************************
1438 SUBROUTINE read_tersoff_section(nonbonded, section, start, tersoff_section)
1439 TYPE(pair_potential_p_type), POINTER :: nonbonded
1440 TYPE(section_vals_type), POINTER :: section
1441 INTEGER, INTENT(IN) :: start
1442 TYPE(section_vals_type), POINTER :: tersoff_section
1443
1444 CHARACTER(LEN=default_string_length), &
1445 DIMENSION(:), POINTER :: atm_names
1446 INTEGER :: isec, n_items, n_rep
1447 REAL(kind=dp) :: rcut, rcutsq
1448
1449 CALL section_vals_get(section, n_repetition=n_items)
1450 DO isec = 1, n_items
1451 CALL cite_reference(tersoff1988)
1452 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1453
1454 nonbonded%pot(start + isec)%pot%type = tersoff_type
1455 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1456 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1457 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1458 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1459
1460 CALL section_vals_val_get(tersoff_section, "A", i_rep_section=isec, &
1461 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%A)
1462 CALL section_vals_val_get(tersoff_section, "B", i_rep_section=isec, &
1463 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%B)
1464 CALL section_vals_val_get(tersoff_section, "lambda1", i_rep_section=isec, &
1465 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%lambda1)
1466 CALL section_vals_val_get(tersoff_section, "lambda2", i_rep_section=isec, &
1467 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%lambda2)
1468 CALL section_vals_val_get(tersoff_section, "alpha", i_rep_section=isec, &
1469 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%alpha)
1470 CALL section_vals_val_get(tersoff_section, "beta", i_rep_section=isec, &
1471 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%beta)
1472 CALL section_vals_val_get(tersoff_section, "n", i_rep_section=isec, &
1473 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%n)
1474 CALL section_vals_val_get(tersoff_section, "c", i_rep_section=isec, &
1475 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%c)
1476 CALL section_vals_val_get(tersoff_section, "d", i_rep_section=isec, &
1477 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%d)
1478 CALL section_vals_val_get(tersoff_section, "h", i_rep_section=isec, &
1479 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%h)
1480 CALL section_vals_val_get(tersoff_section, "lambda3", i_rep_section=isec, &
1481 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%lambda3)
1482 CALL section_vals_val_get(tersoff_section, "bigR", i_rep_section=isec, &
1483 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%bigR)
1484 CALL section_vals_val_get(tersoff_section, "bigD", i_rep_section=isec, &
1485 r_val=nonbonded%pot(start + isec)%pot%set(1)%tersoff%bigD)
1486
1487 rcutsq = (nonbonded%pot(start + isec)%pot%set(1)%tersoff%bigR + &
1488 nonbonded%pot(start + isec)%pot%set(1)%tersoff%bigD)**2
1489 nonbonded%pot(start + isec)%pot%set(1)%tersoff%rcutsq = rcutsq
1490 nonbonded%pot(start + isec)%pot%rcutsq = rcutsq
1491
1492 ! In case it is defined override the standard specification of RCUT
1493 CALL section_vals_val_get(tersoff_section, "RCUT", i_rep_section=isec, n_rep_val=n_rep)
1494 IF (n_rep == 1) THEN
1495 CALL section_vals_val_get(tersoff_section, "RCUT", i_rep_section=isec, r_val=rcut)
1496 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1497 END IF
1498 END DO
1499 END SUBROUTINE read_tersoff_section
1500
1501! **************************************************************************************************
1502!> \brief Reads the gal19 section
1503!> \param nonbonded ...
1504!> \param section ...
1505!> \param start ...
1506!> \param gal_section ...
1507!> \author Clabaut Paul
1508! **************************************************************************************************
1509 SUBROUTINE read_gal_section(nonbonded, section, start, gal_section)
1510 TYPE(pair_potential_p_type), POINTER :: nonbonded
1511 TYPE(section_vals_type), POINTER :: section
1512 INTEGER, INTENT(IN) :: start
1513 TYPE(section_vals_type), POINTER :: gal_section
1514
1515 CHARACTER(LEN=default_string_length), &
1516 DIMENSION(:), POINTER :: atm_names
1517 INTEGER :: iatom, isec, n_items, n_rep, nval
1518 LOGICAL :: is_ok
1519 REAL(kind=dp) :: rcut, rval
1520 REAL(kind=dp), DIMENSION(:), POINTER :: rvalues
1521 TYPE(cp_sll_val_type), POINTER :: list
1522 TYPE(section_vals_type), POINTER :: subsection
1523 TYPE(val_type), POINTER :: val
1524
1525 CALL section_vals_get(section, n_repetition=n_items)
1526 DO isec = 1, n_items
1527 CALL cite_reference(clabaut2020)
1528 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1529
1530 nonbonded%pot(start + isec)%pot%type = gal_type
1531 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1532 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1533 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1534 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1535
1536 CALL section_vals_val_get(section, "METALS", i_rep_section=isec, c_vals=atm_names)
1537 IF (any(len_trim(atm_names(:)) > 2)) THEN
1538 cpwarn("The atom name will be truncated.")
1539 END IF
1540 nonbonded%pot(start + isec)%pot%set(1)%gal%met1 = trim(atm_names(1))
1541 nonbonded%pot(start + isec)%pot%set(1)%gal%met2 = trim(atm_names(2))
1542
1543 CALL section_vals_val_get(gal_section, "epsilon", i_rep_section=isec, &
1544 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%epsilon)
1545 CALL section_vals_val_get(gal_section, "bxy", i_rep_section=isec, &
1546 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%bxy)
1547 CALL section_vals_val_get(gal_section, "bz", i_rep_section=isec, &
1548 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%bz)
1549
1550 CALL section_vals_val_get(gal_section, "r", i_rep_section=isec, r_vals=rvalues)
1551 nonbonded%pot(start + isec)%pot%set(1)%gal%r1 = rvalues(1)
1552 nonbonded%pot(start + isec)%pot%set(1)%gal%r2 = rvalues(2)
1553
1554 CALL section_vals_val_get(gal_section, "a1", i_rep_section=isec, &
1555 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%a1)
1556 CALL section_vals_val_get(gal_section, "a2", i_rep_section=isec, &
1557 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%a2)
1558 CALL section_vals_val_get(gal_section, "a3", i_rep_section=isec, &
1559 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%a3)
1560 CALL section_vals_val_get(gal_section, "a4", i_rep_section=isec, &
1561 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%a4)
1562 CALL section_vals_val_get(gal_section, "A", i_rep_section=isec, &
1563 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%a)
1564 CALL section_vals_val_get(gal_section, "B", i_rep_section=isec, &
1565 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%b)
1566 CALL section_vals_val_get(gal_section, "C", i_rep_section=isec, &
1567 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal%c)
1568 NULLIFY (list)
1569 subsection => section_vals_get_subs_vals(section, "GCN", i_rep_section=isec)
1570 CALL section_vals_val_get(subsection, "_DEFAULT_KEYWORD_", n_rep_val=nval)
1571 ALLOCATE (nonbonded%pot(start + isec)%pot%set(1)%gal%gcn(nval))
1572 CALL section_vals_list_get(subsection, "_DEFAULT_KEYWORD_", list=list)
1573 DO iatom = 1, nval
1574 ! we use only the first default_string_length characters of each line
1575 is_ok = cp_sll_val_next(list, val)
1576 CALL val_get(val, r_val=rval)
1577 ! assign values
1578 nonbonded%pot(start + isec)%pot%set(1)%gal%gcn(iatom) = rval
1579 END DO
1580
1581 CALL section_vals_val_get(gal_section, "Fit_express", i_rep_section=isec, &
1582 l_val=nonbonded%pot(start + isec)%pot%set(1)%gal%express)
1583
1584 ! ! In case it is defined override the standard specification of RCUT
1585 CALL section_vals_val_get(gal_section, "RCUT", i_rep_section=isec, n_rep_val=n_rep)
1586 IF (n_rep == 1) THEN
1587 CALL section_vals_val_get(gal_section, "RCUT", i_rep_section=isec, r_val=rcut)
1588 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1589 nonbonded%pot(start + isec)%pot%set(1)%gal%rcutsq = rcut**2
1590 END IF
1591 END DO
1592 END SUBROUTINE read_gal_section
1593
1594! **************************************************************************************************
1595!> \brief Reads the gal21 section
1596!> \param nonbonded ...
1597!> \param section ...
1598!> \param start ...
1599!> \param gal21_section ...
1600!> \author Clabaut Paul
1601! **************************************************************************************************
1602 SUBROUTINE read_gal21_section(nonbonded, section, start, gal21_section)
1603 TYPE(pair_potential_p_type), POINTER :: nonbonded
1604 TYPE(section_vals_type), POINTER :: section
1605 INTEGER, INTENT(IN) :: start
1606 TYPE(section_vals_type), POINTER :: gal21_section
1607
1608 CHARACTER(LEN=default_string_length), &
1609 DIMENSION(:), POINTER :: atm_names
1610 INTEGER :: iatom, isec, n_items, n_rep, nval
1611 LOGICAL :: is_ok
1612 REAL(kind=dp) :: rcut, rval
1613 REAL(kind=dp), DIMENSION(:), POINTER :: rvalues
1614 TYPE(cp_sll_val_type), POINTER :: list
1615 TYPE(section_vals_type), POINTER :: subsection
1616 TYPE(val_type), POINTER :: val
1617
1618 CALL section_vals_get(section, n_repetition=n_items)
1619 DO isec = 1, n_items
1620 CALL cite_reference(clabaut2021)
1621 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1622
1623 nonbonded%pot(start + isec)%pot%type = gal21_type
1624 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1625 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1626 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1627 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1628
1629 CALL section_vals_val_get(section, "METALS", i_rep_section=isec, c_vals=atm_names)
1630 IF (any(len_trim(atm_names(:)) > 2)) THEN
1631 cpwarn("The atom name will be truncated.")
1632 END IF
1633 nonbonded%pot(start + isec)%pot%set(1)%gal21%met1 = trim(atm_names(1))
1634 nonbonded%pot(start + isec)%pot%set(1)%gal21%met2 = trim(atm_names(2))
1635
1636 CALL section_vals_val_get(gal21_section, "epsilon", i_rep_section=isec, r_vals=rvalues)
1637 nonbonded%pot(start + isec)%pot%set(1)%gal21%epsilon1 = rvalues(1)
1638 nonbonded%pot(start + isec)%pot%set(1)%gal21%epsilon2 = rvalues(2)
1639 nonbonded%pot(start + isec)%pot%set(1)%gal21%epsilon3 = rvalues(3)
1640
1641 CALL section_vals_val_get(gal21_section, "bxy", i_rep_section=isec, r_vals=rvalues)
1642 nonbonded%pot(start + isec)%pot%set(1)%gal21%bxy1 = rvalues(1)
1643 nonbonded%pot(start + isec)%pot%set(1)%gal21%bxy2 = rvalues(2)
1644
1645 CALL section_vals_val_get(gal21_section, "bz", i_rep_section=isec, r_vals=rvalues)
1646 nonbonded%pot(start + isec)%pot%set(1)%gal21%bz1 = rvalues(1)
1647 nonbonded%pot(start + isec)%pot%set(1)%gal21%bz2 = rvalues(2)
1648
1649 CALL section_vals_val_get(gal21_section, "r", i_rep_section=isec, r_vals=rvalues)
1650 nonbonded%pot(start + isec)%pot%set(1)%gal21%r1 = rvalues(1)
1651 nonbonded%pot(start + isec)%pot%set(1)%gal21%r2 = rvalues(2)
1652
1653 CALL section_vals_val_get(gal21_section, "a1", i_rep_section=isec, r_vals=rvalues)
1654 nonbonded%pot(start + isec)%pot%set(1)%gal21%a11 = rvalues(1)
1655 nonbonded%pot(start + isec)%pot%set(1)%gal21%a12 = rvalues(2)
1656 nonbonded%pot(start + isec)%pot%set(1)%gal21%a13 = rvalues(3)
1657
1658 CALL section_vals_val_get(gal21_section, "a2", i_rep_section=isec, r_vals=rvalues)
1659 nonbonded%pot(start + isec)%pot%set(1)%gal21%a21 = rvalues(1)
1660 nonbonded%pot(start + isec)%pot%set(1)%gal21%a22 = rvalues(2)
1661 nonbonded%pot(start + isec)%pot%set(1)%gal21%a23 = rvalues(3)
1662
1663 CALL section_vals_val_get(gal21_section, "a3", i_rep_section=isec, r_vals=rvalues)
1664 nonbonded%pot(start + isec)%pot%set(1)%gal21%a31 = rvalues(1)
1665 nonbonded%pot(start + isec)%pot%set(1)%gal21%a32 = rvalues(2)
1666 nonbonded%pot(start + isec)%pot%set(1)%gal21%a33 = rvalues(3)
1667
1668 CALL section_vals_val_get(gal21_section, "a4", i_rep_section=isec, r_vals=rvalues)
1669 nonbonded%pot(start + isec)%pot%set(1)%gal21%a41 = rvalues(1)
1670 nonbonded%pot(start + isec)%pot%set(1)%gal21%a42 = rvalues(2)
1671 nonbonded%pot(start + isec)%pot%set(1)%gal21%a43 = rvalues(3)
1672
1673 CALL section_vals_val_get(gal21_section, "A", i_rep_section=isec, r_vals=rvalues)
1674 nonbonded%pot(start + isec)%pot%set(1)%gal21%AO1 = rvalues(1)
1675 nonbonded%pot(start + isec)%pot%set(1)%gal21%AO2 = rvalues(2)
1676
1677 CALL section_vals_val_get(gal21_section, "B", i_rep_section=isec, r_vals=rvalues)
1678 nonbonded%pot(start + isec)%pot%set(1)%gal21%BO1 = rvalues(1)
1679 nonbonded%pot(start + isec)%pot%set(1)%gal21%BO2 = rvalues(2)
1680
1681 CALL section_vals_val_get(gal21_section, "C", i_rep_section=isec, &
1682 r_val=nonbonded%pot(start + isec)%pot%set(1)%gal21%c)
1683
1684 CALL section_vals_val_get(gal21_section, "AH", i_rep_section=isec, r_vals=rvalues)
1685 nonbonded%pot(start + isec)%pot%set(1)%gal21%AH1 = rvalues(1)
1686 nonbonded%pot(start + isec)%pot%set(1)%gal21%AH2 = rvalues(2)
1687
1688 CALL section_vals_val_get(gal21_section, "BH", i_rep_section=isec, r_vals=rvalues)
1689 nonbonded%pot(start + isec)%pot%set(1)%gal21%BH1 = rvalues(1)
1690 nonbonded%pot(start + isec)%pot%set(1)%gal21%BH2 = rvalues(2)
1691
1692 NULLIFY (list)
1693 subsection => section_vals_get_subs_vals(section, "GCN", i_rep_section=isec)
1694 CALL section_vals_val_get(subsection, "_DEFAULT_KEYWORD_", n_rep_val=nval)
1695 ALLOCATE (nonbonded%pot(start + isec)%pot%set(1)%gal21%gcn(nval))
1696 CALL section_vals_list_get(subsection, "_DEFAULT_KEYWORD_", list=list)
1697 DO iatom = 1, nval
1698 ! we use only the first default_string_length characters of each line
1699 is_ok = cp_sll_val_next(list, val)
1700 CALL val_get(val, r_val=rval)
1701 ! assign values
1702 nonbonded%pot(start + isec)%pot%set(1)%gal21%gcn(iatom) = rval
1703 END DO
1704
1705 CALL section_vals_val_get(gal21_section, "Fit_express", i_rep_section=isec, &
1706 l_val=nonbonded%pot(start + isec)%pot%set(1)%gal21%express)
1707
1708 ! ! In case it is defined override the standard specification of RCUT
1709 CALL section_vals_val_get(gal21_section, "RCUT", i_rep_section=isec, n_rep_val=n_rep)
1710 IF (n_rep == 1) THEN
1711 CALL section_vals_val_get(gal21_section, "RCUT", i_rep_section=isec, r_val=rcut)
1712 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1713 nonbonded%pot(start + isec)%pot%set(1)%gal21%rcutsq = rcut**2
1714 END IF
1715 END DO
1716 END SUBROUTINE read_gal21_section
1717
1718! **************************************************************************************************
1719!> \brief Reads the siepmann section
1720!> \param nonbonded ...
1721!> \param section ...
1722!> \param start ...
1723!> \param siepmann_section ...
1724!> \author Dorothea Golze
1725! **************************************************************************************************
1726 SUBROUTINE read_siepmann_section(nonbonded, section, start, siepmann_section)
1727 TYPE(pair_potential_p_type), POINTER :: nonbonded
1728 TYPE(section_vals_type), POINTER :: section
1729 INTEGER, INTENT(IN) :: start
1730 TYPE(section_vals_type), POINTER :: siepmann_section
1731
1732 CHARACTER(LEN=default_string_length), &
1733 DIMENSION(:), POINTER :: atm_names
1734 INTEGER :: isec, n_items, n_rep
1735 REAL(kind=dp) :: rcut
1736
1737 CALL section_vals_get(section, n_repetition=n_items)
1738 DO isec = 1, n_items
1739 CALL cite_reference(siepmann1995)
1740 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1741
1742 nonbonded%pot(start + isec)%pot%type = siepmann_type
1743 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1744 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1745 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1746 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1747
1748 CALL section_vals_val_get(siepmann_section, "B", i_rep_section=isec, &
1749 r_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%B)
1750 CALL section_vals_val_get(siepmann_section, "D", i_rep_section=isec, &
1751 r_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%D)
1752 CALL section_vals_val_get(siepmann_section, "E", i_rep_section=isec, &
1753 r_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%E)
1754 CALL section_vals_val_get(siepmann_section, "F", i_rep_section=isec, &
1755 r_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%F)
1756 CALL section_vals_val_get(siepmann_section, "beta", i_rep_section=isec, &
1757 r_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%beta)
1758 CALL section_vals_val_get(siepmann_section, "ALLOW_OH_FORMATION", i_rep_section=isec, &
1759 l_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%allow_oh_formation)
1760 CALL section_vals_val_get(siepmann_section, "ALLOW_H3O_FORMATION", i_rep_section=isec, &
1761 l_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%allow_h3o_formation)
1762 CALL section_vals_val_get(siepmann_section, "ALLOW_O_FORMATION", i_rep_section=isec, &
1763 l_val=nonbonded%pot(start + isec)%pot%set(1)%siepmann%allow_o_formation)
1764
1765 ! ! In case it is defined override the standard specification of RCUT
1766 CALL section_vals_val_get(siepmann_section, "RCUT", i_rep_section=isec, n_rep_val=n_rep)
1767 IF (n_rep == 1) THEN
1768 CALL section_vals_val_get(siepmann_section, "RCUT", i_rep_section=isec, r_val=rcut)
1769 nonbonded%pot(start + isec)%pot%rcutsq = rcut**2
1770 nonbonded%pot(start + isec)%pot%set(1)%siepmann%rcutsq = rcut**2
1771 END IF
1772 END DO
1773 END SUBROUTINE read_siepmann_section
1774
1775! **************************************************************************************************
1776!> \brief Reads the Buckingham plus Morse potential section
1777!> \param nonbonded ...
1778!> \param section ...
1779!> \param start ...
1780!> \author MI
1781! **************************************************************************************************
1782 SUBROUTINE read_bm_section(nonbonded, section, start)
1783 TYPE(pair_potential_p_type), POINTER :: nonbonded
1784 TYPE(section_vals_type), POINTER :: section
1785 INTEGER, INTENT(IN) :: start
1786
1787 CHARACTER(LEN=default_string_length), &
1788 DIMENSION(:), POINTER :: atm_names
1789 INTEGER :: isec, n_items, n_rep
1790 REAL(kind=dp) :: a1, a2, b1, b2, beta, c, d, f0, r0, rcut
1791
1792 CALL section_vals_get(section, n_repetition=n_items)
1793 DO isec = 1, n_items
1794 CALL cite_reference(yamada2000)
1795 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1796 CALL section_vals_val_get(section, "F0", i_rep_section=isec, r_val=f0)
1797 CALL section_vals_val_get(section, "A1", i_rep_section=isec, r_val=a1)
1798 CALL section_vals_val_get(section, "A2", i_rep_section=isec, r_val=a2)
1799 CALL section_vals_val_get(section, "B1", i_rep_section=isec, r_val=b1)
1800 CALL section_vals_val_get(section, "B2", i_rep_section=isec, r_val=b2)
1801 CALL section_vals_val_get(section, "C", i_rep_section=isec, r_val=c)
1802 CALL section_vals_val_get(section, "D", i_rep_section=isec, r_val=d)
1803 CALL section_vals_val_get(section, "R0", i_rep_section=isec, r_val=r0)
1804 CALL section_vals_val_get(section, "Beta", i_rep_section=isec, r_val=beta)
1805 CALL section_vals_val_get(section, "RCUT", i_rep_section=isec, r_val=rcut)
1806
1807 nonbonded%pot(start + isec)%pot%type = bm_type
1808 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1809 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1810 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1811 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1812 nonbonded%pot(start + isec)%pot%set(1)%buckmo%f0 = f0
1813 nonbonded%pot(start + isec)%pot%set(1)%buckmo%a1 = a1
1814 nonbonded%pot(start + isec)%pot%set(1)%buckmo%a2 = a2
1815 nonbonded%pot(start + isec)%pot%set(1)%buckmo%b1 = b1
1816 nonbonded%pot(start + isec)%pot%set(1)%buckmo%b2 = b2
1817 nonbonded%pot(start + isec)%pot%set(1)%buckmo%c = c
1818 nonbonded%pot(start + isec)%pot%set(1)%buckmo%d = d
1819 nonbonded%pot(start + isec)%pot%set(1)%buckmo%r0 = r0
1820 nonbonded%pot(start + isec)%pot%set(1)%buckmo%beta = beta
1821 nonbonded%pot(start + isec)%pot%rcutsq = rcut*rcut
1822 !
1823 CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, n_rep_val=n_rep)
1824 IF (n_rep == 1) CALL section_vals_val_get(section, "RMIN", i_rep_section=isec, &
1825 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmin)
1826 CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, n_rep_val=n_rep)
1827 IF (n_rep == 1) CALL section_vals_val_get(section, "RMAX", i_rep_section=isec, &
1828 r_val=nonbonded%pot(start + isec)%pot%set(1)%rmax)
1829 END DO
1830 END SUBROUTINE read_bm_section
1831
1832! **************************************************************************************************
1833!> \brief Reads the TABPOT section
1834!> \param nonbonded ...
1835!> \param section ...
1836!> \param start ...
1837!> \param para_env ...
1838!> \param mm_section ...
1839!> \author Alex Mironenko, Da Teng
1840! **************************************************************************************************
1841 SUBROUTINE read_tabpot_section(nonbonded, section, start, para_env, mm_section)
1842 TYPE(pair_potential_p_type), POINTER :: nonbonded
1843 TYPE(section_vals_type), POINTER :: section
1844 INTEGER, INTENT(IN) :: start
1845 TYPE(mp_para_env_type), POINTER :: para_env
1846 TYPE(section_vals_type), POINTER :: mm_section
1847
1848 CHARACTER(LEN=default_string_length), &
1849 DIMENSION(:), POINTER :: atm_names
1850 INTEGER :: isec, n_items
1851
1852 CALL section_vals_get(section, n_repetition=n_items)
1853 DO isec = 1, n_items
1854 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
1855 nonbonded%pot(start + isec)%pot%type = tab_type
1856 nonbonded%pot(start + isec)%pot%at1 = atm_names(1)
1857 nonbonded%pot(start + isec)%pot%at2 = atm_names(2)
1858 CALL uppercase(nonbonded%pot(start + isec)%pot%at1)
1859 CALL uppercase(nonbonded%pot(start + isec)%pot%at2)
1860 CALL section_vals_val_get(section, "PARM_FILE_NAME", i_rep_section=isec, &
1861 c_val=nonbonded%pot(start + isec)%pot%set(1)%tab%tabpot_file_name)
1862 CALL read_tabpot_data(nonbonded%pot(start + isec)%pot%set(1)%tab, para_env, mm_section)
1863 nonbonded%pot(start + isec)%pot%set(1)%tab%index = isec
1864 END DO
1865 END SUBROUTINE read_tabpot_section
1866
1867! **************************************************************************************************
1868!> \brief Reads the CHARGE section
1869!> \param charge_atm ...
1870!> \param charge ...
1871!> \param section ...
1872!> \param start ...
1873!> \author teo
1874! **************************************************************************************************
1875 SUBROUTINE read_chrg_section(charge_atm, charge, section, start)
1876 CHARACTER(LEN=default_string_length), &
1877 DIMENSION(:), POINTER :: charge_atm
1878 REAL(kind=dp), DIMENSION(:), POINTER :: charge
1879 TYPE(section_vals_type), POINTER :: section
1880 INTEGER, INTENT(IN) :: start
1881
1882 CHARACTER(LEN=default_string_length) :: atm_name
1883 INTEGER :: isec, n_items
1884
1885 CALL section_vals_get(section, n_repetition=n_items)
1886 DO isec = 1, n_items
1887 CALL section_vals_val_get(section, "ATOM", i_rep_section=isec, c_val=atm_name)
1888 charge_atm(start + isec) = atm_name
1889 CALL uppercase(charge_atm(start + isec))
1890 CALL section_vals_val_get(section, "CHARGE", i_rep_section=isec, r_val=charge(start + isec))
1891 END DO
1892 END SUBROUTINE read_chrg_section
1893
1894! **************************************************************************************************
1895!> \brief Reads the POLARIZABILITY section
1896!> \param apol_atm ...
1897!> \param apol ...
1898!> \param damping_list ...
1899!> \param section ...
1900!> \param start ...
1901!> \author Marcel Baer
1902! **************************************************************************************************
1903 SUBROUTINE read_apol_section(apol_atm, apol, damping_list, section, &
1904 start)
1905 CHARACTER(LEN=default_string_length), &
1906 DIMENSION(:), POINTER :: apol_atm
1907 REAL(kind=dp), DIMENSION(:), POINTER :: apol
1908 TYPE(damping_info_type), DIMENSION(:), POINTER :: damping_list
1909 TYPE(section_vals_type), POINTER :: section
1910 INTEGER, INTENT(IN) :: start
1911
1912 CHARACTER(LEN=default_string_length) :: atm_name
1913 INTEGER :: isec, isec_damp, n_damp, n_items, &
1914 start_damp, tmp_damp
1915 TYPE(section_vals_type), POINTER :: tmp_section
1916
1917 CALL section_vals_get(section, n_repetition=n_items)
1918 NULLIFY (tmp_section)
1919 n_damp = 0
1920! *** Counts number of DIPOLE%DAMPING sections ****
1921 DO isec = 1, n_items
1922 tmp_section => section_vals_get_subs_vals(section, "DAMPING", &
1923 i_rep_section=isec)
1924 CALL section_vals_get(tmp_section, n_repetition=tmp_damp)
1925 n_damp = n_damp + tmp_damp
1926
1927 END DO
1928
1929 IF (n_damp > 0) THEN
1930 ALLOCATE (damping_list(1:n_damp))
1931 END IF
1932
1933! *** Reads DIPOLE sections *****
1934 start_damp = 0
1935 DO isec = 1, n_items
1936 CALL section_vals_val_get(section, "ATOM", i_rep_section=isec, c_val=atm_name)
1937 apol_atm(start + isec) = atm_name
1938 CALL uppercase(apol_atm(start + isec))
1939 CALL section_vals_val_get(section, "APOL", i_rep_section=isec, r_val=apol(start + isec))
1940
1941 tmp_section => section_vals_get_subs_vals(section, "DAMPING", &
1942 i_rep_section=isec)
1943 CALL section_vals_get(tmp_section, n_repetition=tmp_damp)
1944 DO isec_damp = 1, tmp_damp
1945 damping_list(start_damp + isec_damp)%atm_name1 = apol_atm(start + isec)
1946 CALL section_vals_val_get(tmp_section, "ATOM", i_rep_section=isec_damp, &
1947 c_val=atm_name)
1948 damping_list(start_damp + isec_damp)%atm_name2 = atm_name
1949 CALL uppercase(damping_list(start_damp + isec_damp)%atm_name2)
1950 CALL section_vals_val_get(tmp_section, "TYPE", i_rep_section=isec_damp, &
1951 c_val=atm_name)
1952 damping_list(start_damp + isec_damp)%dtype = atm_name
1953 CALL uppercase(damping_list(start_damp + isec_damp)%dtype)
1954
1955 CALL section_vals_val_get(tmp_section, "ORDER", i_rep_section=isec_damp, &
1956 i_val=damping_list(start_damp + isec_damp)%order)
1957 CALL section_vals_val_get(tmp_section, "BIJ", i_rep_section=isec_damp, &
1958 r_val=damping_list(start_damp + isec_damp)%bij)
1959 CALL section_vals_val_get(tmp_section, "CIJ", i_rep_section=isec_damp, &
1960 r_val=damping_list(start_damp + isec_damp)%cij)
1961 END DO
1962 start_damp = start_damp + tmp_damp
1963
1964 END DO
1965
1966 END SUBROUTINE read_apol_section
1967
1968! **************************************************************************************************
1969!> \brief Reads the QUADRUPOLE POLARIZABILITY section
1970!> \param cpol_atm ...
1971!> \param cpol ...
1972!> \param section ...
1973!> \param start ...
1974!> \author Marcel Baer
1975! **************************************************************************************************
1976 SUBROUTINE read_cpol_section(cpol_atm, cpol, section, start)
1977 CHARACTER(LEN=default_string_length), &
1978 DIMENSION(:), POINTER :: cpol_atm
1979 REAL(kind=dp), DIMENSION(:), POINTER :: cpol
1980 TYPE(section_vals_type), POINTER :: section
1981 INTEGER, INTENT(IN) :: start
1982
1983 CHARACTER(LEN=default_string_length) :: atm_name
1984 INTEGER :: isec, n_items
1985
1986 CALL section_vals_get(section, n_repetition=n_items)
1987 DO isec = 1, n_items
1988 CALL section_vals_val_get(section, "ATOM", i_rep_section=isec, c_val=atm_name)
1989 cpol_atm(start + isec) = atm_name
1990 CALL uppercase(cpol_atm(start + isec))
1991 CALL section_vals_val_get(section, "CPOL", i_rep_section=isec, r_val=cpol(start + isec))
1992 END DO
1993 END SUBROUTINE read_cpol_section
1994
1995! **************************************************************************************************
1996!> \brief Reads the SHELL section
1997!> \param shell_list ...
1998!> \param section ...
1999!> \param start ...
2000!> \author Marcella Iannuzzi
2001! **************************************************************************************************
2002 SUBROUTINE read_shell_section(shell_list, section, start)
2003
2004 TYPE(shell_p_type), DIMENSION(:), POINTER :: shell_list
2005 TYPE(section_vals_type), POINTER :: section
2006 INTEGER, INTENT(IN) :: start
2007
2008 CHARACTER(LEN=default_string_length) :: atm_name
2009 INTEGER :: i_rep, n_rep
2010 REAL(dp) :: ccharge, cutoff, k, maxdist, mfrac, &
2011 scharge
2012
2013 CALL section_vals_get(section, n_repetition=n_rep)
2014
2015 DO i_rep = 1, n_rep
2016 CALL section_vals_val_get(section, "_SECTION_PARAMETERS_", &
2017 c_val=atm_name, i_rep_section=i_rep)
2018 CALL uppercase(atm_name)
2019 shell_list(start + i_rep)%atm_name = atm_name
2020 CALL section_vals_val_get(section, "CORE_CHARGE", i_rep_section=i_rep, r_val=ccharge)
2021 shell_list(start + i_rep)%shell%charge_core = ccharge
2022 CALL section_vals_val_get(section, "SHELL_CHARGE", i_rep_section=i_rep, r_val=scharge)
2023 shell_list(start + i_rep)%shell%charge_shell = scharge
2024 CALL section_vals_val_get(section, "MASS_FRACTION", i_rep_section=i_rep, r_val=mfrac)
2025 shell_list(start + i_rep)%shell%massfrac = mfrac
2026 CALL section_vals_val_get(section, "K2_SPRING", i_rep_section=i_rep, r_val=k)
2027 IF (k < 0.0_dp) THEN
2028 CALL cp_abort(__location__, &
2029 "An invalid value was specified for the force constant k2 of the core-shell "// &
2030 "spring potential")
2031 END IF
2032 shell_list(start + i_rep)%shell%k2_spring = k
2033 CALL section_vals_val_get(section, "K4_SPRING", i_rep_section=i_rep, r_val=k)
2034 IF (k < 0.0_dp) THEN
2035 CALL cp_abort(__location__, &
2036 "An invalid value was specified for the force constant k4 of the core-shell "// &
2037 "spring potential")
2038 END IF
2039 shell_list(start + i_rep)%shell%k4_spring = k
2040 CALL section_vals_val_get(section, "MAX_DISTANCE", i_rep_section=i_rep, r_val=maxdist)
2041 shell_list(start + i_rep)%shell%max_dist = maxdist
2042 CALL section_vals_val_get(section, "SHELL_CUTOFF", i_rep_section=i_rep, r_val=cutoff)
2043 shell_list(start + i_rep)%shell%shell_cutoff = cutoff
2044 END DO
2045
2046 END SUBROUTINE read_shell_section
2047
2048! **************************************************************************************************
2049!> \brief Reads the BONDS section
2050!> \param bond_kind ...
2051!> \param bond_a ...
2052!> \param bond_b ...
2053!> \param bond_k ...
2054!> \param bond_r0 ...
2055!> \param bond_cs ...
2056!> \param section ...
2057!> \param start ...
2058!> \author teo
2059! **************************************************************************************************
2060 SUBROUTINE read_bonds_section(bond_kind, bond_a, bond_b, bond_k, bond_r0, bond_cs, section, start)
2061 INTEGER, DIMENSION(:), POINTER :: bond_kind
2062 CHARACTER(LEN=default_string_length), &
2063 DIMENSION(:), POINTER :: bond_a, bond_b
2064 REAL(kind=dp), DIMENSION(:, :), POINTER :: bond_k
2065 REAL(kind=dp), DIMENSION(:), POINTER :: bond_r0, bond_cs
2066 TYPE(section_vals_type), POINTER :: section
2067 INTEGER, INTENT(IN) :: start
2068
2069 CHARACTER(LEN=default_string_length), &
2070 DIMENSION(:), POINTER :: atm_names
2071 INTEGER :: isec, k, n_items
2072 REAL(kind=dp), DIMENSION(:), POINTER :: kvals
2073
2074 NULLIFY (kvals, atm_names)
2075 CALL section_vals_get(section, n_repetition=n_items)
2076 DO isec = 1, n_items
2077 CALL section_vals_val_get(section, "KIND", i_rep_section=isec, i_val=bond_kind(start + isec))
2078 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2079 bond_a(start + isec) = atm_names(1)
2080 bond_b(start + isec) = atm_names(2)
2081 CALL uppercase(bond_a(start + isec))
2082 CALL uppercase(bond_b(start + isec))
2083 CALL section_vals_val_get(section, "K", i_rep_section=isec, r_vals=kvals)
2084 cpassert(SIZE(kvals) <= 3)
2085 bond_k(:, start + isec) = 0.0_dp
2086 DO k = 1, SIZE(kvals)
2087 bond_k(k, start + isec) = kvals(k)
2088 END DO
2089 CALL section_vals_val_get(section, "R0", i_rep_section=isec, r_val=bond_r0(start + isec))
2090 CALL section_vals_val_get(section, "CS", i_rep_section=isec, r_val=bond_cs(start + isec))
2091 END DO
2092 END SUBROUTINE read_bonds_section
2093
2094! **************************************************************************************************
2095!> \brief Reads the BENDS section
2096!> \param bend_kind ...
2097!> \param bend_a ...
2098!> \param bend_b ...
2099!> \param bend_c ...
2100!> \param bend_k ...
2101!> \param bend_theta0 ...
2102!> \param bend_cb ...
2103!> \param bend_r012 ...
2104!> \param bend_r032 ...
2105!> \param bend_kbs12 ...
2106!> \param bend_kbs32 ...
2107!> \param bend_kss ...
2108!> \param bend_legendre ...
2109!> \param section ...
2110!> \param start ...
2111!> \author teo
2112! **************************************************************************************************
2113 SUBROUTINE read_bends_section(bend_kind, bend_a, bend_b, bend_c, bend_k, bend_theta0, bend_cb, &
2114 bend_r012, bend_r032, bend_kbs12, bend_kbs32, bend_kss, bend_legendre, &
2115 section, start)
2116 INTEGER, DIMENSION(:), POINTER :: bend_kind
2117 CHARACTER(LEN=default_string_length), &
2118 DIMENSION(:), POINTER :: bend_a, bend_b, bend_c
2119 REAL(kind=dp), DIMENSION(:), POINTER :: bend_k, bend_theta0, bend_cb, bend_r012, &
2120 bend_r032, bend_kbs12, bend_kbs32, &
2121 bend_kss
2122 TYPE(legendre_data_type), DIMENSION(:), POINTER :: bend_legendre
2123 TYPE(section_vals_type), POINTER :: section
2124 INTEGER, INTENT(IN) :: start
2125
2126 CHARACTER(LEN=default_string_length), &
2127 DIMENSION(:), POINTER :: atm_names
2128 INTEGER :: isec, k, n_items, n_rep
2129 REAL(kind=dp), DIMENSION(:), POINTER :: kvals, r_values
2130
2131 NULLIFY (kvals, atm_names)
2132 CALL section_vals_get(section, n_repetition=n_items)
2133 bend_legendre%order = 0
2134 DO isec = 1, n_items
2135 CALL section_vals_val_get(section, "KIND", i_rep_section=isec, i_val=bend_kind(start + isec))
2136 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2137 bend_a(start + isec) = atm_names(1)
2138 bend_b(start + isec) = atm_names(2)
2139 bend_c(start + isec) = atm_names(3)
2140 CALL uppercase(bend_a(start + isec))
2141 CALL uppercase(bend_b(start + isec))
2142 CALL uppercase(bend_c(start + isec))
2143 CALL section_vals_val_get(section, "K", i_rep_section=isec, r_vals=kvals)
2144 cpassert(SIZE(kvals) == 1)
2145 bend_k(start + isec) = kvals(1)
2146 CALL section_vals_val_get(section, "THETA0", i_rep_section=isec, r_val=bend_theta0(start + isec))
2147 CALL section_vals_val_get(section, "CB", i_rep_section=isec, r_val=bend_cb(start + isec))
2148 CALL section_vals_val_get(section, "R012", i_rep_section=isec, r_val=bend_r012(start + isec))
2149 CALL section_vals_val_get(section, "R032", i_rep_section=isec, r_val=bend_r032(start + isec))
2150 CALL section_vals_val_get(section, "KBS12", i_rep_section=isec, r_val=bend_kbs12(start + isec))
2151 CALL section_vals_val_get(section, "KBS32", i_rep_section=isec, r_val=bend_kbs32(start + isec))
2152 CALL section_vals_val_get(section, "KSS", i_rep_section=isec, r_val=bend_kss(start + isec))
2153 ! get legendre based data
2154 CALL section_vals_val_get(section, "LEGENDRE", i_rep_section=isec, n_rep_val=n_rep)
2155 DO k = 1, n_rep
2156 CALL section_vals_val_get(section, "LEGENDRE", i_rep_val=k, r_vals=r_values, i_rep_section=isec)
2157 bend_legendre(start + isec)%order = SIZE(r_values)
2158 IF (ASSOCIATED(bend_legendre(start + isec)%coeffs)) THEN
2159 DEALLOCATE (bend_legendre(start + isec)%coeffs)
2160 END IF
2161 ALLOCATE (bend_legendre(start + isec)%coeffs(bend_legendre(start + isec)%order))
2162 bend_legendre(start + isec)%coeffs = r_values
2163 END DO
2164 END DO
2165 END SUBROUTINE read_bends_section
2166
2167! **************************************************************************************************
2168!> \brief ...
2169!> \param ub_kind ...
2170!> \param ub_a ...
2171!> \param ub_b ...
2172!> \param ub_c ...
2173!> \param ub_k ...
2174!> \param ub_r0 ...
2175!> \param section ...
2176!> \param start ...
2177! **************************************************************************************************
2178 SUBROUTINE read_ubs_section(ub_kind, ub_a, ub_b, ub_c, ub_k, ub_r0, section, start)
2179 INTEGER, DIMENSION(:), POINTER :: ub_kind
2180 CHARACTER(LEN=default_string_length), &
2181 DIMENSION(:), POINTER :: ub_a, ub_b, ub_c
2182 REAL(kind=dp), DIMENSION(:, :), POINTER :: ub_k
2183 REAL(kind=dp), DIMENSION(:), POINTER :: ub_r0
2184 TYPE(section_vals_type), POINTER :: section
2185 INTEGER, INTENT(IN) :: start
2186
2187 CHARACTER(LEN=default_string_length), &
2188 DIMENSION(:), POINTER :: atm_names
2189 INTEGER :: isec, k, n_items
2190 LOGICAL :: explicit
2191 REAL(kind=dp), DIMENSION(:), POINTER :: kvals
2192 TYPE(section_vals_type), POINTER :: subsection
2193
2194 NULLIFY (atm_names)
2195 CALL section_vals_get(section, n_repetition=n_items)
2196 DO isec = 1, n_items
2197 subsection => section_vals_get_subs_vals(section, "UB", i_rep_section=isec)
2198 CALL section_vals_get(subsection, explicit=explicit)
2199 IF (explicit) THEN
2200 CALL section_vals_val_get(subsection, "KIND", i_val=ub_kind(start + isec))
2201 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2202 ub_a(start + isec) = atm_names(1)
2203 ub_b(start + isec) = atm_names(2)
2204 ub_c(start + isec) = atm_names(3)
2205 CALL uppercase(ub_a(start + isec))
2206 CALL uppercase(ub_b(start + isec))
2207 CALL uppercase(ub_c(start + isec))
2208 CALL section_vals_val_get(subsection, "K", r_vals=kvals)
2209 cpassert(SIZE(kvals) <= 3)
2210 ub_k(:, start + isec) = 0.0_dp
2211 DO k = 1, SIZE(kvals)
2212 ub_k(k, start + isec) = kvals(k)
2213 END DO
2214 CALL section_vals_val_get(subsection, "R0", r_val=ub_r0(start + isec))
2215 END IF
2216 END DO
2217 END SUBROUTINE read_ubs_section
2218
2219! **************************************************************************************************
2220!> \brief Reads the TORSIONS section
2221!> \param torsion_kind ...
2222!> \param torsion_a ...
2223!> \param torsion_b ...
2224!> \param torsion_c ...
2225!> \param torsion_d ...
2226!> \param torsion_k ...
2227!> \param torsion_phi0 ...
2228!> \param torsion_m ...
2229!> \param section ...
2230!> \param start ...
2231!> \author teo
2232! **************************************************************************************************
2233 SUBROUTINE read_torsions_section(torsion_kind, torsion_a, torsion_b, torsion_c, torsion_d, torsion_k, &
2234 torsion_phi0, torsion_m, section, start)
2235 INTEGER, DIMENSION(:), POINTER :: torsion_kind
2236 CHARACTER(LEN=default_string_length), &
2237 DIMENSION(:), POINTER :: torsion_a, torsion_b, torsion_c, &
2238 torsion_d
2239 REAL(kind=dp), DIMENSION(:), POINTER :: torsion_k, torsion_phi0
2240 INTEGER, DIMENSION(:), POINTER :: torsion_m
2241 TYPE(section_vals_type), POINTER :: section
2242 INTEGER, INTENT(IN) :: start
2243
2244 CHARACTER(LEN=default_string_length), &
2245 DIMENSION(:), POINTER :: atm_names
2246 INTEGER :: isec, n_items
2247
2248 NULLIFY (atm_names)
2249 CALL section_vals_get(section, n_repetition=n_items)
2250 DO isec = 1, n_items
2251 CALL section_vals_val_get(section, "KIND", i_rep_section=isec, i_val=torsion_kind(start + isec))
2252 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2253 torsion_a(start + isec) = atm_names(1)
2254 torsion_b(start + isec) = atm_names(2)
2255 torsion_c(start + isec) = atm_names(3)
2256 torsion_d(start + isec) = atm_names(4)
2257 CALL uppercase(torsion_a(start + isec))
2258 CALL uppercase(torsion_b(start + isec))
2259 CALL uppercase(torsion_c(start + isec))
2260 CALL uppercase(torsion_d(start + isec))
2261 CALL section_vals_val_get(section, "K", i_rep_section=isec, r_val=torsion_k(start + isec))
2262 CALL section_vals_val_get(section, "PHI0", i_rep_section=isec, r_val=torsion_phi0(start + isec))
2263 CALL section_vals_val_get(section, "M", i_rep_section=isec, i_val=torsion_m(start + isec))
2264 ! Modify parameterisation for OPLS case
2265 IF (torsion_kind(start + isec) == do_ff_opls) THEN
2266 IF (torsion_phi0(start + isec) /= 0.0_dp) THEN
2267 CALL cp_warn(__location__, "PHI0 parameter was non-zero "// &
2268 "for an OPLS-type TORSION. It will be ignored.")
2269 END IF
2270 IF (modulo(torsion_m(start + isec), 2) == 0) THEN
2271 ! For even M, negate the cosine using a Pi phase factor
2272 torsion_phi0(start + isec) = pi
2273 END IF
2274 ! the K parameter appears as K/2 in the OPLS parameterisation
2275 torsion_k(start + isec) = torsion_k(start + isec)*0.5_dp
2276 END IF
2277 END DO
2278 END SUBROUTINE read_torsions_section
2279
2280! **************************************************************************************************
2281!> \brief Reads the IMPROPER section
2282!> \param impr_kind ...
2283!> \param impr_a ...
2284!> \param impr_b ...
2285!> \param impr_c ...
2286!> \param impr_d ...
2287!> \param impr_k ...
2288!> \param impr_phi0 ...
2289!> \param section ...
2290!> \param start ...
2291!> \author louis vanduyfhuys
2292! **************************************************************************************************
2293 SUBROUTINE read_improper_section(impr_kind, impr_a, impr_b, impr_c, impr_d, impr_k, &
2294 impr_phi0, section, start)
2295 INTEGER, DIMENSION(:), POINTER :: impr_kind
2296 CHARACTER(LEN=default_string_length), &
2297 DIMENSION(:), POINTER :: impr_a, impr_b, impr_c, impr_d
2298 REAL(kind=dp), DIMENSION(:), POINTER :: impr_k, impr_phi0
2299 TYPE(section_vals_type), POINTER :: section
2300 INTEGER, INTENT(IN) :: start
2301
2302 CHARACTER(LEN=default_string_length), &
2303 DIMENSION(:), POINTER :: atm_names
2304 INTEGER :: isec, n_items
2305
2306 NULLIFY (atm_names)
2307 CALL section_vals_get(section, n_repetition=n_items)
2308 DO isec = 1, n_items
2309 CALL section_vals_val_get(section, "KIND", i_rep_section=isec, i_val=impr_kind(start + isec))
2310 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2311 impr_a(start + isec) = atm_names(1)
2312 impr_b(start + isec) = atm_names(2)
2313 impr_c(start + isec) = atm_names(3)
2314 impr_d(start + isec) = atm_names(4)
2315 CALL uppercase(impr_a(start + isec))
2316 CALL uppercase(impr_b(start + isec))
2317 CALL uppercase(impr_c(start + isec))
2318 CALL uppercase(impr_d(start + isec))
2319 CALL section_vals_val_get(section, "K", i_rep_section=isec, r_val=impr_k(start + isec))
2320 CALL section_vals_val_get(section, "PHI0", i_rep_section=isec, r_val=impr_phi0(start + isec))
2321 END DO
2322 END SUBROUTINE read_improper_section
2323
2324! **************************************************************************************************
2325!> \brief Reads the OPBEND section
2326!> \param opbend_kind ...
2327!> \param opbend_a ...
2328!> \param opbend_b ...
2329!> \param opbend_c ...
2330!> \param opbend_d ...
2331!> \param opbend_k ...
2332!> \param opbend_phi0 ...
2333!> \param section ...
2334!> \param start ...
2335!> \author louis vanduyfhuys
2336! **************************************************************************************************
2337 SUBROUTINE read_opbend_section(opbend_kind, opbend_a, opbend_b, opbend_c, opbend_d, opbend_k, &
2338 opbend_phi0, section, start)
2339 INTEGER, DIMENSION(:), POINTER :: opbend_kind
2340 CHARACTER(LEN=default_string_length), &
2341 DIMENSION(:), POINTER :: opbend_a, opbend_b, opbend_c, opbend_d
2342 REAL(kind=dp), DIMENSION(:), POINTER :: opbend_k, opbend_phi0
2343 TYPE(section_vals_type), POINTER :: section
2344 INTEGER, INTENT(IN) :: start
2345
2346 CHARACTER(LEN=default_string_length), &
2347 DIMENSION(:), POINTER :: atm_names
2348 INTEGER :: isec, n_items
2349
2350 NULLIFY (atm_names)
2351 CALL section_vals_get(section, n_repetition=n_items)
2352 DO isec = 1, n_items
2353 CALL section_vals_val_get(section, "KIND", i_rep_section=isec, i_val=opbend_kind(start + isec))
2354 CALL section_vals_val_get(section, "ATOMS", i_rep_section=isec, c_vals=atm_names)
2355 opbend_a(start + isec) = atm_names(1)
2356 opbend_b(start + isec) = atm_names(2)
2357 opbend_c(start + isec) = atm_names(3)
2358 opbend_d(start + isec) = atm_names(4)
2359 CALL uppercase(opbend_a(start + isec))
2360 CALL uppercase(opbend_b(start + isec))
2361 CALL uppercase(opbend_c(start + isec))
2362 CALL uppercase(opbend_d(start + isec))
2363 CALL section_vals_val_get(section, "K", i_rep_section=isec, r_val=opbend_k(start + isec))
2364 CALL section_vals_val_get(section, "PHI0", i_rep_section=isec, r_val=opbend_phi0(start + isec))
2365 END DO
2366 END SUBROUTINE read_opbend_section
2367
2368! **************************************************************************************************
2369!> \brief Reads the force_field input section
2370!> \param ff_type ...
2371!> \param para_env ...
2372!> \param mm_section ...
2373!> \par History
2374!> JGH (30.11.2001) : moved determination of setup variables to
2375!> molecule_input
2376!> \author CJM
2377! **************************************************************************************************
2378 SUBROUTINE read_force_field_section(ff_type, para_env, mm_section)
2379 TYPE(force_field_type), INTENT(INOUT) :: ff_type
2380 TYPE(mp_para_env_type), POINTER :: para_env
2381 TYPE(section_vals_type), POINTER :: mm_section
2382
2383 TYPE(section_vals_type), POINTER :: ff_section
2384
2385 NULLIFY (ff_section)
2386 ff_section => section_vals_get_subs_vals(mm_section, "FORCEFIELD")
2387 CALL read_force_field_section1(ff_section, mm_section, ff_type, para_env)
2388 END SUBROUTINE read_force_field_section
2389
2390! **************************************************************************************************
2391!> \brief reads EAM potential from library
2392!> \param eam ...
2393!> \param para_env ...
2394!> \param mm_section ...
2395! **************************************************************************************************
2396 SUBROUTINE read_eam_data(eam, para_env, mm_section)
2397 TYPE(eam_pot_type), POINTER :: eam
2398 TYPE(mp_para_env_type), POINTER :: para_env
2399 TYPE(section_vals_type), POINTER :: mm_section
2400
2401 CHARACTER(len=*), PARAMETER :: routinen = 'read_eam_data'
2402
2403 INTEGER :: handle, i, iw
2404 TYPE(cp_logger_type), POINTER :: logger
2405 TYPE(cp_parser_type) :: parser
2406
2407 CALL timeset(routinen, handle)
2408 NULLIFY (logger)
2409 logger => cp_get_default_logger()
2410 iw = cp_print_key_unit_nr(logger, mm_section, "PRINT%FF_INFO", &
2411 extension=".mmLog")
2412 IF (iw > 0) WRITE (iw, *) "Reading EAM data from: ", trim(eam%eam_file_name)
2413 CALL parser_create(parser, trim(eam%eam_file_name), para_env=para_env)
2414
2415 CALL parser_get_next_line(parser, 1)
2416 IF (iw > 0) WRITE (iw, *) "Title: ", parser%input_line
2417
2418 CALL parser_get_next_line(parser, 2)
2419 READ (parser%input_line, *) eam%drar, eam%drhoar, eam%acutal, eam%npoints
2420 eam%drar = cp_unit_to_cp2k(eam%drar, "angstrom")
2421 eam%acutal = cp_unit_to_cp2k(eam%acutal, "angstrom")
2422 ! Relocating arrays with the right size
2423 CALL reallocate(eam%rho, 1, eam%npoints)
2424 CALL reallocate(eam%rhop, 1, eam%npoints)
2425 CALL reallocate(eam%rval, 1, eam%npoints)
2426 CALL reallocate(eam%rhoval, 1, eam%npoints)
2427 CALL reallocate(eam%phi, 1, eam%npoints)
2428 CALL reallocate(eam%phip, 1, eam%npoints)
2429 CALL reallocate(eam%frho, 1, eam%npoints)
2430 CALL reallocate(eam%frhop, 1, eam%npoints)
2431 ! Reading density and derivative of density (with respect to r)
2432 DO i = 1, eam%npoints
2433 CALL parser_get_next_line(parser, 1)
2434 READ (parser%input_line, *) eam%rho(i), eam%rhop(i)
2435 eam%rhop(i) = cp_unit_to_cp2k(eam%rhop(i), "angstrom^-1")
2436 eam%rval(i) = real(i - 1, kind=dp)*eam%drar
2437 eam%rhoval(i) = real(i - 1, kind=dp)*eam%drhoar
2438 END DO
2439 ! Reading pair potential PHI and its derivative (with respect to r)
2440 DO i = 1, eam%npoints
2441 CALL parser_get_next_line(parser, 1)
2442 READ (parser%input_line, *) eam%phi(i), eam%phip(i)
2443 eam%phi(i) = cp_unit_to_cp2k(eam%phi(i), "eV")
2444 eam%phip(i) = cp_unit_to_cp2k(eam%phip(i), "eV*angstrom^-1")
2445 END DO
2446 ! Reading embedded function and its derivative (with respect to density)
2447 DO i = 1, eam%npoints
2448 CALL parser_get_next_line(parser, 1)
2449 READ (parser%input_line, *) eam%frho(i), eam%frhop(i)
2450 eam%frho(i) = cp_unit_to_cp2k(eam%frho(i), "eV")
2451 eam%frhop(i) = cp_unit_to_cp2k(eam%frhop(i), "eV")
2452 END DO
2453
2454 IF (iw > 0) WRITE (iw, *) "Finished EAM data"
2455 CALL parser_release(parser)
2456 CALL cp_print_key_finished_output(iw, logger, mm_section, "PRINT%FF_INFO")
2457 CALL timestop(handle)
2458
2459 END SUBROUTINE read_eam_data
2460
2461! **************************************************************************************************
2462!> \brief reads NequIP potential from .pth file
2463!> \param nequip ...
2464!> \author Gabriele Tocci
2465! **************************************************************************************************
2466 SUBROUTINE read_nequip_data(nequip)
2467 TYPE(nequip_pot_type) :: nequip
2468
2469 CHARACTER(LEN=*), PARAMETER :: routinen = 'read_nequip_data'
2470 CHARACTER(LEN=1), PARAMETER :: delimiter = ' '
2471
2472 CHARACTER(LEN=100), ALLOCATABLE, DIMENSION(:) :: tokenized_string
2473 CHARACTER(LEN=default_path_length) :: allow_tf32_str, cutoff_str, &
2474 default_dtype, model_dtype, types_str
2475 INTEGER :: handle
2476 LOGICAL :: allow_tf32, found
2477
2478 CALL timeset(routinen, handle)
2479
2480 INQUIRE (file=nequip%nequip_file_name, exist=found)
2481 IF (.NOT. found) THEN
2482 CALL cp_abort(__location__, &
2483 "Nequip model file <"//trim(nequip%nequip_file_name)// &
2484 "> not found.")
2485 END IF
2486
2487 nequip%nequip_version = torch_model_read_metadata(nequip%nequip_file_name, "nequip_version")
2488 cutoff_str = torch_model_read_metadata(nequip%nequip_file_name, "r_max")
2489 types_str = torch_model_read_metadata(nequip%nequip_file_name, "type_names")
2490 CALL tokenize_string(trim(types_str), delimiter, tokenized_string)
2491
2492 IF (ALLOCATED(nequip%type_names_torch)) THEN
2493 DEALLOCATE (nequip%type_names_torch)
2494 END IF
2495 ALLOCATE (nequip%type_names_torch(SIZE(tokenized_string)))
2496
2497 nequip%type_names_torch(:) = tokenized_string(:)
2498
2499 READ (cutoff_str, *) nequip%rcutsq
2500 nequip%rcutsq = cp_unit_to_cp2k(nequip%rcutsq, nequip%unit_coords)
2501 nequip%rcutsq = nequip%rcutsq*nequip%rcutsq
2502 nequip%unit_coords_val = cp_unit_to_cp2k(nequip%unit_coords_val, nequip%unit_coords)
2503 nequip%unit_forces_val = cp_unit_to_cp2k(nequip%unit_forces_val, nequip%unit_forces)
2504 nequip%unit_energy_val = cp_unit_to_cp2k(nequip%unit_energy_val, nequip%unit_energy)
2505 nequip%unit_cell_val = cp_unit_to_cp2k(nequip%unit_cell_val, nequip%unit_cell)
2506
2507 default_dtype = torch_model_read_metadata(nequip%nequip_file_name, "default_dtype")
2508 model_dtype = torch_model_read_metadata(nequip%nequip_file_name, "model_dtype")
2509 IF (trim(default_dtype) == "float32" .AND. trim(model_dtype) == "float32") THEN
2510 nequip%do_nequip_sp = .true.
2511 ELSE IF (trim(default_dtype) == "float64" .AND. trim(model_dtype) == "float64") THEN
2512 nequip%do_nequip_sp = .false.
2513 ELSE
2514 CALL cp_abort(__location__, &
2515 "Both default_dtype and model_dtype should be either float32 or float64. Currently, default_dtype is <"// &
2516 trim(default_dtype)//"> and model_dtype is <"//trim(model_dtype)//">.")
2517 END IF
2518
2519 allow_tf32_str = torch_model_read_metadata(nequip%nequip_file_name, "allow_tf32")
2520 allow_tf32 = (trim(allow_tf32_str) == "1")
2521 IF (trim(allow_tf32_str) /= "1" .AND. trim(allow_tf32_str) /= "0") THEN
2522 CALL cp_abort(__location__, &
2523 "The value for allow_tf32 <"//trim(allow_tf32_str)// &
2524 "> is not supported. Check the .yaml and .pth files.")
2525 END IF
2526 CALL torch_allow_tf32(allow_tf32)
2527
2528 CALL timestop(handle)
2529 END SUBROUTINE read_nequip_data
2530
2531! **************************************************************************************************
2532!> \brief reads ALLEGRO potential from .pth file
2533!> \param allegro ...
2534!> \author Gabriele Tocci
2535! **************************************************************************************************
2536 SUBROUTINE read_allegro_data(allegro)
2537 TYPE(allegro_pot_type) :: allegro
2538
2539 CHARACTER(len=*), PARAMETER :: routinen = 'read_allegro_data'
2540 CHARACTER(LEN=1), PARAMETER :: delimiter = ' '
2541
2542 CHARACTER(LEN=100), ALLOCATABLE, DIMENSION(:) :: tokenized_string
2543 CHARACTER(LEN=default_path_length) :: allow_tf32_str, cutoff_str, &
2544 default_dtype, model_dtype, types_str
2545 INTEGER :: handle
2546 LOGICAL :: allow_tf32, found
2547
2548 CALL timeset(routinen, handle)
2549
2550 INQUIRE (file=allegro%allegro_file_name, exist=found)
2551 IF (.NOT. found) THEN
2552 CALL cp_abort(__location__, &
2553 "Allegro model file <"//trim(allegro%allegro_file_name)// &
2554 "> not found.")
2555 END IF
2556
2557 allegro%nequip_version = torch_model_read_metadata(allegro%allegro_file_name, "nequip_version")
2558 IF (allegro%nequip_version == "") THEN
2559 CALL cp_abort(__location__, &
2560 "Allegro model file <"//trim(allegro%allegro_file_name)// &
2561 "> has not been deployed; did you forget to run `nequip-deploy`?")
2562 END IF
2563 cutoff_str = torch_model_read_metadata(allegro%allegro_file_name, "r_max")
2564 types_str = torch_model_read_metadata(allegro%allegro_file_name, "type_names")
2565 CALL tokenize_string(trim(types_str), delimiter, tokenized_string)
2566
2567 IF (ALLOCATED(allegro%type_names_torch)) THEN
2568 DEALLOCATE (allegro%type_names_torch)
2569 END IF
2570 ALLOCATE (allegro%type_names_torch(SIZE(tokenized_string)))
2571 allegro%type_names_torch(:) = tokenized_string(:)
2572
2573 READ (cutoff_str, *) allegro%rcutsq
2574 allegro%rcutsq = cp_unit_to_cp2k(allegro%rcutsq, allegro%unit_coords)
2575 allegro%rcutsq = allegro%rcutsq*allegro%rcutsq
2576 allegro%unit_coords_val = cp_unit_to_cp2k(allegro%unit_coords_val, allegro%unit_coords)
2577 allegro%unit_forces_val = cp_unit_to_cp2k(allegro%unit_forces_val, allegro%unit_forces)
2578 allegro%unit_energy_val = cp_unit_to_cp2k(allegro%unit_energy_val, allegro%unit_energy)
2579 allegro%unit_cell_val = cp_unit_to_cp2k(allegro%unit_cell_val, allegro%unit_cell)
2580
2581 default_dtype = torch_model_read_metadata(allegro%allegro_file_name, "default_dtype")
2582 model_dtype = torch_model_read_metadata(allegro%allegro_file_name, "model_dtype")
2583 IF (trim(default_dtype) == "float32" .AND. trim(model_dtype) == "float32") THEN
2584 allegro%do_allegro_sp = .true.
2585 ELSE IF (trim(default_dtype) == "float64" .AND. trim(model_dtype) == "float64") THEN
2586 allegro%do_allegro_sp = .false.
2587 ELSE
2588 CALL cp_abort(__location__, &
2589 "Both default_dtype and model_dtype should be either float32 or float64. Currently, default_dtype is <"// &
2590 trim(default_dtype)//"> and model_dtype is <"//trim(model_dtype)//">.")
2591 END IF
2592
2593 allow_tf32_str = torch_model_read_metadata(allegro%allegro_file_name, "allow_tf32")
2594 allow_tf32 = (trim(allow_tf32_str) == "1")
2595 IF (trim(allow_tf32_str) /= "1" .AND. trim(allow_tf32_str) /= "0") THEN
2596 CALL cp_abort(__location__, &
2597 "The value for allow_tf32 <"//trim(allow_tf32_str)// &
2598 "> is not supported. Check the .yaml and .pth files.")
2599 END IF
2600 CALL torch_allow_tf32(allow_tf32)
2601
2602 CALL timestop(handle)
2603 END SUBROUTINE read_allegro_data
2604
2605! **************************************************************************************************
2606!> \brief returns tokenized string of kinds from .pth file
2607!> \param element ...
2608!> \param delimiter ...
2609!> \param tokenized_array ...
2610!> \author Maria Bilichenko
2611! **************************************************************************************************
2612 SUBROUTINE tokenize_string(element, delimiter, tokenized_array)
2613 CHARACTER(LEN=*), INTENT(IN) :: element
2614 CHARACTER(LEN=1), INTENT(IN) :: delimiter
2615 CHARACTER(LEN=100), ALLOCATABLE, DIMENSION(:), &
2616 INTENT(OUT) :: tokenized_array
2617
2618 CHARACTER(LEN=100) :: temp_kinds
2619 INTEGER :: end_pos, i, num_elements, start
2620 LOGICAL, ALLOCATABLE, DIMENSION(:) :: delim_positions
2621
2622 ! Find positions of delimiter within element
2623 ALLOCATE (delim_positions(len(element)))
2624 delim_positions = .false.
2625
2626 DO i = 1, len(element)
2627 IF (element(i:i) == delimiter) delim_positions(i) = .true.
2628 END DO
2629
2630 num_elements = count(delim_positions) + 1
2631
2632 ALLOCATE (tokenized_array(num_elements))
2633
2634 start = 1
2635 DO i = 1, num_elements
2636 IF (len(element) < 3 .AND. count(delim_positions) == 0) THEN ! if there is 1 kind only and it has one or two
2637 !characters (C or Cl) the end_pos will be the index of the last character (1 or 2)
2638 end_pos = len(element)
2639 ELSE ! else, the end_pos is determined by the index of the space - 1
2640 end_pos = find_end_pos(start, delim_positions)
2641 END IF
2642 temp_kinds = element(start:end_pos)
2643 IF (trim(temp_kinds) /= '') THEN
2644 tokenized_array(i) = temp_kinds
2645 END IF
2646 start = end_pos + 2
2647 END DO
2648 DEALLOCATE (delim_positions)
2649 END SUBROUTINE tokenize_string
2650
2651! **************************************************************************************************
2652!> \brief finds the position of the atom by the spacing
2653!> \param start ...
2654!> \param delim_positions ...
2655!> \return ...
2656!> \author Maria Bilichenko
2657! **************************************************************************************************
2658 INTEGER FUNCTION find_end_pos(start, delim_positions)
2659 INTEGER, INTENT(IN) :: start
2660 LOGICAL, DIMENSION(:), INTENT(IN) :: delim_positions
2661
2662 INTEGER :: end_pos, i
2663
2664 end_pos = start
2665 DO i = start, SIZE(delim_positions)
2666 IF (delim_positions(i)) THEN
2667 end_pos = i - 1
2668 EXIT
2669 END IF
2670 END DO
2671
2672 find_end_pos = end_pos
2673 END FUNCTION find_end_pos
2674
2675! **************************************************************************************************
2676!> \brief checks if all the ATOMS from *.inp file are available in *.pth file
2677!> \param cp2k_inp_atom_types ...
2678!> \param torch_atom_types ...
2679!> \author Maria Bilichenko
2680! **************************************************************************************************
2681 SUBROUTINE check_cp2k_atom_names_in_torch(cp2k_inp_atom_types, torch_atom_types)
2682 CHARACTER(LEN=*), DIMENSION(:), INTENT(IN) :: cp2k_inp_atom_types, torch_atom_types
2683
2684 INTEGER :: i, j
2685 LOGICAL :: found
2686
2687 DO i = 1, SIZE(cp2k_inp_atom_types)
2688 found = .false.
2689 DO j = 1, SIZE(torch_atom_types)
2690 IF (trim(cp2k_inp_atom_types(i)) == trim(torch_atom_types(j))) THEN
2691 found = .true.
2692 EXIT
2693 END IF
2694 END DO
2695 IF (.NOT. found) THEN
2696 CALL cp_abort(__location__, &
2697 "Atom "//trim(cp2k_inp_atom_types(i))// &
2698 " is defined in the CP2K input file but is missing in the torch model file")
2699 END IF
2700 END DO
2701 END SUBROUTINE check_cp2k_atom_names_in_torch
2702
2703! **************************************************************************************************
2704!> \brief reads TABPOT potential from file
2705!> \param tab ...
2706!> \param para_env ...
2707!> \param mm_section ...
2708!> \author Da Teng, Alex Mironenko
2709! **************************************************************************************************
2710 SUBROUTINE read_tabpot_data(tab, para_env, mm_section)
2711 TYPE(tab_pot_type), POINTER :: tab
2712 TYPE(mp_para_env_type), POINTER :: para_env
2713 TYPE(section_vals_type), POINTER :: mm_section
2714
2715 CHARACTER(len=*), PARAMETER :: routinen = 'read_tabpot_data'
2716
2717 CHARACTER :: d1, d2
2718 INTEGER :: d, handle, i, iw
2719 TYPE(cp_logger_type), POINTER :: logger
2720 TYPE(cp_parser_type) :: parser
2721
2722 CALL timeset(routinen, handle)
2723 NULLIFY (logger)
2724 logger => cp_get_default_logger()
2725 iw = cp_print_key_unit_nr(logger, mm_section, "PRINT%FF_INFO", &
2726 extension=".mmLog")
2727 IF (iw > 0) WRITE (iw, *) "Reading TABPOT data from: ", trim(tab%tabpot_file_name)
2728 CALL parser_create(parser, trim(tab%tabpot_file_name), para_env=para_env)
2729 CALL parser_get_next_line(parser, 1)
2730 IF (iw > 0) WRITE (iw, *) "Title: ", parser%input_line
2731 CALL parser_get_next_line(parser, 1)
2732
2733 ! example format: N 1000 R 1.00 20.0
2734 ! Assume the data is evenly spaced
2735 READ (parser%input_line, *) d1, tab%npoints, d2, tab%dr, tab%rcut
2736
2737 ! Relocating arrays with the right size
2738 CALL reallocate(tab%r, 1, tab%npoints)
2739 CALL reallocate(tab%e, 1, tab%npoints)
2740 CALL reallocate(tab%f, 1, tab%npoints)
2741
2742 ! Reading r, e, f
2743 DO i = 1, tab%npoints
2744 CALL parser_get_next_line(parser, 1)
2745 READ (parser%input_line, *) d, tab%r(i), tab%e(i), tab%f(i)
2746 tab%r(i) = cp_unit_to_cp2k(tab%r(i), "angstrom")
2747 tab%e(i) = cp_unit_to_cp2k(tab%e(i), "kcalmol")
2748 tab%f(i) = cp_unit_to_cp2k(tab%f(i), "kcalmol*angstrom^-1")
2749 END DO
2750
2751 tab%dr = tab%r(2) - tab%r(1)
2752 tab%rcut = cp_unit_to_cp2k(tab%rcut, "angstrom")
2753
2754 IF (iw > 0) WRITE (iw, *) "Finished TABPOT data"
2755 CALL parser_release(parser)
2756 CALL cp_print_key_finished_output(iw, logger, mm_section, "PRINT%FF_INFO")
2757 CALL timestop(handle)
2758 END SUBROUTINE read_tabpot_data
2759END MODULE force_fields_input
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
Interface to ACE C wrapper.
Definition ace_wrapper.F:12
subroutine, public ace_model_initialize(ntypec, symbolc, fname, rcutc, model)
...
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public tosi1964b
integer, save, public tersoff1988
integer, save, public tosi1964a
integer, save, public siepmann1995
integer, save, public yamada2000
integer, save, public clabaut2021
integer, save, public clabaut2020
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
character(len=default_path_length) function, public discover_file(file_name)
Checks various locations for a file name.
Definition cp_files.F:521
logical function, public cp_sll_val_next(iterator, el_att)
returns true if the actual element is valid (i.e. iterator ont at end) moves the iterator to the next...
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_get_next_line(parser, nline, at_end)
Read the next input line and broadcast the input information. Skip (nline-1) lines and skip also all ...
Utility routines to read data from files. Kept as close as possible to the old parser because.
subroutine, public parser_release(parser)
releases the parser
subroutine, public parser_create(parser, file_name, unit_nr, para_env, end_section_label, separator_chars, comment_char, continuation_char, quote_char, section_char, parse_white_lines, initial_variables, apply_preprocessing)
Start a parser run. Initial variables allow to @SET stuff before opening the file.
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
Define all structure types related to force field kinds.
integer, parameter, public do_ff_undef
integer, parameter, public do_ff_charmm
integer, parameter, public do_ff_g87
integer, parameter, public do_ff_g96
integer, parameter, public do_ff_amber
integer, parameter, public do_ff_opls
Define all structures types related to force_fields.
subroutine, public read_gd_section(nonbonded, section, start)
Reads the GOODWIN section.
subroutine, public read_force_field_section(ff_type, para_env, mm_section)
Reads the force_field input section.
subroutine, public read_gp_section(nonbonded, section, start)
Reads the GENPOT - generic potential section.
subroutine, public read_wl_section(nonbonded, section, start)
Reads the WILLIAMS section.
subroutine, public read_chrg_section(charge_atm, charge, section, start)
Reads the CHARGE section.
subroutine, public read_lj_section(nonbonded, section, start)
Reads the LJ section.
subroutine, public get_generic_info(gen_section, func_name, xfunction, parameters, values, var_values, size_variables, i_rep_sec, input_variables)
Reads from the input structure all information for generic functions.
objects that represent the structure of input sections and the data contained in an input section
subroutine, public section_vals_list_get(section_vals, keyword_name, i_rep_section, list)
returns the requested list
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
a wrapper for basic fortran types.
subroutine, public val_get(val, has_l, has_i, has_r, has_lc, has_c, l_val, l_vals, i_val, i_vals, r_val, r_vals, c_val, c_vals, len_c, type_of_var, enum)
returns the stored values
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
integer, parameter, public default_path_length
Definition kinds.F:58
An array-based list which grows on demand. When the internal array is full, a new array of twice the ...
Definition list.F:24
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Collection of simple mathematical functions and subroutines.
Definition mathlib.F:15
Utility routines for the memory handling.
Interface to the message passing library MPI.
integer, parameter, public lj_charmm_type
integer, parameter, public allegro_type
integer, parameter, public bm_type
integer, parameter, public gal_type
integer, parameter, public nequip_type
integer, parameter, public wl_type
integer, parameter, public ft_type
integer, parameter, public tab_type
subroutine, public pair_potential_reallocate(p, lb1_new, ub1_new, lj, lj_charmm, williams, goodwin, eam, nequip, allegro, bmhft, bmhftd, ipbv, buck4r, buckmo, gp, tersoff, siepmann, gal, gal21, tab, deepmd, ace)
Cleans the potential parameter type.
integer, parameter, public ftd_type
integer, parameter, public ip_type
integer, parameter, public deepmd_type
integer, parameter, public gp_type
integer, parameter, public siepmann_type
integer, parameter, public ace_type
integer, dimension(2), parameter, public do_potential_single_allocation
integer, parameter, public gw_type
integer, dimension(2), parameter, public no_potential_single_allocation
integer, parameter, public b4_type
integer, parameter, public gal21_type
integer, dimension(2), public potential_single_allocation
integer, parameter, public ea_type
integer, parameter, public tersoff_type
subroutine, public shell_p_create(shell_list, ndim)
...
Utilities for string manipulations.
elemental subroutine, public uppercase(string)
Convert all lower case characters in a string to upper case.
subroutine, public torch_allow_tf32(allow_tf32)
Set whether to allow the use of TF32. Needed due to changes in defaults from pytorch 1....
Definition torch_api.F:1599
character(:) function, allocatable, public torch_model_read_metadata(filename, key)
Reads metadata entry from given "*.pth" file. (In Torch lingo they are called extra files)
Definition torch_api.F:1518
represent a single linked list that stores pointers to the elements
type of a logger, at the moment it contains just a print level starting at which level it should be l...
a type to have a wrapper that stores any basic fortran type
stores all the informations relevant to an mpi environment