c8cfb79c161814627f4f164b890e606ffd9f147b
[alexxy/gromacs.git] / src / gromacs / mdlib / force.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, 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_MDLIB_FORCE_H
39 #define GMX_MDLIB_FORCE_H
40
41 #include <cstdio>
42
43 #include "gromacs/math/vectypes.h"
44
45 class DDBalanceRegionHandler;
46 struct gmx_edsam;
47 struct gmx_enerdata_t;
48 struct gmx_enfrot;
49 struct SimulationGroups;
50 struct gmx_localtop_t;
51 struct gmx_multisim_t;
52 struct gmx_wallcycle;
53 class history_t;
54 class InteractionDefinitions;
55 struct pull_t;
56 struct t_commrec;
57 struct t_forcerec;
58 struct t_inputrec;
59 struct t_lambda;
60 struct t_mdatoms;
61 struct t_nrnb;
62
63 namespace gmx
64 {
65 template<typename>
66 class ArrayRef;
67 template<typename>
68 class ArrayRefWithPadding;
69 class Awh;
70 class ForceBuffersView;
71 class ForceWithVirial;
72 class ImdSession;
73 class MdrunScheduleWorkload;
74 class MDLogger;
75 class StepWorkload;
76 class VirtualSitesHandler;
77 } // namespace gmx
78
79 void do_force(FILE*                               log,
80               const t_commrec*                    cr,
81               const gmx_multisim_t*               ms,
82               const t_inputrec*                   inputrec,
83               gmx::Awh*                           awh,
84               gmx_enfrot*                         enforcedRotation,
85               gmx::ImdSession*                    imdSession,
86               pull_t*                             pull_work,
87               int64_t                             step,
88               t_nrnb*                             nrnb,
89               gmx_wallcycle*                      wcycle,
90               const gmx_localtop_t*               top,
91               const matrix                        box,
92               gmx::ArrayRefWithPadding<gmx::RVec> coordinates,
93               history_t*                          hist,
94               gmx::ForceBuffersView*              force,
95               tensor                              vir_force,
96               const t_mdatoms*                    mdatoms,
97               gmx_enerdata_t*                     enerd,
98               gmx::ArrayRef<const real>           lambda,
99               t_forcerec*                         fr,
100               gmx::MdrunScheduleWorkload*         runScheduleWork,
101               gmx::VirtualSitesHandler*           vsite,
102               rvec                                mu_tot,
103               double                              t,
104               gmx_edsam*                          ed,
105               int                                 legacyFlags,
106               const DDBalanceRegionHandler&       ddBalanceRegionHandler);
107
108 /* Communicate coordinates (if parallel).
109  * Do neighbor searching (if necessary).
110  * Calculate forces.
111  * Communicate forces (if parallel).
112  * Spread forces for vsites (if present).
113  *
114  * f is always required.
115  */
116
117
118 /* Calculate CPU Ewald or PME-mesh forces when done on this rank and Ewald corrections, when used
119  *
120  * Note that Ewald dipole and net charge corrections are always computed here, independently
121  * on whether the PME-mesh contribution is computed on a separate PME rank or on a GPU.
122  */
123 void calculateLongRangeNonbondeds(t_forcerec*                    fr,
124                                   const t_inputrec*              ir,
125                                   const t_commrec*               cr,
126                                   t_nrnb*                        nrnb,
127                                   gmx_wallcycle*                 wcycle,
128                                   const t_mdatoms*               md,
129                                   gmx::ArrayRef<const gmx::RVec> coordinates,
130                                   gmx::ForceWithVirial*          forceWithVirial,
131                                   gmx_enerdata_t*                enerd,
132                                   const matrix                   box,
133                                   const real*                    lambda,
134                                   const rvec*                    mu_tot,
135                                   const gmx::StepWorkload&       stepWork,
136                                   const DDBalanceRegionHandler&  ddBalanceRegionHandler);
137
138 #endif