Fix remaining copyright headers
[alexxy/gromacs.git] / src / gromacs / legacyheaders / gmx_ga2la.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2010, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef _gmx_ga2la_h
38 #define _gmx_ga2la_h
39
40 #include "typedefs.h"
41 #include "smalloc.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47 typedef struct {
48     int  la;
49     int  cell;
50 } gmx_laa_t;
51
52 typedef struct {
53     int  ga;
54     int  la;
55     int  cell;
56     int  next;
57 } gmx_lal_t;
58
59 typedef struct gmx_ga2la {
60     gmx_bool      bAll;
61     int           mod;
62     int           nalloc;
63     gmx_laa_t    *laa;
64     gmx_lal_t    *lal;
65     int           start_space_search;
66 } t_gmx_ga2la;
67
68 /* Clear all the entries in the ga2la list */
69 static void ga2la_clear(gmx_ga2la_t ga2la)
70 {
71     int i;
72
73     if (ga2la->bAll)
74     {
75         for (i = 0; i < ga2la->nalloc; i++)
76         {
77             ga2la->laa[i].cell = -1;
78         }
79     }
80     else
81     {
82         for (i = 0; i < ga2la->nalloc; i++)
83         {
84             ga2la->lal[i].ga   = -1;
85             ga2la->lal[i].next = -1;
86         }
87         ga2la->start_space_search = ga2la->mod;
88     }
89 }
90
91 static gmx_ga2la_t ga2la_init(int nat_tot, int nat_loc)
92 {
93     gmx_ga2la_t ga2la;
94
95     snew(ga2la, 1);
96
97     /* There are two methods implemented for finding the local atom number
98      * belonging to a global atom number:
99      * 1) a simple, direct array
100      * 2) a list of linked lists indexed with the global number modulo mod.
101      * Memory requirements:
102      * 1) nat_tot*2 ints
103      * 2) nat_loc*(2+1-2(1-e^-1/2))*4 ints
104      * where nat_loc is the number of atoms in the home + communicated zones.
105      * Method 1 is faster for low parallelization, 2 for high parallelization.
106      * We switch to method 2 when it uses less than half the memory method 1.
107      */
108     ga2la->bAll = (nat_tot < 1024 || 9*nat_loc >= nat_tot);
109     if (ga2la->bAll)
110     {
111         ga2la->nalloc = nat_tot;
112         snew(ga2la->laa, ga2la->nalloc);
113     }
114     else
115     {
116         /* Make the direct list twice as long as the number of local atoms.
117          * The fraction of entries in the list with:
118          * 0   size lists: e^-1/f
119          * >=1 size lists: 1 - e^-1/f
120          * where f is: the direct list length / #local atoms
121          * The fraction of atoms not in the direct list is: 1-f(1-e^-1/f).
122          */
123         ga2la->mod    = 2*nat_loc;
124         ga2la->nalloc = over_alloc_dd(ga2la->mod);
125         snew(ga2la->lal, ga2la->nalloc);
126     }
127
128     ga2la_clear(ga2la);
129
130     return ga2la;
131 }
132
133 /* Set the ga2la entry for global atom a_gl to local atom a_loc and cell. */
134 static void ga2la_set(gmx_ga2la_t ga2la, int a_gl, int a_loc, int cell)
135 {
136     int ind, ind_prev, i;
137
138     if (ga2la->bAll)
139     {
140         ga2la->laa[a_gl].la   = a_loc;
141         ga2la->laa[a_gl].cell = cell;
142
143         return;
144     }
145
146     ind = a_gl % ga2la->mod;
147
148     if (ga2la->lal[ind].ga >= 0)
149     {
150         /* Search the last entry in the linked list for this index */
151         ind_prev = ind;
152         while (ga2la->lal[ind_prev].next >= 0)
153         {
154             ind_prev = ga2la->lal[ind_prev].next;
155         }
156         /* Search for space in the array */
157         ind = ga2la->start_space_search;
158         while (ind < ga2la->nalloc && ga2la->lal[ind].ga >= 0)
159         {
160             ind++;
161         }
162         /* If we are at the end of the list we need to increase the size */
163         if (ind == ga2la->nalloc)
164         {
165             ga2la->nalloc = over_alloc_dd(ind+1);
166             srenew(ga2la->lal, ga2la->nalloc);
167             for (i = ind; i < ga2la->nalloc; i++)
168             {
169                 ga2la->lal[i].ga   = -1;
170                 ga2la->lal[i].next = -1;
171             }
172         }
173         ga2la->lal[ind_prev].next = ind;
174
175         ga2la->start_space_search = ind + 1;
176     }
177     ga2la->lal[ind].ga   = a_gl;
178     ga2la->lal[ind].la   = a_loc;
179     ga2la->lal[ind].cell = cell;
180 }
181
182 /* Delete the ga2la entry for global atom a_gl */
183 static void ga2la_del(gmx_ga2la_t ga2la, int a_gl)
184 {
185     int ind, ind_prev;
186
187     if (ga2la->bAll)
188     {
189         ga2la->laa[a_gl].cell = -1;
190
191         return;
192     }
193
194     ind_prev = -1;
195     ind      = a_gl % ga2la->mod;
196     do
197     {
198         if (ga2la->lal[ind].ga == a_gl)
199         {
200             if (ind_prev >= 0)
201             {
202                 ga2la->lal[ind_prev].next = ga2la->lal[ind].next;
203
204                 /* This index is a linked entry, so we free an entry.
205                  * Check if we are creating the first empty space.
206                  */
207                 if (ind < ga2la->start_space_search)
208                 {
209                     ga2la->start_space_search = ind;
210                 }
211             }
212             ga2la->lal[ind].ga   = -1;
213             ga2la->lal[ind].cell = -1;
214             ga2la->lal[ind].next = -1;
215
216             return;
217         }
218         ind_prev = ind;
219         ind      = ga2la->lal[ind].next;
220     }
221     while (ind >= 0);
222
223     return;
224 }
225
226 /* Change the local atom for present ga2la entry for global atom a_gl */
227 static void ga2la_change_la(gmx_ga2la_t ga2la, int a_gl, int a_loc)
228 {
229     int ind;
230
231     if (ga2la->bAll)
232     {
233         ga2la->laa[a_gl].la = a_loc;
234
235         return;
236     }
237
238     ind = a_gl % ga2la->mod;
239     do
240     {
241         if (ga2la->lal[ind].ga == a_gl)
242         {
243             ga2la->lal[ind].la = a_loc;
244
245             return;
246         }
247         ind = ga2la->lal[ind].next;
248     }
249     while (ind >= 0);
250
251     return;
252 }
253
254 /* Returns if the global atom a_gl available locally.
255  * Sets the local atom and cell,
256  * cell can be larger than the number of zones,
257  * in which case it indicates that it is more than one cell away
258  * in zone cell - #zones.
259  */
260 static gmx_bool ga2la_get(const gmx_ga2la_t ga2la, int a_gl, int *a_loc, int *cell)
261 {
262     int ind;
263
264     if (ga2la->bAll)
265     {
266         *a_loc = ga2la->laa[a_gl].la;
267         *cell  = ga2la->laa[a_gl].cell;
268
269         return (ga2la->laa[a_gl].cell >= 0);
270     }
271
272     ind = a_gl % ga2la->mod;
273     do
274     {
275         if (ga2la->lal[ind].ga == a_gl)
276         {
277             *a_loc = ga2la->lal[ind].la;
278             *cell  = ga2la->lal[ind].cell;
279
280             return TRUE;
281         }
282         ind = ga2la->lal[ind].next;
283     }
284     while (ind >= 0);
285
286     return FALSE;
287 }
288
289 /* Returns if the global atom a_gl is a home atom.
290  * Sets the local atom.
291  */
292 static gmx_bool ga2la_get_home(const gmx_ga2la_t ga2la, int a_gl, int *a_loc)
293 {
294     int ind;
295
296     if (ga2la->bAll)
297     {
298         *a_loc = ga2la->laa[a_gl].la;
299
300         return (ga2la->laa[a_gl].cell == 0);
301     }
302
303     ind = a_gl % ga2la->mod;
304     do
305     {
306         if (ga2la->lal[ind].ga == a_gl)
307         {
308             if (ga2la->lal[ind].cell == 0)
309             {
310                 *a_loc = ga2la->lal[ind].la;
311
312                 return TRUE;
313             }
314             else
315             {
316                 return FALSE;
317             }
318         }
319         ind = ga2la->lal[ind].next;
320     }
321     while (ind >= 0);
322
323     return FALSE;
324 }
325
326 /* Returns if the global atom a_gl is a home atom.
327  */
328 static gmx_bool ga2la_is_home(const gmx_ga2la_t ga2la, int a_gl)
329 {
330     int ind;
331
332     if (ga2la->bAll)
333     {
334         return (ga2la->laa[a_gl].cell == 0);
335     }
336
337     ind = a_gl % ga2la->mod;
338     do
339     {
340         if (ga2la->lal[ind].ga == a_gl)
341         {
342             return (ga2la->lal[ind].cell == 0);
343         }
344         ind = ga2la->lal[ind].next;
345     }
346     while (ind >= 0);
347
348     return FALSE;
349 }
350
351 #ifdef __cplusplus
352 }
353 #endif
354
355 #endif /* _gmx_ga2la_h */