fixed rest of include typos in headers
[gromacs/adressmacs.git] / include / gmx_ga2la.h
blob44cd4a5ad8c71f4412a76cbe22931af8ed6d1b06
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 #include "typedefs.h"
44 #include "smalloc.h"
46 #ifdef __cplusplus
47 extern "C" {
48 #endif
50 typedef struct {
51 int la;
52 int cell;
53 } gmx_laa_t;
55 typedef struct {
56 int ga;
57 int la;
58 int cell;
59 int next;
60 } gmx_lal_t;
62 typedef struct gmx_ga2la {
63 bool bAll;
64 int mod;
65 int nalloc;
66 gmx_laa_t *laa;
67 gmx_lal_t *lal;
68 int start_space_search;
69 } t_gmx_ga2la;
71 /* Clear all the entries in the ga2la list */
72 static void ga2la_clear(gmx_ga2la_t ga2la)
74 int i;
76 if (ga2la->bAll)
78 for(i=0; i<ga2la->nalloc; i++)
80 ga2la->laa[i].cell = -1;
83 else
85 for(i=0; i<ga2la->nalloc; i++)
87 ga2la->lal[i].ga = -1;
88 ga2la->lal[i].next = -1;
90 ga2la->start_space_search = ga2la->mod;
94 static gmx_ga2la_t ga2la_init(int nat_tot,int nat_loc)
96 gmx_ga2la_t ga2la;
98 snew(ga2la,1);
100 /* There are two methods implemented for finding the local atom number
101 * belonging to a global atom number:
102 * 1) a simple, direct array
103 * 2) a list of linked lists indexed with the global number modulo mod.
104 * Memory requirements:
105 * 1) nat_tot*2 ints
106 * 2) nat_loc*(2+1-2(1-e^-1/2))*4 ints
107 * where nat_loc is the number of atoms in the home + communicated zones.
108 * Method 1 is faster for low parallelization, 2 for high parallelization.
109 * We switch to method 2 when it uses less than half the memory method 1.
111 ga2la->bAll = (nat_tot < 1024 || 9*nat_loc >= nat_tot);
112 if (ga2la->bAll)
114 ga2la->nalloc = nat_tot;
115 snew(ga2la->laa,ga2la->nalloc);
117 else
119 /* Make the direct list twice as long as the number of local atoms.
120 * The fraction of entries in the list with:
121 * 0 size lists: e^-1/f
122 * >=1 size lists: 1 - e^-1/f
123 * where f is: the direct list length / #local atoms
124 * The fraction of atoms not in the direct list is: 1-f(1-e^-1/f).
126 ga2la->mod = 2*nat_loc;
127 ga2la->nalloc = over_alloc_dd(ga2la->mod);
128 snew(ga2la->lal,ga2la->nalloc);
131 ga2la_clear(ga2la);
133 return ga2la;
136 /* Set the ga2la entry for global atom a_gl to local atom a_loc and cell. */
137 static void ga2la_set(gmx_ga2la_t ga2la,int a_gl,int a_loc,int cell)
139 int ind,ind_prev,i;
141 if (ga2la->bAll)
143 ga2la->laa[a_gl].la = a_loc;
144 ga2la->laa[a_gl].cell = cell;
146 return;
149 ind = a_gl % ga2la->mod;
151 if (ga2la->lal[ind].ga >= 0)
153 /* Search the last entry in the linked list for this index */
154 ind_prev = ind;
155 while(ga2la->lal[ind_prev].next >= 0)
157 ind_prev = ga2la->lal[ind_prev].next;
159 /* Search for space in the array */
160 ind = ga2la->start_space_search;
161 while (ind < ga2la->nalloc && ga2la->lal[ind].ga >= 0)
163 ind++;
165 /* If we are at the end of the list we need to increase the size */
166 if (ind == ga2la->nalloc)
168 ga2la->nalloc = over_alloc_dd(ind+1);
169 srenew(ga2la->lal,ga2la->nalloc);
170 for(i=ind; i<ga2la->nalloc; i++)
172 ga2la->lal[i].ga = -1;
173 ga2la->lal[i].next = -1;
176 ga2la->lal[ind_prev].next = ind;
178 ga2la->start_space_search = ind + 1;
180 ga2la->lal[ind].ga = a_gl;
181 ga2la->lal[ind].la = a_loc;
182 ga2la->lal[ind].cell = cell;
185 /* Delete the ga2la entry for global atom a_gl */
186 static void ga2la_del(gmx_ga2la_t ga2la,int a_gl)
188 int ind,ind_prev;
190 if (ga2la->bAll)
192 ga2la->laa[a_gl].cell = -1;
194 return;
197 ind_prev = -1;
198 ind = a_gl % ga2la->mod;
201 if (ga2la->lal[ind].ga == a_gl)
203 if (ind_prev >= 0)
205 ga2la->lal[ind_prev].next = ga2la->lal[ind].next;
207 /* This index is a linked entry, so we free an entry.
208 * Check if we are creating the first empty space.
210 if (ind < ga2la->start_space_search)
212 ga2la->start_space_search = ind;
215 ga2la->lal[ind].ga = -1;
216 ga2la->lal[ind].cell = -1;
217 ga2la->lal[ind].next = -1;
219 return;
221 ind_prev = ind;
222 ind = ga2la->lal[ind].next;
224 while (ind >= 0);
226 return;
229 /* Change the local atom for present ga2la entry for global atom a_gl */
230 static void ga2la_change_la(gmx_ga2la_t ga2la,int a_gl,int a_loc)
232 int ind;
234 if (ga2la->bAll)
236 ga2la->laa[a_gl].la = a_loc;
238 return;
241 ind = a_gl % ga2la->mod;
244 if (ga2la->lal[ind].ga == a_gl)
246 ga2la->lal[ind].la = a_loc;
248 return;
250 ind = ga2la->lal[ind].next;
252 while (ind >= 0);
254 return;
257 /* Returns if the global atom a_gl available locally.
258 * Sets the local atom and cell,
259 * cell can be larger than the number of zones,
260 * in which case it indicates that it is more than one cell away
261 * in zone cell - #zones.
263 static bool ga2la_get(const gmx_ga2la_t ga2la,int a_gl,int *a_loc,int *cell)
265 int ind;
267 if (ga2la->bAll)
269 *a_loc = ga2la->laa[a_gl].la;
270 *cell = ga2la->laa[a_gl].cell;
272 return (ga2la->laa[a_gl].cell >= 0);
275 ind = a_gl % ga2la->mod;
278 if (ga2la->lal[ind].ga == a_gl)
280 *a_loc = ga2la->lal[ind].la;
281 *cell = ga2la->lal[ind].cell;
283 return TRUE;
285 ind = ga2la->lal[ind].next;
287 while (ind >= 0);
289 return FALSE;
292 /* Returns if the global atom a_gl is a home atom.
293 * Sets the local atom.
295 static bool ga2la_get_home(const gmx_ga2la_t ga2la,int a_gl,int *a_loc)
297 int ind;
299 if (ga2la->bAll)
301 *a_loc = ga2la->laa[a_gl].la;
303 return (ga2la->laa[a_gl].cell == 0);
306 ind = a_gl % ga2la->mod;
309 if (ga2la->lal[ind].ga == a_gl)
311 if (ga2la->lal[ind].cell == 0)
313 *a_loc = ga2la->lal[ind].la;
315 return TRUE;
317 else
319 return FALSE;
322 ind = ga2la->lal[ind].next;
324 while (ind >= 0);
326 return FALSE;
329 /* Returns if the global atom a_gl is a home atom.
331 static bool ga2la_is_home(const gmx_ga2la_t ga2la,int a_gl)
333 int ind;
335 if (ga2la->bAll)
337 return (ga2la->laa[a_gl].cell == 0);
340 ind = a_gl % ga2la->mod;
343 if (ga2la->lal[ind].ga == a_gl)
345 return (ga2la->lal[ind].cell == 0);
347 ind = ga2la->lal[ind].next;
349 while (ind >= 0);
351 return FALSE;
354 #ifdef __cplusplus
356 #endif
358 #endif /* _gmx_ga2la_h */