2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 #ifndef GMX_MDLIB_FORCE_H
39 #define GMX_MDLIB_FORCE_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/arrayref.h"
48 #include "gromacs/utility/enumerationhelpers.h"
50 class DDBalanceRegionHandler;
52 struct gmx_enerdata_t;
54 struct SimulationGroups;
55 struct gmx_localtop_t;
56 struct gmx_multisim_t;
60 class InteractionDefinitions;
68 struct gmx_ewald_tab_t;
69 class CpuPpLongRangeNonbondeds;
74 class ArrayRefWithPadding;
76 class ForceBuffersView;
77 class ForceWithVirial;
79 class MdrunScheduleWorkload;
82 class VirtualSitesHandler;
85 struct ewald_corr_thread_t
89 gmx::EnumerationArray<FreeEnergyPerturbationCouplingType, real> dvdl;
94 /* Perform the force and, if requested, energy computation
96 * Without multiple time stepping the force is returned in force->force().
98 * With multiple time stepping the behavior depends on the integration step.
99 * At fast steps (step % mtsFactor != 0), the fast force is returned in
100 * force->force(). The force->forceMtsCombined() buffer is unused.
101 * At slow steps, the normal force is returned in force->force(),
102 * unless the \p GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE is set in \p legacyFlags.
103 * A MTS-combined force, F_fast + mtsFactor*F_slow, is always returned in
104 * force->forceMtsCombined(). This forceMts can be used directly in a standard
105 * leap-frog integrator to do multiple time stepping.
107 void do_force(FILE* log,
109 const gmx_multisim_t* ms,
110 const t_inputrec& inputrec,
112 gmx_enfrot* enforcedRotation,
113 gmx::ImdSession* imdSession,
117 gmx_wallcycle* wcycle,
118 const gmx_localtop_t* top,
120 gmx::ArrayRefWithPadding<gmx::RVec> coordinates,
121 const history_t* hist,
122 gmx::ForceBuffersView* force,
124 const t_mdatoms* mdatoms,
125 gmx_enerdata_t* enerd,
126 gmx::ArrayRef<const real> lambda,
128 gmx::MdrunScheduleWorkload* runScheduleWork,
129 gmx::VirtualSitesHandler* vsite,
133 CpuPpLongRangeNonbondeds* longRangeNonbondeds,
135 const DDBalanceRegionHandler& ddBalanceRegionHandler);
137 /* Communicate coordinates (if parallel).
138 * Do neighbor searching (if necessary).
140 * Communicate forces (if parallel).
141 * Spread forces for vsites (if present).
143 * f is always required.
147 class CpuPpLongRangeNonbondeds
150 /* \brief Constructor
152 * Should be called after init_forcerec if params come from a populated forcerec
154 CpuPpLongRangeNonbondeds(int numberOfTestPaticles,
157 gmx::ArrayRef<const double> chargeC6Sum,
158 CoulombInteractionType eeltype,
159 VanDerWaalsType vdwtype,
160 const t_inputrec& inputrec,
162 gmx_wallcycle* wcycle,
165 ~CpuPpLongRangeNonbondeds();
167 void updateAfterPartition(const t_mdatoms& md);
169 /* Calculate CPU Ewald or PME-mesh forces when done on this rank and Ewald corrections, when used
171 * Note that Ewald dipole and net charge corrections are always computed here, independently
172 * of whether the PME-mesh contribution is computed on a separate PME rank or on a GPU.
174 void calculate(gmx_pme_t* pmedata,
175 const t_commrec* commrec,
176 gmx::ArrayRef<const gmx::RVec> coordinates,
177 gmx::ForceWithVirial* forceWithVirial,
178 gmx_enerdata_t* enerd,
180 gmx::ArrayRef<const real> lambda,
181 gmx::ArrayRef<const gmx::RVec> mu_tot,
182 const gmx::StepWorkload& stepWork,
183 const DDBalanceRegionHandler& ddBalanceRegionHandler);
186 //! Number of particles for test particle insertion
188 //! Ewald charge coefficient
190 //! Dielectric constant
192 //! [0]: sum of charges; [1]: sum of C6's
193 gmx::ArrayRef<const double> chargeC6Sum_;
194 //! Cut-off treatment for Coulomb
195 CoulombInteractionType coulombInteractionType_;
196 //! Van der Waals interaction treatment
197 VanDerWaalsType vanDerWaalsType_;
199 EwaldGeometry ewaldGeometry_;
200 //! Epsilon for PME dipole correction
201 real epsilonSurface_;
202 //! Whether a long range correction is used
203 bool haveEwaldSurfaceTerm_;
204 //! Scaling factor for the box for Ewald
206 //! Whether the simulation is 2D periodic with two walls
207 bool havePbcXY2Walls_;
208 //! Free energy perturbation type
209 FreeEnergyPerturbationType freeEnergyPerturbationType_;
210 //! Number of atoms on this node
212 //! Whether there are perturbed interactions
215 gmx::ArrayRef<const real> chargeA_;
217 gmx::ArrayRef<const real> chargeB_;
219 gmx::ArrayRef<const real> sqrt_c6A_;
221 gmx::ArrayRef<const real> sqrt_c6B_;
223 gmx::ArrayRef<const real> sigmaA_;
225 gmx::ArrayRef<const real> sigmaB_;
226 //! Ewald correction thread local virial and energy data
227 std::vector<ewald_corr_thread_t> outputPerThread_;
229 std::unique_ptr<gmx_ewald_tab_t> ewaldTable_;
230 //! Non bonded kernel flop counters
232 //! Wall cycle counters
233 gmx_wallcycle* wcycle_;