Update instructions in containers.rst
[gromacs.git] / src / gromacs / domdec / domdec_specatomcomm.h
blobab7d723c819ed7e8de4a4fc3d28f3162bf36be54
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005,2006,2007,2008,2009 by the GROMACS development team.
5 * Copyright (c) 2010,2012,2013,2014,2015 by the GROMACS development team.
6 * Copyright (c) 2017,2018,2019,2020, 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.
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.
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.
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.
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.
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.
37 /*! \internal \file
39 * \brief This file declares functions for domdec to use
40 * while managing communication of atoms required for special purposes
42 * \author Berk Hess <hess@kth.se>
43 * \ingroup module_domdec
46 #ifndef GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
47 #define GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
49 #include <vector>
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/utility/basedefinitions.h"
54 struct gmx_domdec_t;
55 struct t_commrec;
57 namespace gmx
59 template<typename T>
60 class HashedMap;
61 } // namespace gmx
63 /*! \internal \brief The communication setup along a single dimension */
64 struct gmx_specatsend_t
66 std::vector<int> a; /**< The indices of atoms to send */
67 int nrecv; /**< The number of atoms to receive */
70 /*! \internal \brief Struct with setup and buffers for special atom communication */
71 struct gmx_domdec_specat_comm_t
73 /* The number of indices to receive during the setup */
74 int nreq[DIM][2][2] = { { { 0 } } }; /**< The nr. of atoms requested, per DIM, direction and direct/indirect */
75 /* The atoms to send */
76 gmx_specatsend_t spas[DIM][2]; /**< The communication setup per DIM, direction */
77 std::vector<bool> sendAtom; /**< Work buffer that tells if spec.atoms should be sent */
79 /* Send buffers */
80 std::vector<int> ibuf; /**< Integer send buffer */
81 std::vector<gmx::RVec> vbuf; /**< rvec send buffer */
82 std::vector<gmx::RVec> vbuf2; /**< rvec send buffer */
83 /* The range in the local buffer(s) for received atoms */
84 int at_start; /**< Start index of received atoms */
85 int at_end; /**< End index of received atoms */
88 /*! \brief Communicates the force for special atoms, the shift forces are reduced with \p fshift != NULL */
89 void dd_move_f_specat(const gmx_domdec_t* dd, gmx_domdec_specat_comm_t* spac, rvec* f, rvec* fshift);
91 /*! \brief Communicates the coordinates for special atoms
93 * \param[in] dd Domain decomposition struct
94 * \param[in] spac Special atom communication struct
95 * \param[in] box Box, used for pbc
96 * \param[in,out] x0 Vector to communicate
97 * \param[in,out] x1 Vector to communicate, when != NULL
98 * \param[in] bX1IsCoord Tells is \p x1 is a coordinate vector (needs pbc)
100 void dd_move_x_specat(const gmx_domdec_t* dd,
101 gmx_domdec_specat_comm_t* spac,
102 const matrix box,
103 rvec* x0,
104 rvec* x1,
105 gmx_bool bX1IsCoord);
107 /*! \brief Sets up the communication for special atoms
109 * \param[in] dd Domain decomposition struct
110 * \param[in,out] ireq List of requested atom indices, updated due to aggregation
111 * \param[in,out] spac Special atom communication struct
112 * \param[out] ga2la_specat Global to local special atom index
113 * \param[in] at_start Index in local state where to start storing communicated atoms
114 * \param[in] vbuf_fac Buffer factor, 1 or 2 for communicating 1 or 2 vectors
115 * \param[in] specat_type Name of the special atom, used for error message
116 * \param[in] add_err Text to add at the end of error message when atoms can't be found
118 int setup_specat_communication(gmx_domdec_t* dd,
119 std::vector<int>* ireq,
120 gmx_domdec_specat_comm_t* spac,
121 gmx::HashedMap<int>* ga2la_specat,
122 int at_start,
123 int vbuf_fac,
124 const char* specat_type,
125 const char* add_err);
127 #endif