a013cf2bf9aef0bf49cd541145ebf763b9d23331
[alexxy/gromacs.git] / src / gromacs / mdlib / calc_verletbuf.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_MDLIB_CALC_VERLETBUF_H
38 #define GMX_MDLIB_CALC_VERLETBUF_H
39
40 #include "gromacs/utility/basedefinitions.h"
41 #include "gromacs/utility/real.h"
42
43 struct gmx_mtop_t;
44 struct t_inputrec;
45
46 namespace gmx
47 {
48 template<typename>
49 class ArrayRef;
50 class RangePartitioning;
51 } // namespace gmx
52
53 namespace Nbnxm
54 {
55 enum class KernelType;
56 } // namespace Nbnxm
57
58
59 struct VerletbufListSetup
60 {
61     int cluster_size_i; /* Cluster pair-list i-cluster size atom count */
62     int cluster_size_j; /* Cluster pair-list j-cluster size atom count */
63 };
64
65
66 /* Add a 5% and 10% rlist buffer for simulations without dynamics (EM, NM, ...)
67  * and NVE simulations with zero initial temperature, respectively.
68  * 10% should be enough for any NVE simulation with PME and nstlist=10,
69  * for other settings it might not be enough, but then it's difficult
70  * to come up with any reasonable (not crazily expensive) value
71  * and grompp will notify the user when using the 10% buffer.
72  */
73 static const real verlet_buffer_ratio_nodynamics = 0.05;
74 static const real verlet_buffer_ratio_NVE_T0     = 0.10;
75
76
77 /* Returns the pair-list setup for the given nbnxn kernel type.
78  */
79 VerletbufListSetup verletbufGetListSetup(Nbnxm::KernelType nbnxnKernelType);
80
81 /* Enum for choosing the list type for verletbufGetSafeListSetup() */
82 enum class ListSetupType
83 {
84     CpuNoSimd,            /* CPU Plain-C 4x4 list */
85     CpuSimdWhenSupported, /* CPU 4xN list, where N=4 when the binary doesn't support SIMD or the smallest N supported by SIMD in this binary */
86     Gpu                   /* GPU (8x2x)8x4 list */
87 };
88
89 /* Returns the pair-list setup assumed for the current Gromacs configuration.
90  * The setup with smallest cluster sizes is returned, such that the Verlet
91  * buffer size estimated with this setup will be conservative.
92  */
93 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType);
94
95 /* Returns the non-bonded pair-list radius including computed buffer
96  *
97  * Calculate the non-bonded pair-list buffer size for the Verlet list
98  * based on the particle masses, temperature, LJ types, charges
99  * and constraints as well as the non-bonded force behavior at the cut-off.
100  * The pair list update frequency and the list lifetime, which is nstlist-1
101  * for normal pair-list buffering, are passed separately, as in some cases
102  * we want an estimate for different values than the ones set in the inputrec.
103  * If referenceTemperature < 0, the maximum coupling temperature will be used.
104  * The target is a maximum average energy jump per atom of
105  * inputrec.verletbuf_tol*nstlist*inputrec.delta_t over the lifetime of the list.
106  *
107  * \note For non-linear virtual sites it can be problematic to determine their
108  *       contribution to the drift exaclty, so we approximate.
109  *
110  * \param[in] mtop          The system topology
111  * \param[in] boxVolume     The volume of the unit cell
112  * \param[in] inputrec      The input record
113  * \param[in] nstlist       The pair list update frequency in steps (is not taken from \p inputrec)
114  * \param[in] listLifetime  The lifetime of the pair-list, usually nstlist-1, but could be different
115  * for dynamic pruning \param[in] referenceTemperature  The reference temperature for the ensemble
116  * \param[in] listSetup     The pair-list setup
117  * \returns The computed pair-list radius including buffer
118  */
119 real calcVerletBufferSize(const gmx_mtop_t&         mtop,
120                           real                      boxVolume,
121                           const t_inputrec&         inputrec,
122                           int                       nstlist,
123                           int                       listLifetime,
124                           real                      referenceTemperature,
125                           const VerletbufListSetup& listSetup);
126
127 /* Convenience type */
128 using PartitioningPerMoltype = gmx::ArrayRef<const gmx::RangePartitioning>;
129
130 /* Determines the mininum cell size based on atom displacement
131  *
132  * The value returned is the minimum size for which the chance that
133  * an atom or update group crosses to non nearest-neighbor cells
134  * is <= chanceRequested within ir.nstlist steps.
135  * Update groups are used when !updateGrouping.empty().
136  * Without T-coupling, SD or BD, we can not estimate atom displacements
137  * and fall back to the, crude, estimate of using the pairlist buffer size.
138  *
139  * Note: Like the Verlet buffer estimate, this estimate is based on
140  *       non-interacting atoms and constrained atom-pairs. Therefore for
141  *       any system that is not an ideal gas, this will be an overestimate.
142  *
143  * Note: This size increases (very slowly) with system size.
144  */
145 real minCellSizeForAtomDisplacement(const gmx_mtop_t&      mtop,
146                                     const t_inputrec&      ir,
147                                     PartitioningPerMoltype updateGrouping,
148                                     real                   chanceRequested);
149
150 /* Struct for unique atom type for calculating the energy drift.
151  * The atom displacement depends on mass and constraints.
152  * The energy jump for given distance depend on LJ type and q.
153  */
154 struct atom_nonbonded_kinetic_prop_t
155 {
156     real mass     = 0;     /* mass */
157     int  type     = 0;     /* type (used for LJ parameters) */
158     real q        = 0;     /* charge */
159     bool bConstr  = false; /* constrained, if TRUE, use #DOF=2 iso 3 */
160     real con_mass = 0;     /* mass of heaviest atom connected by constraints */
161     real con_len  = 0;     /* constraint length to the heaviest atom */
162 };
163
164 /* This function computes two components of the estimate of the variance
165  * in the displacement of one atom in a system of two constrained atoms.
166  * Returns in sigma2_2d the variance due to rotation of the constrained
167  * atom around the atom to which it constrained.
168  * Returns in sigma2_3d the variance due to displacement of the COM
169  * of the whole system of the two constrained atoms.
170  *
171  * Only exposed here for testing purposes.
172  */
173 void constrained_atom_sigma2(real kT_fac, const atom_nonbonded_kinetic_prop_t* prop, real* sigma2_2d, real* sigma2_3d);
174
175 #endif