Improve Verlet buffer constraint estimate
[gromacs.git] / src / gromacs / mdlib / calc_verletbuf.h
blob097e900d788855d500e6b61e3bf2a4706e16c6f0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 #ifndef GMX_MDLIB_CALC_VERLETBUF_H
37 #define GMX_MDLIB_CALC_VERLETBUF_H
39 #include "gromacs/utility/basedefinitions.h"
40 #include "gromacs/utility/real.h"
42 struct gmx_mtop_t;
43 struct t_inputrec;
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
49 typedef struct
51 int cluster_size_i; /* Cluster pair-list i-cluster size atom count */
52 int cluster_size_j; /* Cluster pair-list j-cluster size atom count */
53 } verletbuf_list_setup_t;
56 /* Add a 5% and 10% rlist buffer for simulations without dynamics (EM, NM, ...)
57 * and NVE simulations with zero initial temperature, respectively.
58 * 10% should be enough for any NVE simulation with PME and nstlist=10,
59 * for other settings it might not be enough, but then it's difficult
60 * to come up with any reasonable (not crazily expensive) value
61 * and grompp will notify the user when using the 10% buffer.
63 static const real verlet_buffer_ratio_nodynamics = 0.05;
64 static const real verlet_buffer_ratio_NVE_T0 = 0.10;
67 /* Sets the pair-list setup assumed for the current Gromacs configuration.
68 * The setup with smallest cluster sizes is return, such that the Verlet
69 * buffer size estimated with this setup will be conservative.
70 * bSIMD tells if to take into account SIMD, when supported.
71 * bGPU tells to estimate for GPU kernels (bSIMD is ignored with bGPU=TRUE)
73 void verletbuf_get_list_setup(gmx_bool bSIMD,
74 gmx_bool bGPU,
75 verletbuf_list_setup_t *list_setup);
78 /* Calculate the non-bonded pair-list buffer size for the Verlet list
79 * based on the particle masses, temperature, LJ types, charges
80 * and constraints as well as the non-bonded force behavior at the cut-off.
81 * If reference_temperature < 0, the maximum coupling temperature will be used.
82 * The target is a maximum energy drift of ir->verletbuf_tol.
83 * Returns the number of non-linear virtual sites. For these it's difficult
84 * to determine their contribution to the drift exaclty, so we approximate.
85 * Returns the pair-list cut-off.
87 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
88 const t_inputrec *ir,
89 real reference_temperature,
90 const verletbuf_list_setup_t *list_setup,
91 int *n_nonlin_vsite,
92 real *rlist);
94 /* Struct for unique atom type for calculating the energy drift.
95 * The atom displacement depends on mass and constraints.
96 * The energy jump for given distance depend on LJ type and q.
98 struct atom_nonbonded_kinetic_prop_t
100 real mass; /* mass */
101 int type; /* type (used for LJ parameters) */
102 real q; /* charge */
103 gmx_bool bConstr; /* constrained, if TRUE, use #DOF=2 iso 3 */
104 real con_mass; /* mass of heaviest atom connected by constraints */
105 real con_len; /* constraint length to the heaviest atom */
108 /* This function computes two components of the estimate of the variance
109 * in the displacement of one atom in a system of two constrained atoms.
110 * Returns in sigma2_2d the variance due to rotation of the constrained
111 * atom around the atom to which it constrained.
112 * Returns in sigma2_3d the variance due to displacement of the COM
113 * of the whole system of the two constrained atoms.
115 * Only exposed here for testing purposes.
117 void constrained_atom_sigma2(real kT_fac,
118 const atom_nonbonded_kinetic_prop_t *prop,
119 real *sigma2_2d,
120 real *sigma2_3d);
122 #ifdef __cplusplus
124 #endif
126 #endif