Make domdec a proper module
[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, 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 "gromacs/legacyheaders/types/commrec.h"
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/utility/basedefinitions.h"
50
51 typedef struct {
52     int  nsend;
53     int *a;
54     int  a_nalloc;
55     int  nrecv;
56 } gmx_specatsend_t;
57
58 typedef struct {
59     int *ind;
60     int  nalloc;
61     int  n;
62 } ind_req_t;
63
64 /*! \internal \brief Struct with setup and buffers for special atom communication */
65 typedef struct gmx_domdec_specat_comm {
66     /* The number of indices to receive during the setup */
67     int              nreq[DIM][2][2];  /**< The nr. of atoms requested, per DIM, direction and direct/indirect */
68     /* The atoms to send */
69     gmx_specatsend_t spas[DIM][2];     /**< The communication setup per DIM, direction */
70     gmx_bool        *bSendAtom;        /**< Work buffer that tells if spec.atoms should be sent */
71     int              bSendAtom_nalloc; /**< Allocation size of \p bSendAtom */
72     /* Send buffers */
73     int             *ibuf;             /**< Integer send buffer */
74     int              ibuf_nalloc;      /**< Allocation size of \p ibuf */
75     rvec            *vbuf;             /**< rvec send buffer */
76     int              vbuf_nalloc;      /**< Allocation size of \p vbuf */
77     rvec            *vbuf2;            /**< rvec send buffer */
78     int              vbuf2_nalloc;     /**< Allocation size of \p vbuf2 */
79     /* The range in the local buffer(s) for received atoms */
80     int              at_start;         /**< Start index of received atoms */
81     int              at_end;           /**< End index of received atoms */
82
83     /* The atom indices we need from the surrounding cells.
84      * We can gather the indices over nthread threads.
85      */
86     int        nthread;                /**< Number of threads used for spec.atom communication */
87     ind_req_t *ireq;                   /**< Index request buffer per thread, allocation size \p nthread */
88 } gmx_domdec_specat_comm_t;
89
90 /*! \brief Communicates the force for special atoms, the shift forces are reduced with \p fshift != NULL */
91 void dd_move_f_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
92                       rvec *f, rvec *fshift);
93
94 /*! \brief Communicates the coordinates for special atoms
95  *
96  * \param[in]     dd         Domain decomposition struct
97  * \param[in]     spac       Special atom communication struct
98  * \param[in]     box        Box, used for pbc
99  * \param[in,out] x0         Vector to communicate
100  * \param[in,out] x1         Vector to communicate, when != NULL
101  * \param[in]     bX1IsCoord Tells is \p x1 is a coordinate vector (needs pbc)
102  */
103 void dd_move_x_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
104                       matrix box,
105                       rvec *x0,
106                       rvec *x1, gmx_bool bX1IsCoord);
107
108 /*! \brief Sets up the communication for special atoms
109  *
110  * \param[in]     dd         Domain decomposition struct
111  * \param[in]     ireq List of requested atom indices
112  * \param[in,out] spac   Special atom communication struct
113  * \param[out]    ga2la_specat Global to local special atom index
114  * \param[in]     at_start     Index in local state where to start storing communicated atoms
115  * \param[in]     vbuf_fac     Buffer factor, 1 or 2 for communicating 1 or 2 vectors
116  * \param[in]     specat_type  Name of the special atom, used for error message
117  * \param[in]     add_err      Text to add at the end of error message when atoms can't be found
118  */
119 int setup_specat_communication(gmx_domdec_t             *dd,
120                                ind_req_t                *ireq,
121                                gmx_domdec_specat_comm_t *spac,
122                                gmx_hash_t                ga2la_specat,
123                                int                       at_start,
124                                int                       vbuf_fac,
125                                const char               *specat_type,
126                                const char               *add_err);
127
128 /*! \brief Initialize a special communication struct */
129 gmx_domdec_specat_comm_t *specat_comm_init(int nthread);
130
131 #endif