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