Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / gpp_nextnb.cpp
blob704ff8a894ee0d8c7d98ee530073d9022a05f12c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2014,2015, 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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "gpp_nextnb.h"
42 #include <stdlib.h>
44 #include "gromacs/gmxpreprocess/toputil.h"
45 #include "gromacs/topology/ifunc.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
49 /* #define DEBUG_NNB */
51 typedef struct {
52 int ai, aj;
53 } sortable;
55 static int
56 bond_sort (const void *a, const void *b)
58 sortable *sa, *sb;
60 sa = (sortable *) a;
61 sb = (sortable *) b;
63 if (sa->ai == sb->ai)
65 return (sa->aj-sb->aj);
67 else
69 return (sa->ai-sb->ai);
73 static int
74 compare_int (const void * a, const void * b)
76 return ( *(int*)a - *(int*)b );
80 #ifdef DEBUG
81 #define prints(str, n, s) __prints(str, n, s)
82 static void __prints(char *str, int n, sortable *s)
84 int i;
86 if (debug)
88 fprintf(debug, "%s\n", str);
89 fprintf(debug, "Sortables \n");
90 for (i = 0; (i < n); i++)
92 fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
95 fflush(debug);
98 #else
99 #define prints(str, n, s)
100 #endif
102 void init_nnb(t_nextnb *nnb, int nr, int nrex)
104 int i;
106 /* initiate nnb */
107 nnb->nr = nr;
108 nnb->nrex = nrex;
110 snew(nnb->a, nr);
111 snew(nnb->nrexcl, nr);
112 for (i = 0; (i < nr); i++)
114 snew(nnb->a[i], nrex+1);
115 snew(nnb->nrexcl[i], nrex+1);
119 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
121 srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
122 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
123 nnb->nrexcl[i][nre]++;
126 void done_nnb (t_nextnb *nnb)
128 int i, nre;
130 for (i = 0; (i < nnb->nr); i++)
132 for (nre = 0; (nre <= nnb->nrex); nre++)
134 if (nnb->nrexcl[i][nre] > 0)
136 sfree (nnb->a[i][nre]);
139 sfree (nnb->nrexcl[i]);
140 sfree (nnb->a[i]);
142 sfree (nnb->a);
143 sfree (nnb->nrexcl);
144 nnb->nr = 0;
145 nnb->nrex = 0;
148 #ifdef DEBUG_NNB
149 void __print_nnb(t_nextnb *nnb, char *s)
151 int i, j, k;
153 if (debug)
155 fprintf(debug, "%s\n", s);
156 fprintf(debug, "nnb->nr: %d\n", nnb->nr);
157 fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
158 for (i = 0; (i < nnb->nr); i++)
160 for (j = 0; (j <= nnb->nrex); j++)
162 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
163 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
165 fprintf(debug, "%d, ", nnb->a[i][j][k]);
167 fprintf(debug, "\n");
172 #endif
174 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
176 int i, j, j_index;
177 int nre, nrx, nrs, nr_of_sortables;
178 sortable *s;
180 srenew(excl->index, nnb->nr+1);
181 excl->index[0] = 0;
182 for (i = 0; (i < nnb->nr); i++)
184 /* calculate the total number of exclusions for atom i */
185 nr_of_sortables = 0;
186 for (nre = 0; (nre <= nnb->nrex); nre++)
188 nr_of_sortables += nnb->nrexcl[i][nre];
191 if (debug)
193 fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
195 /* make space for sortable array */
196 snew(s, nr_of_sortables);
198 /* fill the sortable array and sort it */
199 nrs = 0;
200 for (nre = 0; (nre <= nnb->nrex); nre++)
202 for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
204 s[nrs].ai = i;
205 s[nrs].aj = nnb->a[i][nre][nrx];
206 nrs++;
209 if (nrs != nr_of_sortables)
211 gmx_incons("Generating exclusions");
213 prints("nnb2excl before qsort", nr_of_sortables, s);
214 if (nr_of_sortables > 1)
216 qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
217 prints("nnb2excl after qsort", nr_of_sortables, s);
220 /* remove duplicate entries from the list */
221 j_index = 0;
222 if (nr_of_sortables > 0)
224 for (j = 1; (j < nr_of_sortables); j++)
226 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
228 s[j_index++] = s[j-1];
231 s[j_index++] = s[j-1];
233 nr_of_sortables = j_index;
234 prints("after rm-double", j_index, s);
236 /* make space for arrays */
237 srenew(excl->a, excl->nra+nr_of_sortables);
239 /* put the sorted exclusions in the target list */
240 for (nrs = 0; (nrs < nr_of_sortables); nrs++)
242 excl->a[excl->nra+nrs] = s[nrs].aj;
244 excl->nra += nr_of_sortables;
245 excl->index[i+1] = excl->nra;
247 /* cleanup temporary space */
248 sfree (s);
250 if (debug)
252 print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
256 static void do_gen(int nrbonds, /* total number of bonds in s */
257 sortable *s, /* bidirectional list of bonds */
258 t_nextnb *nnb) /* the tmp storage for excl */
259 /* Assume excl is initalised and s[] contains all bonds bidirectional */
261 int i, j, k, n, nb;
263 /* exclude self */
264 for (i = 0; (i < nnb->nr); i++)
266 add_nnb(nnb, 0, i, i);
268 print_nnb(nnb, "After exclude self");
270 /* exclude all the bonded atoms */
271 if (nnb->nrex > 0)
273 for (i = 0; (i < nrbonds); i++)
275 add_nnb(nnb, 1, s[i].ai, s[i].aj);
278 print_nnb(nnb, "After exclude bonds");
280 /* for the nr of exclusions per atom */
281 for (n = 1; (n < nnb->nrex); n++)
283 /* now for all atoms */
284 for (i = 0; (i < nnb->nr); i++)
286 /* for all directly bonded atoms of atom i */
287 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
290 /* store the 1st neighbour in nb */
291 nb = nnb->a[i][1][j];
293 /* store all atoms in nb's n-th list into i's n+1-th list */
294 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
296 if (i != nnb->a[nb][n][k])
298 add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
304 print_nnb(nnb, "After exclude rest");
308 static void add_b(t_params *bonds, int *nrf, sortable *s)
310 int i;
311 int ai, aj;
313 for (i = 0; (i < bonds->nr); i++)
315 ai = bonds->param[i].ai();
316 aj = bonds->param[i].aj();
317 if ((ai < 0) || (aj < 0))
319 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
320 i, ai, aj);
322 /* Add every bond twice */
323 s[(*nrf)].ai = ai;
324 s[(*nrf)++].aj = aj;
325 s[(*nrf)].aj = ai;
326 s[(*nrf)++].ai = aj;
330 void gen_nnb(t_nextnb *nnb, t_params plist[])
332 sortable *s;
333 int i, nrbonds, nrf;
335 nrbonds = 0;
336 for (i = 0; (i < F_NRE); i++)
338 if (IS_CHEMBOND(i))
340 /* we need every bond twice (bidirectional) */
341 nrbonds += 2*plist[i].nr;
345 snew(s, nrbonds);
347 nrf = 0;
348 for (i = 0; (i < F_NRE); i++)
350 if (IS_CHEMBOND(i))
352 add_b(&plist[i], &nrf, s);
356 /* now sort the bonds */
357 prints("gen_excl before qsort", nrbonds, s);
358 if (nrbonds > 1)
360 qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
361 prints("gen_excl after qsort", nrbonds, s);
364 do_gen(nrbonds, s, nnb);
365 sfree(s);
368 static void
369 sort_and_purge_nnb(t_nextnb *nnb)
371 int i, j, k, m, n, cnt, found, prev, idx;
373 for (i = 0; (i < nnb->nr); i++)
375 for (n = 0; (n <= nnb->nrex); n++)
377 /* Sort atoms in this list */
378 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
380 cnt = 0;
381 prev = -1;
382 for (j = 0; j < nnb->nrexcl[i][n]; j++)
384 idx = nnb->a[i][n][j];
386 found = 0;
387 for (m = 0; m < n && !found; m++)
389 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
391 found = (idx == nnb->a[i][m][k]);
395 if (!found && nnb->a[i][n][j] != prev)
397 nnb->a[i][n][cnt] = nnb->a[i][n][j];
398 prev = nnb->a[i][n][cnt];
399 cnt++;
402 nnb->nrexcl[i][n] = cnt;
408 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
410 if (nrexcl < 0)
412 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
414 init_nnb(nnb, nratoms, nrexcl);
415 gen_nnb(nnb, plist);
416 excl->nr = nratoms;
417 sort_and_purge_nnb(nnb);
418 nnb2excl (nnb, excl);