2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2010,2014,2015, 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 /*! \libinternal \file
39 * Defines structures and functions for mapping from global to local atom
40 * indices. The functions are performance critical and should be inlined.
43 * \ingroup module_domdec
45 * \author Berk Hess <hess@kth.se>
48 #ifndef GMX_DOMDEC_GA2LA_H
49 #define GMX_DOMDEC_GA2LA_H
51 #include "gromacs/legacyheaders/types/commrec.h"
52 #include "gromacs/utility/basedefinitions.h"
53 #include "gromacs/utility/smalloc.h"
55 /*! \libinternal \brief Structure for the local atom info for a plain list */
57 int la
; /**< The local atom index */
58 int cell
; /**< The DD zone index for neighboring domains, zone+zone otherwise */
61 /*! \libinternal \brief Structure for the local atom info for a hash table */
63 int ga
; /**< The global atom index */
64 int la
; /**< The local atom index */
65 int cell
; /**< The DD zone index for neighboring domains, zone+zone otherwise */
66 int next
; /**< Index in the list of the next element with the same hash, -1 if none */
69 /*! \libinternal \brief Structure for all global to local mapping information */
71 gmx_bool bDirectList
; /**< Use a direct list */
72 int mod
; /**< The hash size */
73 int nalloc
; /**< The alloction size of laa or la1 */
74 gmx_laa_t
*laa
; /**< The direct list */
75 gmx_lal_t
*lal
; /**< The hash table list */
76 int start_space_search
; /**< Index in lal at which to start looking for empty space */
79 /*! \brief Clear all the entries in the ga2la list
81 * \param[in,out] ga2la The global to local atom struct
83 static void ga2la_clear(gmx_ga2la_t
*ga2la
)
87 if (ga2la
->bDirectList
)
89 for (i
= 0; i
< ga2la
->nalloc
; i
++)
91 ga2la
->laa
[i
].cell
= -1;
96 for (i
= 0; i
< ga2la
->nalloc
; i
++)
98 ga2la
->lal
[i
].ga
= -1;
99 ga2la
->lal
[i
].next
= -1;
101 ga2la
->start_space_search
= ga2la
->mod
;
105 /*! \brief Initializes and returns a pointer to a gmx_ga2la_t structure
107 * \param[in] natoms_total The total number of atoms in the system
108 * \param[in] natoms_local An estimate of the number of home+communicated atoms
109 * \return a pointer to an initialized gmx_ga2la_t struct
111 static gmx_ga2la_t
*ga2la_init(int natoms_total
, int natoms_local
)
117 /* There are two methods implemented for finding the local atom number
118 * belonging to a global atom number:
119 * 1) a simple, direct array
120 * 2) a hash table consisting of list of linked lists indexed with
121 * the global number modulo mod.
122 * Memory requirements:
124 * 2) nat_loc*(2+1-2(1-e^-1/2))*4 ints
125 * where nat_loc is the number of atoms in the home + communicated zones.
126 * Method 1 is faster for low parallelization, 2 for high parallelization.
127 * We switch to method 2 when it uses less than half the memory method 1.
129 ga2la
->bDirectList
= (natoms_total
<= 1024 ||
130 natoms_total
<= natoms_local
*9);
132 if (ga2la
->bDirectList
)
134 ga2la
->nalloc
= natoms_total
;
135 snew(ga2la
->laa
, ga2la
->nalloc
);
139 /* Make the direct list twice as long as the number of local atoms.
140 * The fraction of entries in the list with:
141 * 0 size lists: e^-1/f
142 * >=1 size lists: 1 - e^-1/f
143 * where f is: the direct list length / #local atoms
144 * The fraction of atoms not in the direct list is: 1-f(1-e^-1/f).
146 ga2la
->mod
= 2*natoms_local
;
147 ga2la
->nalloc
= over_alloc_dd(ga2la
->mod
);
148 snew(ga2la
->lal
, ga2la
->nalloc
);
156 /*! \brief Sets the ga2la entry for global atom a_gl
158 * \param[in,out] ga2la The global to local atom struct
159 * \param[in] a_gl The global atom index
160 * \param[in] a_loc The local atom index
161 * \param[in] cell The cell index
163 static void ga2la_set(gmx_ga2la_t
*ga2la
, int a_gl
, int a_loc
, int cell
)
165 int ind
, ind_prev
, i
;
167 if (ga2la
->bDirectList
)
169 ga2la
->laa
[a_gl
].la
= a_loc
;
170 ga2la
->laa
[a_gl
].cell
= cell
;
175 ind
= a_gl
% ga2la
->mod
;
177 if (ga2la
->lal
[ind
].ga
>= 0)
179 /* Search the last entry in the linked list for this index */
181 while (ga2la
->lal
[ind_prev
].next
>= 0)
183 ind_prev
= ga2la
->lal
[ind_prev
].next
;
185 /* Search for space in the array */
186 ind
= ga2la
->start_space_search
;
187 while (ind
< ga2la
->nalloc
&& ga2la
->lal
[ind
].ga
>= 0)
191 /* If we are at the end of the list we need to increase the size */
192 if (ind
== ga2la
->nalloc
)
194 ga2la
->nalloc
= over_alloc_dd(ind
+1);
195 srenew(ga2la
->lal
, ga2la
->nalloc
);
196 for (i
= ind
; i
< ga2la
->nalloc
; i
++)
198 ga2la
->lal
[i
].ga
= -1;
199 ga2la
->lal
[i
].next
= -1;
202 ga2la
->lal
[ind_prev
].next
= ind
;
204 ga2la
->start_space_search
= ind
+ 1;
206 ga2la
->lal
[ind
].ga
= a_gl
;
207 ga2la
->lal
[ind
].la
= a_loc
;
208 ga2la
->lal
[ind
].cell
= cell
;
211 /*! \brief Delete the ga2la entry for global atom a_gl
213 * \param[in,out] ga2la The global to local atom struct
214 * \param[in] a_gl The global atom index
216 static void ga2la_del(gmx_ga2la_t
*ga2la
, int a_gl
)
220 if (ga2la
->bDirectList
)
222 ga2la
->laa
[a_gl
].cell
= -1;
228 ind
= a_gl
% ga2la
->mod
;
231 if (ga2la
->lal
[ind
].ga
== a_gl
)
235 ga2la
->lal
[ind_prev
].next
= ga2la
->lal
[ind
].next
;
237 /* This index is a linked entry, so we free an entry.
238 * Check if we are creating the first empty space.
240 if (ind
< ga2la
->start_space_search
)
242 ga2la
->start_space_search
= ind
;
245 ga2la
->lal
[ind
].ga
= -1;
246 ga2la
->lal
[ind
].cell
= -1;
247 ga2la
->lal
[ind
].next
= -1;
252 ind
= ga2la
->lal
[ind
].next
;
259 /*! \brief Change the local atom for present ga2la entry for global atom a_gl
261 * \param[in,out] ga2la The global to local atom struct
262 * \param[in] a_gl The global atom index
263 * \param[in] a_loc The new local atom index
265 static void ga2la_change_la(gmx_ga2la_t
*ga2la
, int a_gl
, int a_loc
)
269 if (ga2la
->bDirectList
)
271 ga2la
->laa
[a_gl
].la
= a_loc
;
276 ind
= a_gl
% ga2la
->mod
;
279 if (ga2la
->lal
[ind
].ga
== a_gl
)
281 ga2la
->lal
[ind
].la
= a_loc
;
285 ind
= ga2la
->lal
[ind
].next
;
292 /*! \brief Returns if the global atom a_gl available locally
294 * \param[in] ga2la The global to local atom struct
295 * \param[in] a_gl The global atom index
296 * \param[out] a_loc If the return value is TRUE, the local atom index
297 * \param[out] cell If the return value is TRUE, the zone or for atoms more than one cell away zone+nzone
298 * \return if the global atom a_gl available locally
300 static gmx_bool
ga2la_get(const gmx_ga2la_t
*ga2la
, int a_gl
, int *a_loc
, int *cell
)
304 if (ga2la
->bDirectList
)
306 *a_loc
= ga2la
->laa
[a_gl
].la
;
307 *cell
= ga2la
->laa
[a_gl
].cell
;
309 return (ga2la
->laa
[a_gl
].cell
>= 0);
312 ind
= a_gl
% ga2la
->mod
;
315 if (ga2la
->lal
[ind
].ga
== a_gl
)
317 *a_loc
= ga2la
->lal
[ind
].la
;
318 *cell
= ga2la
->lal
[ind
].cell
;
322 ind
= ga2la
->lal
[ind
].next
;
329 /*! \brief Returns if the global atom a_gl is a home atom
331 * \param[in] ga2la The global to local atom struct
332 * \param[in] a_gl The global atom index
333 * \param[out] a_loc If the return value is TRUE, the local atom index
334 * \return if the global atom a_gl is a home atom
336 static gmx_bool
ga2la_get_home(const gmx_ga2la_t
*ga2la
, int a_gl
, int *a_loc
)
340 if (ga2la
->bDirectList
)
342 *a_loc
= ga2la
->laa
[a_gl
].la
;
344 return (ga2la
->laa
[a_gl
].cell
== 0);
347 ind
= a_gl
% ga2la
->mod
;
350 if (ga2la
->lal
[ind
].ga
== a_gl
)
352 if (ga2la
->lal
[ind
].cell
== 0)
354 *a_loc
= ga2la
->lal
[ind
].la
;
363 ind
= ga2la
->lal
[ind
].next
;
370 /*! \brief Returns if the global atom a_gl is a home atom
372 * \param[in] ga2la The global to local atom struct
373 * \param[in] a_gl The global atom index
374 * \return if the global atom a_gl is a home atom
376 static gmx_bool
ga2la_is_home(const gmx_ga2la_t
*ga2la
, int a_gl
)
380 if (ga2la
->bDirectList
)
382 return (ga2la
->laa
[a_gl
].cell
== 0);
385 ind
= a_gl
% ga2la
->mod
;
388 if (ga2la
->lal
[ind
].ga
== a_gl
)
390 return (ga2la
->lal
[ind
].cell
== 0);
392 ind
= ga2la
->lal
[ind
].next
;