a444b2e097d7c3e5195b372c99b45ee42f218953
[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,2018, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 #ifndef GMX_MDLIB_CALC_VERLETBUF_H
37 #define GMX_MDLIB_CALC_VERLETBUF_H
38
39 #include "gromacs/utility/arrayref.h"
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 class RangePartitioning;
49 } // namespace gmx
50
51 struct VerletbufListSetup
52 {
53     int  cluster_size_i;  /* Cluster pair-list i-cluster size atom count */
54     int  cluster_size_j;  /* Cluster pair-list j-cluster size atom count */
55 };
56
57
58 /* Add a 5% and 10% rlist buffer for simulations without dynamics (EM, NM, ...)
59  * and NVE simulations with zero initial temperature, respectively.
60  * 10% should be enough for any NVE simulation with PME and nstlist=10,
61  * for other settings it might not be enough, but then it's difficult
62  * to come up with any reasonable (not crazily expensive) value
63  * and grompp will notify the user when using the 10% buffer.
64  */
65 static const real verlet_buffer_ratio_nodynamics = 0.05;
66 static const real verlet_buffer_ratio_NVE_T0     = 0.10;
67
68
69 /* Returns the pair-list setup for the given nbnxn kernel type.
70  */
71 VerletbufListSetup verletbufGetListSetup(int nbnxnKernelType);
72
73 /* Enum for choosing the list type for verletbufGetSafeListSetup() */
74 enum class ListSetupType
75 {
76     CpuNoSimd,            /* CPU Plain-C 4x4 list */
77     CpuSimdWhenSupported, /* CPU 4xN list, where N=4 when the binary doesn't support SIMD or the smallest N supported by SIMD in this binary */
78     Gpu                   /* GPU (8x2x)8x4 list */
79 };
80
81 /* Returns the pair-list setup assumed for the current Gromacs configuration.
82  * The setup with smallest cluster sizes is returned, such that the Verlet
83  * buffer size estimated with this setup will be conservative.
84  */
85 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType);
86
87 /* Calculate the non-bonded pair-list buffer size for the Verlet list
88  * based on the particle masses, temperature, LJ types, charges
89  * and constraints as well as the non-bonded force behavior at the cut-off.
90  * The pair list update frequency and the list lifetime, which is nstlist-1
91  * for normal pair-list buffering, are passed separately, as in some cases
92  * we want an estimate for different values than the ones set in the inputrec.
93  * If reference_temperature < 0, the maximum coupling temperature will be used.
94  * The target is a maximum average energy jump per atom of
95  * ir->verletbuf_tol*nstlist*ir->delta_t over the lifetime of the list.
96  * Returns the number of non-linear virtual sites. For these it's difficult
97  * to determine their contribution to the drift exaclty, so we approximate.
98  * Returns the pair-list cut-off.
99  */
100 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
101                              const t_inputrec *ir,
102                              int               nstlist,
103                              int               list_lifetime,
104                              real reference_temperature,
105                              const VerletbufListSetup *list_setup,
106                              int *n_nonlin_vsite,
107                              real *rlist);
108
109 /* Convenience type */
110 using PartitioningPerMoltype = gmx::ArrayRef<const gmx::RangePartitioning>;
111
112 /* Determines the mininum cell size based on atom displacement
113  *
114  * The value returned is the minimum size for which the chance that
115  * an atom or update group crosses to non nearest-neighbor cells
116  * is <= chanceRequested within ir.nstlist steps.
117  * Update groups are used when !updateGrouping.empty().
118  * Without T-coupling, SD or BD, we can not estimate atom displacements
119  * and fall back to the, crude, estimate of using the pairlist buffer size.
120  *
121  * Note: Like the Verlet buffer estimate, this estimate is based on
122  *       non-interacting atoms and constrained atom-pairs. Therefore for
123  *       any system that is not an ideal gas, this will be an overestimate.
124  *
125  * Note: This size increases (very slowly) with system size.
126  */
127 real
128 minCellSizeForAtomDisplacement(const gmx_mtop_t       &mtop,
129                                const t_inputrec       &ir,
130                                PartitioningPerMoltype  updateGrouping,
131                                real                    chanceRequested);
132
133 /* Struct for unique atom type for calculating the energy drift.
134  * The atom displacement depends on mass and constraints.
135  * The energy jump for given distance depend on LJ type and q.
136  */
137 struct atom_nonbonded_kinetic_prop_t
138 {
139     real mass     = 0;     /* mass */
140     int  type     = 0;     /* type (used for LJ parameters) */
141     real q        = 0;     /* charge */
142     bool bConstr  = false; /* constrained, if TRUE, use #DOF=2 iso 3 */
143     real con_mass = 0;     /* mass of heaviest atom connected by constraints */
144     real con_len  = 0;     /* constraint length to the heaviest atom */
145 };
146
147 /* This function computes two components of the estimate of the variance
148  * in the displacement of one atom in a system of two constrained atoms.
149  * Returns in sigma2_2d the variance due to rotation of the constrained
150  * atom around the atom to which it constrained.
151  * Returns in sigma2_3d the variance due to displacement of the COM
152  * of the whole system of the two constrained atoms.
153  *
154  * Only exposed here for testing purposes.
155  */
156 void constrained_atom_sigma2(real                                 kT_fac,
157                              const atom_nonbonded_kinetic_prop_t *prop,
158                              real                                *sigma2_2d,
159                              real                                *sigma2_3d);
160
161 #endif