a051190d754c7f3065934c066d18e63bd8bb5e27
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_specatomcomm.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2005,2006,2007,2008,2009,2010,2012,2013,2014,2015,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *
37  * \brief This file declares functions for domdec to use
38  * while managing communication of atoms required for special purposes
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_domdec
42  */
43
44 #ifndef GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
45 #define GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
46
47 #include <vector>
48
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/utility/basedefinitions.h"
51
52 struct gmx_domdec_t;
53 struct t_commrec;
54
55 namespace gmx
56 {
57 template<typename T>
58 class HashedMap;
59 } // namespace gmx
60
61 /*! \internal \brief The communication setup along a single dimension */
62 struct gmx_specatsend_t
63 {
64     std::vector<int> a;     /**< The indices of atoms to send */
65     int              nrecv; /**< The number of atoms to receive */
66 };
67
68 /*! \internal \brief Struct with setup and buffers for special atom communication */
69 struct gmx_domdec_specat_comm_t
70 {
71     /* The number of indices to receive during the setup */
72     int nreq[DIM][2][2] = { { { 0 } } }; /**< The nr. of atoms requested, per DIM, direction and direct/indirect */
73     /* The atoms to send */
74     gmx_specatsend_t  spas[DIM][2]; /**< The communication setup per DIM, direction */
75     std::vector<bool> sendAtom;     /**< Work buffer that tells if spec.atoms should be sent */
76
77     /* Send buffers */
78     std::vector<int>       ibuf;  /**< Integer send buffer */
79     std::vector<gmx::RVec> vbuf;  /**< rvec send buffer */
80     std::vector<gmx::RVec> vbuf2; /**< rvec send buffer */
81     /* The range in the local buffer(s) for received atoms */
82     int at_start; /**< Start index of received atoms */
83     int at_end;   /**< End index of received atoms */
84 };
85
86 /*! \brief Communicates the force for special atoms, the shift forces are reduced with \p fshift != NULL */
87 void dd_move_f_specat(gmx_domdec_t* dd, gmx_domdec_specat_comm_t* spac, rvec* f, rvec* fshift);
88
89 /*! \brief Communicates the coordinates for special atoms
90  *
91  * \param[in]     dd         Domain decomposition struct
92  * \param[in]     spac       Special atom communication struct
93  * \param[in]     box        Box, used for pbc
94  * \param[in,out] x0         Vector to communicate
95  * \param[in,out] x1         Vector to communicate, when != NULL
96  * \param[in]     bX1IsCoord Tells is \p x1 is a coordinate vector (needs pbc)
97  */
98 void dd_move_x_specat(gmx_domdec_t*             dd,
99                       gmx_domdec_specat_comm_t* spac,
100                       const matrix              box,
101                       rvec*                     x0,
102                       rvec*                     x1,
103                       gmx_bool                  bX1IsCoord);
104
105 /*! \brief Sets up the communication for special atoms
106  *
107  * \param[in]     dd           Domain decomposition struct
108  * \param[in,out] ireq         List of requested atom indices, updated due to aggregation
109  * \param[in,out] spac         Special atom communication struct
110  * \param[out]    ga2la_specat Global to local special atom index
111  * \param[in]     at_start     Index in local state where to start storing communicated atoms
112  * \param[in]     vbuf_fac     Buffer factor, 1 or 2 for communicating 1 or 2 vectors
113  * \param[in]     specat_type  Name of the special atom, used for error message
114  * \param[in]     add_err      Text to add at the end of error message when atoms can't be found
115  */
116 int setup_specat_communication(gmx_domdec_t*             dd,
117                                std::vector<int>*         ireq,
118                                gmx_domdec_specat_comm_t* spac,
119                                gmx::HashedMap<int>*      ga2la_specat,
120                                int                       at_start,
121                                int                       vbuf_fac,
122                                const char*               specat_type,
123                                const char*               add_err);
124
125 #endif