b82751206a59532580233c9675c733033015e33f
[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, 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 /* The start and end values of for the vsite indices in the ftype enum.
60  * The validity of these values is checked in init_vsite.
61  * This is used to avoid loops over all ftypes just to get the vsite entries.
62  * (We should replace the fixed ilist array by only the used entries.)
63  */
64 static constexpr int c_ftypeVsiteStart = F_VSITE2;
65 static constexpr int c_ftypeVsiteEnd   = F_VSITEN + 1;
66
67 /* Type for storing PBC atom information for all vsite types in the system */
68 typedef std::array<std::vector<int>, c_ftypeVsiteEnd - c_ftypeVsiteStart> VsitePbc;
69
70 /* Data for handling vsites, needed with OpenMP threading or with charge-groups and PBC */
71 struct gmx_vsite_t
72 {
73     gmx_vsite_t();
74
75     ~gmx_vsite_t();
76
77     gmx_bool                  bHaveChargeGroups;         /* Do we have charge groups?               */
78     int                       n_intercg_vsite;           /* The number of inter charge group vsites */
79     std::vector<VsitePbc>     vsite_pbc_molt;            /* The pbc atoms for intercg vsites        */
80     std::unique_ptr<VsitePbc> vsite_pbc_loc;             /* The local pbc atoms                     */
81     int                       nthreads;                  /* Number of threads used for vsites       */
82     std::vector < std::unique_ptr < VsiteThread>> tData; /* Thread local vsites and work structs    */
83     std::vector<int>          taskIndex;                 /* Work array                              */
84     bool                      useDomdec;                 /* Tells whether we use domain decomposition with more than 1 DD rank */
85 };
86
87 /*! \brief Create positions of vsite atoms based for the local system
88  *
89  * \param[in]     vsite    The vsite struct, when nullptr is passed, no MPI and no multi-threading is used
90  * \param[in,out] x        The coordinates
91  * \param[in]     dt       The time step
92  * \param[in,out] v        When != nullptr, velocities for vsites are set as displacement/dt
93  * \param[in]     ip       Interaction parameters
94  * \param[in]     ilist    The interaction list
95  * \param[in]     ePBC     The type of periodic boundary conditions
96  * \param[in]     bMolPBC  When true, molecules are broken over PBC
97  * \param[in]     cr       The communication record
98  * \param[in]     box      The box
99  */
100 void construct_vsites(const gmx_vsite_t *vsite,
101                       rvec x[],
102                       real dt, rvec v[],
103                       const t_iparams ip[], const t_ilist ilist[],
104                       int ePBC, gmx_bool bMolPBC,
105                       const t_commrec *cr,
106                       const matrix box);
107
108 /*! \brief Create positions of vsite atoms for the whole system assuming all molecules are wholex
109  *
110  * \param[in]     mtop  The global topology
111  * \param[in,out] x     The global coordinates
112  */
113 void constructVsitesGlobal(const gmx_mtop_t         &mtop,
114                            gmx::ArrayRef<gmx::RVec>  x);
115
116 void spread_vsite_f(const gmx_vsite_t *vsite,
117                     const rvec x[],
118                     rvec f[], rvec *fshift,
119                     gmx_bool VirCorr, matrix vir,
120                     t_nrnb *nrnb, const t_idef *idef,
121                     int ePBC, gmx_bool bMolPBC, const t_graph *g, const matrix box,
122                     const t_commrec *cr, gmx_wallcycle *wcycle);
123 /* Spread the force operating on the vsite atoms on the surrounding atoms.
124  * If fshift!=NULL also update the shift forces.
125  * If VirCorr=TRUE add the virial correction for non-linear vsite constructs
126  * to vir. This correction is required when the virial is not calculated
127  * afterwards from the particle position and forces, but in a different way,
128  * as for instance for the PME mesh contribution.
129  */
130
131 int count_intercg_vsites(const gmx_mtop_t *mtop);
132 /* Returns the number of virtual sites that cross charge groups */
133
134 /* Initialize the virtual site struct,
135  *
136  * \param[in] mtop  The global topology
137  * \param[in] cr    The communication record
138  * \returns A valid vsite struct or nullptr when there are no virtual sites
139  */
140 std::unique_ptr<gmx_vsite_t>
141 initVsite(const gmx_mtop_t &mtop,
142           const t_commrec  *cr);
143
144 void split_vsites_over_threads(const t_ilist   *ilist,
145                                const t_iparams *ip,
146                                const t_mdatoms *mdatoms,
147                                gmx_vsite_t     *vsite);
148 /* Divide the vsite work-load over the threads.
149  * Should be called at the end of the domain decomposition.
150  */
151
152 void set_vsite_top(gmx_vsite_t          *vsite,
153                    const gmx_localtop_t *top,
154                    const t_mdatoms      *md);
155 /* Set some vsite data for runs without domain decomposition.
156  * Should be called once after init_vsite, before calling other routines.
157  */
158
159 #endif