Warn when a value is unused
[gromacs.git] / include / gmx_ga2la.h
blob9e385eed1e3301ac39dae515cf95c0589a46262a
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 #include "typedefs.h"
40 #include "smalloc.h"
42 #ifdef __cplusplus
43 extern "C" {
44 #endif
46 typedef struct {
47 int la;
48 int cell;
49 } gmx_laa_t;
51 typedef struct {
52 int ga;
53 int la;
54 int cell;
55 int next;
56 } gmx_lal_t;
58 typedef struct gmx_ga2la {
59 gmx_bool bAll;
60 int mod;
61 int nalloc;
62 gmx_laa_t *laa;
63 gmx_lal_t *lal;
64 int start_space_search;
65 } t_gmx_ga2la;
67 /* Clear all the entries in the ga2la list */
68 static void ga2la_clear(gmx_ga2la_t ga2la)
70 int i;
72 if (ga2la->bAll)
74 for(i=0; i<ga2la->nalloc; i++)
76 ga2la->laa[i].cell = -1;
79 else
81 for(i=0; i<ga2la->nalloc; i++)
83 ga2la->lal[i].ga = -1;
84 ga2la->lal[i].next = -1;
86 ga2la->start_space_search = ga2la->mod;
90 static gmx_ga2la_t ga2la_init(int nat_tot,int nat_loc)
92 gmx_ga2la_t ga2la;
94 snew(ga2la,1);
96 /* There are two methods implemented for finding the local atom number
97 * belonging to a global atom number:
98 * 1) a simple, direct array
99 * 2) a list of linked lists indexed with the global number modulo mod.
100 * Memory requirements:
101 * 1) nat_tot*2 ints
102 * 2) nat_loc*(2+1-2(1-e^-1/2))*4 ints
103 * where nat_loc is the number of atoms in the home + communicated zones.
104 * Method 1 is faster for low parallelization, 2 for high parallelization.
105 * We switch to method 2 when it uses less than half the memory method 1.
107 ga2la->bAll = (nat_tot < 1024 || 9*nat_loc >= nat_tot);
108 if (ga2la->bAll)
110 ga2la->nalloc = nat_tot;
111 snew(ga2la->laa,ga2la->nalloc);
113 else
115 /* Make the direct list twice as long as the number of local atoms.
116 * The fraction of entries in the list with:
117 * 0 size lists: e^-1/f
118 * >=1 size lists: 1 - e^-1/f
119 * where f is: the direct list length / #local atoms
120 * The fraction of atoms not in the direct list is: 1-f(1-e^-1/f).
122 ga2la->mod = 2*nat_loc;
123 ga2la->nalloc = over_alloc_dd(ga2la->mod);
124 snew(ga2la->lal,ga2la->nalloc);
127 ga2la_clear(ga2la);
129 return ga2la;
132 /* Set the ga2la entry for global atom a_gl to local atom a_loc and cell. */
133 static void ga2la_set(gmx_ga2la_t ga2la,int a_gl,int a_loc,int cell)
135 int ind,ind_prev,i;
137 if (ga2la->bAll)
139 ga2la->laa[a_gl].la = a_loc;
140 ga2la->laa[a_gl].cell = cell;
142 return;
145 ind = a_gl % ga2la->mod;
147 if (ga2la->lal[ind].ga >= 0)
149 /* Search the last entry in the linked list for this index */
150 ind_prev = ind;
151 while(ga2la->lal[ind_prev].next >= 0)
153 ind_prev = ga2la->lal[ind_prev].next;
155 /* Search for space in the array */
156 ind = ga2la->start_space_search;
157 while (ind < ga2la->nalloc && ga2la->lal[ind].ga >= 0)
159 ind++;
161 /* If we are at the end of the list we need to increase the size */
162 if (ind == ga2la->nalloc)
164 ga2la->nalloc = over_alloc_dd(ind+1);
165 srenew(ga2la->lal,ga2la->nalloc);
166 for(i=ind; i<ga2la->nalloc; i++)
168 ga2la->lal[i].ga = -1;
169 ga2la->lal[i].next = -1;
172 ga2la->lal[ind_prev].next = ind;
174 ga2la->start_space_search = ind + 1;
176 ga2la->lal[ind].ga = a_gl;
177 ga2la->lal[ind].la = a_loc;
178 ga2la->lal[ind].cell = cell;
181 /* Delete the ga2la entry for global atom a_gl */
182 static void ga2la_del(gmx_ga2la_t ga2la,int a_gl)
184 int ind,ind_prev;
186 if (ga2la->bAll)
188 ga2la->laa[a_gl].cell = -1;
190 return;
193 ind_prev = -1;
194 ind = a_gl % ga2la->mod;
197 if (ga2la->lal[ind].ga == a_gl)
199 if (ind_prev >= 0)
201 ga2la->lal[ind_prev].next = ga2la->lal[ind].next;
203 /* This index is a linked entry, so we free an entry.
204 * Check if we are creating the first empty space.
206 if (ind < ga2la->start_space_search)
208 ga2la->start_space_search = ind;
211 ga2la->lal[ind].ga = -1;
212 ga2la->lal[ind].cell = -1;
213 ga2la->lal[ind].next = -1;
215 return;
217 ind_prev = ind;
218 ind = ga2la->lal[ind].next;
220 while (ind >= 0);
222 return;
225 /* Change the local atom for present ga2la entry for global atom a_gl */
226 static void ga2la_change_la(gmx_ga2la_t ga2la,int a_gl,int a_loc)
228 int ind;
230 if (ga2la->bAll)
232 ga2la->laa[a_gl].la = a_loc;
234 return;
237 ind = a_gl % ga2la->mod;
240 if (ga2la->lal[ind].ga == a_gl)
242 ga2la->lal[ind].la = a_loc;
244 return;
246 ind = ga2la->lal[ind].next;
248 while (ind >= 0);
250 return;
253 /* Returns if the global atom a_gl available locally.
254 * Sets the local atom and cell,
255 * cell can be larger than the number of zones,
256 * in which case it indicates that it is more than one cell away
257 * in zone cell - #zones.
259 static gmx_bool ga2la_get(const gmx_ga2la_t ga2la,int a_gl,int *a_loc,int *cell)
261 int ind;
263 if (ga2la->bAll)
265 *a_loc = ga2la->laa[a_gl].la;
266 *cell = ga2la->laa[a_gl].cell;
268 return (ga2la->laa[a_gl].cell >= 0);
271 ind = a_gl % ga2la->mod;
274 if (ga2la->lal[ind].ga == a_gl)
276 *a_loc = ga2la->lal[ind].la;
277 *cell = ga2la->lal[ind].cell;
279 return TRUE;
281 ind = ga2la->lal[ind].next;
283 while (ind >= 0);
285 return FALSE;
288 /* Returns if the global atom a_gl is a home atom.
289 * Sets the local atom.
291 static gmx_bool ga2la_get_home(const gmx_ga2la_t ga2la,int a_gl,int *a_loc)
293 int ind;
295 if (ga2la->bAll)
297 *a_loc = ga2la->laa[a_gl].la;
299 return (ga2la->laa[a_gl].cell == 0);
302 ind = a_gl % ga2la->mod;
305 if (ga2la->lal[ind].ga == a_gl)
307 if (ga2la->lal[ind].cell == 0)
309 *a_loc = ga2la->lal[ind].la;
311 return TRUE;
313 else
315 return FALSE;
318 ind = ga2la->lal[ind].next;
320 while (ind >= 0);
322 return FALSE;
325 /* Returns if the global atom a_gl is a home atom.
327 static gmx_bool ga2la_is_home(const gmx_ga2la_t ga2la,int a_gl)
329 int ind;
331 if (ga2la->bAll)
333 return (ga2la->laa[a_gl].cell == 0);
336 ind = a_gl % ga2la->mod;
339 if (ga2la->lal[ind].ga == a_gl)
341 return (ga2la->lal[ind].cell == 0);
343 ind = ga2la->lal[ind].next;
345 while (ind >= 0);
347 return FALSE;
350 #ifdef __cplusplus
352 #endif
354 #endif /* _gmx_ga2la_h */