65287fcb14827696d493f92c5b7c7b9c7a3be7f3
[alexxy/gromacs.git] / src / gromacs / mdlib / vsite.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 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_VSITE_H
39 #define GMX_MDLIB_VSITE_H
40
41 #include <memory>
42
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/pbcutil/ishift.h"
45 #include "gromacs/topology/idef.h"
46 #include "gromacs/utility/arrayref.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/real.h"
49
50 struct gmx_localtop_t;
51 struct gmx_mtop_t;
52 struct t_commrec;
53 struct t_graph;
54 struct t_ilist;
55 struct t_mdatoms;
56 struct t_nrnb;
57 struct gmx_wallcycle;
58 struct VsiteThread;
59 enum class PbcType : int;
60
61 namespace gmx
62 {
63 class RangePartitioning;
64 }
65
66 /* The start and end values of for the vsite indices in the ftype enum.
67  * The validity of these values is checked in init_vsite.
68  * This is used to avoid loops over all ftypes just to get the vsite entries.
69  * (We should replace the fixed ilist array by only the used entries.)
70  */
71 static constexpr int c_ftypeVsiteStart = F_VSITE2;
72 static constexpr int c_ftypeVsiteEnd   = F_VSITEN + 1;
73
74 /* Type for storing PBC atom information for all vsite types in the system */
75 typedef std::array<std::vector<int>, c_ftypeVsiteEnd - c_ftypeVsiteStart> VsitePbc;
76
77 /* Data for handling vsites, needed with OpenMP threading or with charge-groups and PBC */
78 struct gmx_vsite_t
79 {
80     gmx_vsite_t();
81
82     ~gmx_vsite_t();
83
84     /* The number of vsites that cross update groups, when =0 no PBC treatment is needed */
85     int numInterUpdategroupVsites;
86     int nthreads;                                    /* Number of threads used for vsites       */
87     std::vector<std::unique_ptr<VsiteThread>> tData; /* Thread local vsites and work structs    */
88     std::vector<int> taskIndex;                      /* Work array                              */
89     bool useDomdec; /* Tells whether we use domain decomposition with more than 1 DD rank */
90 };
91
92 /*! \brief Create positions of vsite atoms based for the local system
93  *
94  * \param[in]     vsite    The vsite struct, when nullptr is passed, no MPI and no multi-threading
95  *                         is used
96  * \param[in,out] x        The coordinates
97  * \param[in]     dt       The time step
98  * \param[in,out] v        When != nullptr, velocities for vsites are set as displacement/dt
99  * \param[in]     ip       Interaction parameters
100  * \param[in]     ilist    The interaction list
101  * \param[in]     pbcType  The type of periodic boundary conditions
102  * \param[in]     bMolPBC  When true, molecules are broken over PBC
103  * \param[in]     cr       The communication record
104  * \param[in]     box      The box
105  */
106 void construct_vsites(const gmx_vsite_t* vsite,
107                       rvec               x[],
108                       real               dt,
109                       rvec               v[],
110                       const t_iparams    ip[],
111                       const t_ilist      ilist[],
112                       PbcType            pbcType,
113                       gmx_bool           bMolPBC,
114                       const t_commrec*   cr,
115                       const matrix       box);
116
117 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
118  *
119  * \param[in]     mtop  The global topology
120  * \param[in,out] x     The global coordinates
121  */
122 void constructVsitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x);
123
124 void spread_vsite_f(const gmx_vsite_t* vsite,
125                     const rvec         x[],
126                     rvec               f[],
127                     rvec*              fshift,
128                     gmx_bool           VirCorr,
129                     matrix             vir,
130                     t_nrnb*            nrnb,
131                     const t_idef*      idef,
132                     PbcType            pbcType,
133                     gmx_bool           bMolPBC,
134                     const t_graph*     g,
135                     const matrix       box,
136                     const t_commrec*   cr,
137                     gmx_wallcycle*     wcycle);
138 /* Spread the force operating on the vsite atoms on the surrounding atoms.
139  * If fshift!=NULL also update the shift forces.
140  * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
141  * to vir. This correction is required when the virial is not calculated
142  * afterwards from the particle position and forces, but in a different way,
143  * as for instance for the PME mesh contribution.
144  */
145
146 /* Return the number of non-linear virtual site constructions in the system */
147 int countNonlinearVsites(const gmx_mtop_t& mtop);
148
149 /* Return the number of virtual sites that cross update groups
150  *
151  * \param[in] mtop                           The global topology
152  * \param[in] updateGroupingPerMoleculetype  Update grouping per molecule type, pass empty when not using update groups
153  */
154 int countInterUpdategroupVsites(const gmx_mtop_t&                           mtop,
155                                 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype);
156
157 /* Initialize the virtual site struct,
158  *
159  * \param[in] mtop  The global topology
160  * \param[in] cr    The communication record
161  * \returns A valid vsite struct or nullptr when there are no virtual sites
162  */
163 std::unique_ptr<gmx_vsite_t> initVsite(const gmx_mtop_t& mtop, const t_commrec* cr);
164
165 void split_vsites_over_threads(const t_ilist*   ilist,
166                                const t_iparams* ip,
167                                const t_mdatoms* mdatoms,
168                                gmx_vsite_t*     vsite);
169 /* Divide the vsite work-load over the threads.
170  * Should be called at the end of the domain decomposition.
171  */
172
173 void set_vsite_top(gmx_vsite_t* vsite, const gmx_localtop_t* top, const t_mdatoms* md);
174 /* Set some vsite data for runs without domain decomposition.
175  * Should be called once after init_vsite, before calling other routines.
176  */
177
178 #endif