d4bcb4c350add584f95969cebec9eca4bf1c791f
[alexxy/gromacs.git] / src / gromacs / mdtypes / forcerec.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifndef GMX_MDTYPES_TYPES_FORCEREC_H
39 #define GMX_MDTYPES_TYPES_FORCEREC_H
40
41 #include <array>
42 #include <memory>
43 #include <vector>
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/atominfo.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/pbcutil/ishift.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/real.h"
52
53 #include "locality.h"
54
55 /* Abstract type for PME that is defined only in the routine that use them. */
56 struct gmx_pme_t;
57 struct nonbonded_verlet_t;
58 struct bonded_threading_t;
59 class DeviceContext;
60 class DispersionCorrection;
61 class ListedForces;
62 struct t_fcdata;
63 struct t_forcetable;
64 struct interaction_const_t;
65
66 namespace gmx
67 {
68 class DeviceStreamManager;
69 class ListedForcesGpu;
70 class GpuForceReduction;
71 class ForceProviders;
72 class StatePropagatorDataGpu;
73 class PmePpCommGpu;
74 class WholeMoleculeTransform;
75 } // namespace gmx
76
77 /* Value to be used in mdrun for an infinite cut-off.
78  * Since we need to compare with the cut-off squared,
79  * this value should be slighlty smaller than sqrt(GMX_FLOAT_MAX).
80  */
81 #define GMX_CUTOFF_INF 1E+18
82 //! Check the cuttoff
83 real cutoff_inf(real cutoff);
84
85 /* Forward declaration of type for managing Ewald tables */
86 struct gmx_ewald_tab_t;
87
88 struct ewald_corr_thread_t;
89
90 /*! \brief Helper force buffers for ForceOutputs
91  *
92  * This class stores intermediate force buffers that are used
93  * internally in the force calculation and which are reduced into
94  * the output force buffer passed to the force calculation.
95  */
96 class ForceHelperBuffers
97 {
98 public:
99     /*! \brief Constructs helper buffers
100      *
101      * When the forces that will be accumulated with help of these buffers
102      * have direct virial contributions, set the parameter to true, so
103      * an extra force buffer is available for these forces to enable
104      * correct virial computation.
105      */
106     ForceHelperBuffers(bool haveDirectVirialContributions);
107
108     //! Returns whether we have a direct virial contribution force buffer
109     bool haveDirectVirialContributions() const { return haveDirectVirialContributions_; }
110
111     //! Returns the buffer for direct virial contributions
112     gmx::ArrayRef<gmx::RVec> forceBufferForDirectVirialContributions()
113     {
114         GMX_ASSERT(haveDirectVirialContributions_, "Buffer can only be requested when present");
115         return forceBufferForDirectVirialContributions_;
116     }
117
118     //! Returns the buffer for shift forces, size c_numShiftVectors
119     gmx::ArrayRef<gmx::RVec> shiftForces() { return shiftForces_; }
120
121     //! Resizes the direct virial contribution buffer, when present
122     void resize(int numAtoms);
123
124 private:
125     //! True when we have contributions that are directly added to the virial
126     bool haveDirectVirialContributions_ = false;
127     //! Force buffer for force computation with direct virial contributions
128     std::vector<gmx::RVec> forceBufferForDirectVirialContributions_;
129     //! Shift force array for computing the virial, size c_numShiftVectors
130     std::vector<gmx::RVec> shiftForces_;
131 };
132 // NOLINTNEXTLINE (clang-analyzer-optin.performance.Padding)
133 struct t_forcerec
134 {
135     // Declare an explicit constructor and destructor, so they can be
136     // implemented in a single source file, so that not every source
137     // file that includes this one needs to understand how to find the
138     // destructors of the objects pointed to by unique_ptr members.
139     t_forcerec();
140     ~t_forcerec();
141
142     std::unique_ptr<interaction_const_t> ic;
143
144     /* PBC stuff */
145     PbcType pbcType = PbcType::Xyz;
146     //! Tells whether atoms inside a molecule can be in different periodic images,
147     //  i.e. whether we need to take into account PBC when computing distances inside molecules.
148     //  This determines whether PBC must be considered for e.g. bonded interactions.
149     bool            bMolPBC     = false;
150     RefCoordScaling rc_scaling  = RefCoordScaling::No;
151     gmx::RVec       posres_com  = { 0, 0, 0 };
152     gmx::RVec       posres_comB = { 0, 0, 0 };
153
154     bool use_simd_kernels = false;
155
156     /* Interaction for calculated in kernels. In many cases this is similar to
157      * the electrostatics settings in the inputrecord, but the difference is that
158      * these variables always specify the actual interaction in the kernel - if
159      * we are tabulating reaction-field the inputrec will say reaction-field, but
160      * the kernel interaction will say cubic-spline-table. To be safe we also
161      * have a kernel-specific setting for the modifiers - if the interaction is
162      * tabulated we already included the inputrec modification there, so the kernel
163      * modification setting will say 'none' in that case.
164      */
165     NbkernelElecType     nbkernel_elec_interaction = NbkernelElecType::None;
166     NbkernelVdwType      nbkernel_vdw_interaction  = NbkernelVdwType::None;
167     InteractionModifiers nbkernel_elec_modifier    = InteractionModifiers::None;
168     InteractionModifiers nbkernel_vdw_modifier     = InteractionModifiers::None;
169
170     /* Cut-Off stuff.
171      * Infinite cut-off's will be GMX_CUTOFF_INF (unlike in t_inputrec: 0).
172      */
173     real rlist = 0;
174
175     /* Charge sum for topology A/B ([0]/[1]) for Ewald corrections */
176     std::array<double, 2> qsum  = { 0 };
177     std::array<double, 2> q2sum = { 0 };
178     std::array<double, 2> c6sum = { 0 };
179
180     /* Dispersion correction stuff */
181     std::unique_ptr<DispersionCorrection> dispersionCorrection;
182
183     /* Fudge factors */
184     real fudgeQQ = 0;
185
186     std::unique_ptr<t_forcetable> pairsTable; /* for 1-4 interactions, [pairs] and [pairs_nb] */
187
188     /* Free energy */
189     FreeEnergyPerturbationType efep = FreeEnergyPerturbationType::No;
190
191     /* Information about atom properties for the molecule blocks in the global topology */
192     std::vector<gmx::AtomInfoWithinMoleculeBlock> atomInfoForEachMoleculeBlock;
193     /* Information about atom properties for local and non-local atoms */
194     std::vector<int64_t> atomInfo;
195
196     std::vector<gmx::RVec> shift_vec;
197
198     std::unique_ptr<gmx::WholeMoleculeTransform> wholeMoleculeTransform;
199
200     /* The Nbnxm Verlet non-bonded machinery */
201     std::unique_ptr<nonbonded_verlet_t> nbv;
202
203     /* The wall tables (if used) */
204     int                                                     nwall = 0;
205     std::vector<std::vector<std::unique_ptr<t_forcetable>>> wall_tab;
206
207     /* The number of atoms participating in do_force_lowlevel */
208     int natoms_force = 0;
209     /* The number of atoms participating in force calculation and constraints */
210     int natoms_force_constr = 0;
211
212     /* List of helper buffers for ForceOutputs, one for each time step with MTS */
213     std::vector<ForceHelperBuffers> forceHelperBuffers;
214
215     /* Data for PPPM/PME/Ewald */
216     gmx_pme_t*   pmedata                = nullptr;
217     LongRangeVdW ljpme_combination_rule = LongRangeVdW::Geom;
218
219     /* PME/Ewald stuff */
220     std::unique_ptr<gmx_ewald_tab_t> ewald_table;
221
222     /* Non bonded Parameter lists */
223     int               ntype          = 0; /* Number of atom types */
224     bool              haveBuckingham = false;
225     std::vector<real> nbfp;
226     std::vector<real> ljpme_c6grid; /* C6-values used on grid in LJPME */
227
228     /* Energy group pair flags */
229     int* egp_flags = nullptr;
230
231     /* Shell molecular dynamics flexible constraints */
232     real fc_stepsize = 0;
233
234     /* If > 0 signals Test Particle Insertion,
235      * the value is the number of atoms of the molecule to insert
236      * Only the energy difference due to the addition of the last molecule
237      * should be calculated.
238      */
239     int n_tpi = 0;
240
241     /* Limit for printing large forces, negative is don't print */
242     real print_force = 0;
243
244     /* User determined parameters, copied from the inputrec */
245     int  userint1  = 0;
246     int  userint2  = 0;
247     int  userint3  = 0;
248     int  userint4  = 0;
249     real userreal1 = 0;
250     real userreal2 = 0;
251     real userreal3 = 0;
252     real userreal4 = 0;
253
254     /* Data for special listed force calculations */
255     std::unique_ptr<t_fcdata> fcdata;
256
257     // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
258     std::vector<ListedForces> listedForces;
259
260     // The listed forces calculation data for GPU
261     std::unique_ptr<gmx::ListedForcesGpu> listedForcesGpu;
262
263     /* Ewald correction thread local virial and energy data */
264     int                              nthread_ewc = 0;
265     std::vector<ewald_corr_thread_t> ewc_t;
266
267     gmx::ForceProviders* forceProviders = nullptr;
268
269     // The stateGpu object is created in runner, forcerec just keeps the copy of the pointer.
270     // TODO: This is not supposed to be here. StatePropagatorDataGpu should be a part of
271     //       general StatePropagatorData object that is passed around
272     gmx::StatePropagatorDataGpu* stateGpu = nullptr;
273     // TODO: Should not be here. This is here only to pass the pointer around.
274     gmx::DeviceStreamManager* deviceStreamManager = nullptr;
275
276     //! GPU device context
277     DeviceContext* deviceContext = nullptr;
278
279     /* For PME-PP GPU communication */
280     std::unique_ptr<gmx::PmePpCommGpu> pmePpCommGpu;
281
282     /* For GPU force reduction (on both local and non-local atoms) */
283     gmx::EnumerationArray<gmx::AtomLocality, std::unique_ptr<gmx::GpuForceReduction>> gpuForceReduction;
284 };
285
286 /* Important: Starting with Gromacs-4.6, the values of c6 and c12 in the nbfp array have
287  * been scaled by 6.0 or 12.0 to save flops in the kernels. We have corrected this everywhere
288  * in the code, but beware if you are using these macros externally.
289  */
290 #define C6(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj))]
291 #define C12(nbfp, ntp, ai, aj) (nbfp)[2 * ((ntp) * (ai) + (aj)) + 1]
292 #define BHAMC(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj))]
293 #define BHAMA(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 1]
294 #define BHAMB(nbfp, ntp, ai, aj) (nbfp)[3 * ((ntp) * (ai) + (aj)) + 2]
295
296 #endif