Work around memory leak in t_state
[gromacs.git] / src / gromacs / domdec / domdec_specatomcomm.h
blobb08b581a88dfb039d5f6bdad1b1002316539c424
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2005,2006,2007,2008,2009,2010,2012,2013,2014,2015, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
37 * \brief This file declares functions for domdec to use
38 * while managing communication of atoms required for special purposes
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_domdec
44 #ifndef GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
45 #define GMX_DOMDEC_DOMDEC_SPECATOMCOMM_H
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/utility/basedefinitions.h"
51 struct gmx_hash_t;
53 typedef struct {
54 int nsend;
55 int *a;
56 int a_nalloc;
57 int nrecv;
58 } gmx_specatsend_t;
60 typedef struct {
61 int *ind;
62 int nalloc;
63 int n;
64 } ind_req_t;
66 /*! \internal \brief Struct with setup and buffers for special atom communication */
67 struct gmx_domdec_specat_comm_t {
68 /* The number of indices to receive during the setup */
69 int nreq[DIM][2][2]; /**< The nr. of atoms requested, per DIM, direction and direct/indirect */
70 /* The atoms to send */
71 gmx_specatsend_t spas[DIM][2]; /**< The communication setup per DIM, direction */
72 gmx_bool *bSendAtom; /**< Work buffer that tells if spec.atoms should be sent */
73 int bSendAtom_nalloc; /**< Allocation size of \p bSendAtom */
74 /* Send buffers */
75 int *ibuf; /**< Integer send buffer */
76 int ibuf_nalloc; /**< Allocation size of \p ibuf */
77 rvec *vbuf; /**< rvec send buffer */
78 int vbuf_nalloc; /**< Allocation size of \p vbuf */
79 rvec *vbuf2; /**< rvec send buffer */
80 int vbuf2_nalloc; /**< Allocation size of \p vbuf2 */
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 */
85 /* The atom indices we need from the surrounding cells.
86 * We can gather the indices over nthread threads.
88 int nthread; /**< Number of threads used for spec.atom communication */
89 ind_req_t *ireq; /**< Index request buffer per thread, allocation size \p nthread */
92 /*! \brief Communicates the force for special atoms, the shift forces are reduced with \p fshift != NULL */
93 void dd_move_f_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
94 rvec *f, rvec *fshift);
96 /*! \brief Communicates the coordinates for special atoms
98 * \param[in] dd Domain decomposition struct
99 * \param[in] spac Special atom communication struct
100 * \param[in] box Box, used for pbc
101 * \param[in,out] x0 Vector to communicate
102 * \param[in,out] x1 Vector to communicate, when != NULL
103 * \param[in] bX1IsCoord Tells is \p x1 is a coordinate vector (needs pbc)
105 void dd_move_x_specat(gmx_domdec_t *dd, gmx_domdec_specat_comm_t *spac,
106 matrix box,
107 rvec *x0,
108 rvec *x1, gmx_bool bX1IsCoord);
110 /*! \brief Sets up the communication for special atoms
112 * \param[in] dd Domain decomposition struct
113 * \param[in] ireq List of requested atom indices
114 * \param[in,out] spac Special atom communication struct
115 * \param[out] ga2la_specat Global to local special atom index
116 * \param[in] at_start Index in local state where to start storing communicated atoms
117 * \param[in] vbuf_fac Buffer factor, 1 or 2 for communicating 1 or 2 vectors
118 * \param[in] specat_type Name of the special atom, used for error message
119 * \param[in] add_err Text to add at the end of error message when atoms can't be found
121 int setup_specat_communication(gmx_domdec_t *dd,
122 ind_req_t *ireq,
123 gmx_domdec_specat_comm_t *spac,
124 gmx_hash_t *ga2la_specat,
125 int at_start,
126 int vbuf_fac,
127 const char *specat_type,
128 const char *add_err);
130 /*! \brief Initialize a special communication struct */
131 gmx_domdec_specat_comm_t *specat_comm_init(int nthread);
133 #endif