Fixed segfaults with empty selections.
[gromacs/qmmm-gamess-us.git] / include / gmx_ga2la.h
blobc06164a2219ff4ee6eae21003d6009c5fc0b2364
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 typedef struct {
44 int la;
45 int cell;
46 } gmx_laa_t;
48 typedef struct {
49 int ga;
50 int la;
51 int cell;
52 int next;
53 } gmx_lal_t;
55 typedef struct gmx_ga2la {
56 bool bAll;
57 int mod;
58 int nalloc;
59 gmx_laa_t *laa;
60 gmx_lal_t *lal;
61 int start_space_search;
62 } t_gmx_ga2la;
64 /* Clear all the entries in the ga2la list */
65 static void ga2la_clear(gmx_ga2la_t ga2la)
67 int i;
69 if (ga2la->bAll)
71 for(i=0; i<ga2la->nalloc; i++)
73 ga2la->laa[i].cell = -1;
76 else
78 for(i=0; i<ga2la->nalloc; i++)
80 ga2la->lal[i].ga = -1;
81 ga2la->lal[i].next = -1;
83 ga2la->start_space_search = ga2la->mod;
87 static gmx_ga2la_t ga2la_init(int nat_tot,int nat_loc)
89 gmx_ga2la_t ga2la;
91 snew(ga2la,1);
93 /* There are two methods implemented for finding the local atom number
94 * belonging to a global atom number:
95 * 1) a simple, direct arrary
96 * 2) a list of linked lists indexed with the global number modulo mod.
97 * Memory requirements:
98 * 1) nat_tot*2 ints
99 * 2) nat_loc*(2+1-2(1-e^-1/2))*4 ints
100 * where nat_loc is the number of atoms in the home + communicated zones.
101 * Method 1 is faster for low parallelization, 2 for high parallelization.
102 * We switch to method 2 when it uses less than half the memory method 1.
104 ga2la->bAll = (9*nat_loc >= nat_tot);
105 if (ga2la->bAll)
107 ga2la->nalloc = nat_tot;
108 snew(ga2la->laa,ga2la->nalloc);
110 else
112 /* Make the direct list twice as large as the number of local atoms.
113 * The fraction of entries in the list with:
114 * 0 size lists: e^-1/f
115 * >=1 size lists: 1 - e^-1/f
116 * where f is: the direct list length / #local atoms
117 * The fraction of atoms not in the direct list is: 1-f(1-e^-1/f).
119 ga2la->mod = 2*nat_loc;
120 ga2la->nalloc = over_alloc_dd(ga2la->mod);
121 snew(ga2la->lal,ga2la->nalloc);
124 ga2la_clear(ga2la);
126 return ga2la;
129 /* Set the ga2la entry for global atom a_gl to local atom a_loc and cell. */
130 static void ga2la_set(gmx_ga2la_t ga2la,int a_gl,int a_loc,int cell)
132 int ind,ind_prev,i;
134 if (ga2la->bAll)
136 ga2la->laa[a_gl].la = a_loc;
137 ga2la->laa[a_gl].cell = cell;
139 return;
142 ind = a_gl % ga2la->mod;
144 if (ga2la->lal[ind].ga >= 0)
146 /* Search the last entry in the linked list for this index */
147 ind_prev = ind;
148 while(ga2la->lal[ind_prev].next >= 0)
150 ind_prev = ga2la->lal[ind_prev].next;
152 /* Search for space in the array */
153 ind = ga2la->start_space_search;
154 while (ind < ga2la->nalloc && ga2la->lal[ind].ga >= 0)
156 ind++;
158 /* If we are a the end of the list we need to increase the size */
159 if (ind == ga2la->nalloc)
161 ga2la->nalloc = over_alloc_dd(ind+1);
162 srenew(ga2la->lal,ga2la->nalloc);
163 for(i=ind; i<ga2la->nalloc; i++)
165 ga2la->lal[i].ga = -1;
166 ga2la->lal[i].next = -1;
169 ga2la->lal[ind_prev].next = ind;
171 ga2la->start_space_search = ind + 1;
173 ga2la->lal[ind].ga = a_gl;
174 ga2la->lal[ind].la = a_loc;
175 ga2la->lal[ind].cell = cell;
178 /* Delete the ga2la entry for global atom a_gl */
179 static void ga2la_del(gmx_ga2la_t ga2la,int a_gl)
181 int ind,ind_prev;
183 if (ga2la->bAll)
185 ga2la->laa[a_gl].cell = -1;
187 return;
190 ind_prev = -1;
191 ind = a_gl % ga2la->mod;
194 if (ga2la->lal[ind].ga == a_gl)
196 if (ind_prev >= 0)
198 ga2la->lal[ind_prev].next = ga2la->lal[ind].next;
200 /* This index is a linked entry, so we free an entry.
201 * Check if we are creating the first empty space.
203 if (ind < ga2la->start_space_search)
205 ga2la->start_space_search = ind;
208 ga2la->lal[ind].ga = -1;
209 ga2la->lal[ind].cell = -1;
210 ga2la->lal[ind].next = -1;
212 return;
214 ind_prev = ind;
215 ind = ga2la->lal[ind].next;
217 while (ind >= 0);
219 return;
222 /* Change the local atom for present ga2la entry for global atom a_gl */
223 static void ga2la_change_la(gmx_ga2la_t ga2la,int a_gl,int a_loc)
225 int ind;
227 if (ga2la->bAll)
229 ga2la->laa[a_gl].la = a_loc;
231 return;
234 ind = a_gl % ga2la->mod;
237 if (ga2la->lal[ind].ga == a_gl)
239 ga2la->lal[ind].la = a_loc;
241 return;
243 ind = ga2la->lal[ind].next;
245 while (ind >= 0);
247 return;
250 /* Returns if the global atom a_gl available locally.
251 * Sets the local atom and cell,
252 * cell can be larger than the number of zones,
253 * in which case it indicates that it is more than one cell away
254 * in zone cell - #zones.
256 static bool ga2la_get(const gmx_ga2la_t ga2la,int a_gl,int *a_loc,int *cell)
258 int ind;
260 if (ga2la->bAll)
262 *a_loc = ga2la->laa[a_gl].la;
263 *cell = ga2la->laa[a_gl].cell;
265 return (ga2la->laa[a_gl].cell >= 0);
268 ind = a_gl % ga2la->mod;
271 if (ga2la->lal[ind].ga == a_gl)
273 *a_loc = ga2la->lal[ind].la;
274 *cell = ga2la->lal[ind].cell;
276 return TRUE;
278 ind = ga2la->lal[ind].next;
280 while (ind >= 0);
282 return FALSE;
285 /* Returns if the global atom a_gl is a home atom.
286 * Sets the local atom.
288 static bool ga2la_home(const gmx_ga2la_t ga2la,int a_gl,int *a_loc)
290 int ind;
292 if (ga2la->bAll)
294 *a_loc = ga2la->laa[a_gl].la;
296 return (ga2la->laa[a_gl].cell == 0);
299 ind = a_gl % ga2la->mod;
302 if (ga2la->lal[ind].ga == a_gl)
304 if (ga2la->lal[ind].cell == 0)
306 *a_loc = ga2la->lal[ind].la;
308 return TRUE;
310 else
312 return FALSE;
315 ind = ga2la->lal[ind].next;
317 while (ind >= 0);
319 return FALSE;
322 #endif /* _gmx_ga2la_h */