Merge branch 'master' of git@git.gromacs.org:gromacs
[gromacs/rigid-bodies.git] / include / gmx_ga2la.h
blobc872cec00cb04a6a2185fd3a4cbb0e52832b41d0
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _gmx_ga2la_h
37 #define _gmx_ga2la_h
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
47 typedef struct {
48 int la;
49 int cell;
50 } gmx_laa_t;
52 typedef struct {
53 int ga;
54 int la;
55 int cell;
56 int next;
57 } gmx_lal_t;
59 typedef struct gmx_ga2la {
60 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;
68 /* Clear all the entries in the ga2la list */
69 static void ga2la_clear(gmx_ga2la_t ga2la)
71 int i;
73 if (ga2la->bAll)
75 for(i=0; i<ga2la->nalloc; i++)
77 ga2la->laa[i].cell = -1;
80 else
82 for(i=0; i<ga2la->nalloc; i++)
84 ga2la->lal[i].ga = -1;
85 ga2la->lal[i].next = -1;
87 ga2la->start_space_search = ga2la->mod;
91 static gmx_ga2la_t ga2la_init(int nat_tot,int nat_loc)
93 gmx_ga2la_t ga2la;
95 snew(ga2la,1);
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.
108 ga2la->bAll = (nat_tot < 1024 || 9*nat_loc >= nat_tot);
109 if (ga2la->bAll)
111 ga2la->nalloc = nat_tot;
112 snew(ga2la->laa,ga2la->nalloc);
114 else
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).
123 ga2la->mod = 2*nat_loc;
124 ga2la->nalloc = over_alloc_dd(ga2la->mod);
125 snew(ga2la->lal,ga2la->nalloc);
128 ga2la_clear(ga2la);
130 return ga2la;
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)
136 int ind,ind_prev,i;
138 if (ga2la->bAll)
140 ga2la->laa[a_gl].la = a_loc;
141 ga2la->laa[a_gl].cell = cell;
143 return;
146 ind = a_gl % ga2la->mod;
148 if (ga2la->lal[ind].ga >= 0)
150 /* Search the last entry in the linked list for this index */
151 ind_prev = ind;
152 while(ga2la->lal[ind_prev].next >= 0)
154 ind_prev = ga2la->lal[ind_prev].next;
156 /* Search for space in the array */
157 ind = ga2la->start_space_search;
158 while (ind < ga2la->nalloc && ga2la->lal[ind].ga >= 0)
160 ind++;
162 /* If we are at the end of the list we need to increase the size */
163 if (ind == ga2la->nalloc)
165 ga2la->nalloc = over_alloc_dd(ind+1);
166 srenew(ga2la->lal,ga2la->nalloc);
167 for(i=ind; i<ga2la->nalloc; i++)
169 ga2la->lal[i].ga = -1;
170 ga2la->lal[i].next = -1;
173 ga2la->lal[ind_prev].next = ind;
175 ga2la->start_space_search = ind + 1;
177 ga2la->lal[ind].ga = a_gl;
178 ga2la->lal[ind].la = a_loc;
179 ga2la->lal[ind].cell = cell;
182 /* Delete the ga2la entry for global atom a_gl */
183 static void ga2la_del(gmx_ga2la_t ga2la,int a_gl)
185 int ind,ind_prev;
187 if (ga2la->bAll)
189 ga2la->laa[a_gl].cell = -1;
191 return;
194 ind_prev = -1;
195 ind = a_gl % ga2la->mod;
198 if (ga2la->lal[ind].ga == a_gl)
200 if (ind_prev >= 0)
202 ga2la->lal[ind_prev].next = ga2la->lal[ind].next;
204 /* This index is a linked entry, so we free an entry.
205 * Check if we are creating the first empty space.
207 if (ind < ga2la->start_space_search)
209 ga2la->start_space_search = ind;
212 ga2la->lal[ind].ga = -1;
213 ga2la->lal[ind].cell = -1;
214 ga2la->lal[ind].next = -1;
216 return;
218 ind_prev = ind;
219 ind = ga2la->lal[ind].next;
221 while (ind >= 0);
223 return;
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)
229 int ind;
231 if (ga2la->bAll)
233 ga2la->laa[a_gl].la = a_loc;
235 return;
238 ind = a_gl % ga2la->mod;
241 if (ga2la->lal[ind].ga == a_gl)
243 ga2la->lal[ind].la = a_loc;
245 return;
247 ind = ga2la->lal[ind].next;
249 while (ind >= 0);
251 return;
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.
260 static bool ga2la_get(const gmx_ga2la_t ga2la,int a_gl,int *a_loc,int *cell)
262 int ind;
264 if (ga2la->bAll)
266 *a_loc = ga2la->laa[a_gl].la;
267 *cell = ga2la->laa[a_gl].cell;
269 return (ga2la->laa[a_gl].cell >= 0);
272 ind = a_gl % ga2la->mod;
275 if (ga2la->lal[ind].ga == a_gl)
277 *a_loc = ga2la->lal[ind].la;
278 *cell = ga2la->lal[ind].cell;
280 return TRUE;
282 ind = ga2la->lal[ind].next;
284 while (ind >= 0);
286 return FALSE;
289 /* Returns if the global atom a_gl is a home atom.
290 * Sets the local atom.
292 static bool ga2la_get_home(const gmx_ga2la_t ga2la,int a_gl,int *a_loc)
294 int ind;
296 if (ga2la->bAll)
298 *a_loc = ga2la->laa[a_gl].la;
300 return (ga2la->laa[a_gl].cell == 0);
303 ind = a_gl % ga2la->mod;
306 if (ga2la->lal[ind].ga == a_gl)
308 if (ga2la->lal[ind].cell == 0)
310 *a_loc = ga2la->lal[ind].la;
312 return TRUE;
314 else
316 return FALSE;
319 ind = ga2la->lal[ind].next;
321 while (ind >= 0);
323 return FALSE;
326 /* Returns if the global atom a_gl is a home atom.
328 static bool ga2la_is_home(const gmx_ga2la_t ga2la,int a_gl)
330 int ind;
332 if (ga2la->bAll)
334 return (ga2la->laa[a_gl].cell == 0);
337 ind = a_gl % ga2la->mod;
340 if (ga2la->lal[ind].ga == a_gl)
342 return (ga2la->lal[ind].cell == 0);
344 ind = ga2la->lal[ind].next;
346 while (ind >= 0);
348 return FALSE;
351 #ifdef __cplusplus
353 #endif
355 #endif /* _gmx_ga2la_h */