2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, 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/arrayref.h"
40 #include "gromacs/utility/basedefinitions.h"
41 #include "gromacs/utility/real.h"
48 class RangePartitioning;
53 enum class KernelType;
57 struct VerletbufListSetup
59 int cluster_size_i; /* Cluster pair-list i-cluster size atom count */
60 int cluster_size_j; /* Cluster pair-list j-cluster size atom count */
64 /* Add a 5% and 10% rlist buffer for simulations without dynamics (EM, NM, ...)
65 * and NVE simulations with zero initial temperature, respectively.
66 * 10% should be enough for any NVE simulation with PME and nstlist=10,
67 * for other settings it might not be enough, but then it's difficult
68 * to come up with any reasonable (not crazily expensive) value
69 * and grompp will notify the user when using the 10% buffer.
71 static const real verlet_buffer_ratio_nodynamics = 0.05;
72 static const real verlet_buffer_ratio_NVE_T0 = 0.10;
75 /* Returns the pair-list setup for the given nbnxn kernel type.
77 VerletbufListSetup verletbufGetListSetup(Nbnxm::KernelType nbnxnKernelType);
79 /* Enum for choosing the list type for verletbufGetSafeListSetup() */
80 enum class ListSetupType
82 CpuNoSimd, /* CPU Plain-C 4x4 list */
83 CpuSimdWhenSupported, /* CPU 4xN list, where N=4 when the binary doesn't support SIMD or the smallest N supported by SIMD in this binary */
84 Gpu /* GPU (8x2x)8x4 list */
87 /* Returns the pair-list setup assumed for the current Gromacs configuration.
88 * The setup with smallest cluster sizes is returned, such that the Verlet
89 * buffer size estimated with this setup will be conservative.
91 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType);
93 /* Returns the non-bonded pair-list radius including computed buffer
95 * Calculate the non-bonded pair-list buffer size for the Verlet list
96 * based on the particle masses, temperature, LJ types, charges
97 * and constraints as well as the non-bonded force behavior at the cut-off.
98 * The pair list update frequency and the list lifetime, which is nstlist-1
99 * for normal pair-list buffering, are passed separately, as in some cases
100 * we want an estimate for different values than the ones set in the inputrec.
101 * If referenceTemperature < 0, the maximum coupling temperature will be used.
102 * The target is a maximum average energy jump per atom of
103 * inputrec.verletbuf_tol*nstlist*inputrec.delta_t over the lifetime of the list.
105 * \note For non-linear virtual sites it can be problematic to determine their
106 * contribution to the drift exaclty, so we approximate.
108 * \param[in] mtop The system topology
109 * \param[in] boxVolume The volume of the unit cell
110 * \param[in] inputrec The input record
111 * \param[in] nstlist The pair list update frequency in steps (is not taken from \p inputrec)
112 * \param[in] listLifetime The lifetime of the pair-list, usually nstlist-1, but could be different
113 * for dynamic pruning \param[in] referenceTemperature The reference temperature for the ensemble
114 * \param[in] listSetup The pair-list setup
115 * \returns The computed pair-list radius including buffer
117 real calcVerletBufferSize(const gmx_mtop_t& mtop,
119 const t_inputrec& inputrec,
122 real referenceTemperature,
123 const VerletbufListSetup& listSetup);
125 /* Convenience type */
126 using PartitioningPerMoltype = gmx::ArrayRef<const gmx::RangePartitioning>;
128 /* Determines the mininum cell size based on atom displacement
130 * The value returned is the minimum size for which the chance that
131 * an atom or update group crosses to non nearest-neighbor cells
132 * is <= chanceRequested within ir.nstlist steps.
133 * Update groups are used when !updateGrouping.empty().
134 * Without T-coupling, SD or BD, we can not estimate atom displacements
135 * and fall back to the, crude, estimate of using the pairlist buffer size.
137 * Note: Like the Verlet buffer estimate, this estimate is based on
138 * non-interacting atoms and constrained atom-pairs. Therefore for
139 * any system that is not an ideal gas, this will be an overestimate.
141 * Note: This size increases (very slowly) with system size.
143 real minCellSizeForAtomDisplacement(const gmx_mtop_t& mtop,
144 const t_inputrec& ir,
145 PartitioningPerMoltype updateGrouping,
146 real chanceRequested);
148 /* Struct for unique atom type for calculating the energy drift.
149 * The atom displacement depends on mass and constraints.
150 * The energy jump for given distance depend on LJ type and q.
152 struct atom_nonbonded_kinetic_prop_t
154 real mass = 0; /* mass */
155 int type = 0; /* type (used for LJ parameters) */
156 real q = 0; /* charge */
157 bool bConstr = false; /* constrained, if TRUE, use #DOF=2 iso 3 */
158 real con_mass = 0; /* mass of heaviest atom connected by constraints */
159 real con_len = 0; /* constraint length to the heaviest atom */
162 /* This function computes two components of the estimate of the variance
163 * in the displacement of one atom in a system of two constrained atoms.
164 * Returns in sigma2_2d the variance due to rotation of the constrained
165 * atom around the atom to which it constrained.
166 * Returns in sigma2_3d the variance due to displacement of the COM
167 * of the whole system of the two constrained atoms.
169 * Only exposed here for testing purposes.
171 void constrained_atom_sigma2(real kT_fac, const atom_nonbonded_kinetic_prop_t* prop, real* sigma2_2d, real* sigma2_3d);