b224008d1b39960797c9101137c7c8b9ceee0637
[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,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #ifndef GMX_MDLIB_VSITE_H
38 #define GMX_MDLIB_VSITE_H
39
40 #include <memory>
41
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/pbcutil/ishift.h"
44 #include "gromacs/topology/idef.h"
45 #include "gromacs/utility/arrayref.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/real.h"
48
49 struct gmx_localtop_t;
50 struct gmx_mtop_t;
51 struct t_commrec;
52 struct t_graph;
53 struct t_ilist;
54 struct t_mdatoms;
55 struct t_nrnb;
56 struct gmx_wallcycle;
57 struct VsiteThread;
58
59 namespace gmx
60 {
61 class RangePartitioning;
62 }
63
64 /* The start and end values of for the vsite indices in the ftype enum.
65  * The validity of these values is checked in init_vsite.
66  * This is used to avoid loops over all ftypes just to get the vsite entries.
67  * (We should replace the fixed ilist array by only the used entries.)
68  */
69 static constexpr int c_ftypeVsiteStart = F_VSITE2;
70 static constexpr int c_ftypeVsiteEnd   = F_VSITEN + 1;
71
72 /* Type for storing PBC atom information for all vsite types in the system */
73 typedef std::array<std::vector<int>, c_ftypeVsiteEnd - c_ftypeVsiteStart> VsitePbc;
74
75 /* Data for handling vsites, needed with OpenMP threading or with charge-groups and PBC */
76 struct gmx_vsite_t
77 {
78     gmx_vsite_t();
79
80     ~gmx_vsite_t();
81
82     /* The number of vsites that cross update groups, when =0 no PBC treatment is needed */
83     int numInterUpdategroupVsites;
84     int nthreads;                                    /* Number of threads used for vsites       */
85     std::vector<std::unique_ptr<VsiteThread>> tData; /* Thread local vsites and work structs    */
86     std::vector<int> taskIndex;                      /* Work array                              */
87     bool useDomdec; /* Tells whether we use domain decomposition with more than 1 DD rank */
88 };
89
90 /*! \brief Create positions of vsite atoms based for the local system
91  *
92  * \param[in]     vsite    The vsite struct, when nullptr is passed, no MPI and no multi-threading
93  * is used \param[in,out] x        The coordinates \param[in]     dt       The time step \param[in,out]
94  * v        When != nullptr, velocities for vsites are set as displacement/dt \param[in]     ip
95  * Interaction parameters \param[in]     ilist    The interaction list \param[in]     ePBC     The
96  * type of periodic boundary conditions \param[in]     bMolPBC  When true, molecules are broken over
97  * PBC \param[in]     cr       The communication record \param[in]     box      The box
98  */
99 void construct_vsites(const gmx_vsite_t* vsite,
100                       rvec               x[],
101                       real               dt,
102                       rvec               v[],
103                       const t_iparams    ip[],
104                       const t_ilist      ilist[],
105                       int                ePBC,
106                       gmx_bool           bMolPBC,
107                       const t_commrec*   cr,
108                       const matrix       box);
109
110 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
111  *
112  * \param[in]     mtop  The global topology
113  * \param[in,out] x     The global coordinates
114  */
115 void constructVsitesGlobal(const gmx_mtop_t& mtop, gmx::ArrayRef<gmx::RVec> x);
116
117 void spread_vsite_f(const gmx_vsite_t* vsite,
118                     const rvec         x[],
119                     rvec               f[],
120                     rvec*              fshift,
121                     gmx_bool           VirCorr,
122                     matrix             vir,
123                     t_nrnb*            nrnb,
124                     const t_idef*      idef,
125                     int                ePBC,
126                     gmx_bool           bMolPBC,
127                     const t_graph*     g,
128                     const matrix       box,
129                     const t_commrec*   cr,
130                     gmx_wallcycle*     wcycle);
131 /* Spread the force operating on the vsite atoms on the surrounding atoms.
132  * If fshift!=NULL also update the shift forces.
133  * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
134  * to vir. This correction is required when the virial is not calculated
135  * afterwards from the particle position and forces, but in a different way,
136  * as for instance for the PME mesh contribution.
137  */
138
139 /* Return the number of non-linear virtual site constructions in the system */
140 int countNonlinearVsites(const gmx_mtop_t& mtop);
141
142 /* Return the number of virtual sites that cross update groups
143  *
144  * \param[in] mtop                           The global topology
145  * \param[in] updateGroupingPerMoleculetype  Update grouping per molecule type, pass empty when not using update groups
146  */
147 int countInterUpdategroupVsites(const gmx_mtop_t&                           mtop,
148                                 gmx::ArrayRef<const gmx::RangePartitioning> updateGroupingPerMoleculetype);
149
150 /* Initialize the virtual site struct,
151  *
152  * \param[in] mtop  The global topology
153  * \param[in] cr    The communication record
154  * \returns A valid vsite struct or nullptr when there are no virtual sites
155  */
156 std::unique_ptr<gmx_vsite_t> initVsite(const gmx_mtop_t& mtop, const t_commrec* cr);
157
158 void split_vsites_over_threads(const t_ilist*   ilist,
159                                const t_iparams* ip,
160                                const t_mdatoms* mdatoms,
161                                gmx_vsite_t*     vsite);
162 /* Divide the vsite work-load over the threads.
163  * Should be called at the end of the domain decomposition.
164  */
165
166 void set_vsite_top(gmx_vsite_t* vsite, const gmx_localtop_t* top, const t_mdatoms* md);
167 /* Set some vsite data for runs without domain decomposition.
168  * Should be called once after init_vsite, before calling other routines.
169  */
170
171 #endif