1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustr
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifndef _nbnxn_search_h
37 #define _nsnxn_search_h
45 /* Tells if the pair-list corresponding to nb_kernel_type is simple.
46 * Returns FALSE for super-sub type pair-list.
48 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type);
50 /* Due to the cluster size the effective pair-list is longer than
51 * that of a simple atom pair-list. This function gives the extra distance.
53 real nbnxn_get_rlist_effective_inc(int cluster_size,real atom_density);
55 /* Allocates and initializes a pair search data structure */
56 void nbnxn_init_search(nbnxn_search_t * nbs_ptr,
58 gmx_domdec_zones_t *zones,
61 /* Put the atoms on the pair search grid.
62 * Only atoms a0 to a1 in x are put on the grid.
63 * The atom_density is used to determine the grid size.
64 * When atom_density=-1, the density is determined from a1-a0 and the corners.
65 * With domain decomposition part of the n particles might have migrated,
66 * but have not been removed yet. This count is given by nmoved.
67 * When move[i] < 0 particle i has migrated and will not be put on the grid.
68 * Without domain decomposition move will be NULL.
70 void nbnxn_put_on_grid(nbnxn_search_t nbs,
73 rvec corner0,rvec corner1,
80 nbnxn_atomdata_t *nbat);
82 /* As nbnxn_put_on_grid, but for the non-local atoms
83 * with domain decomposition. Should be called after calling
84 * nbnxn_search_put_on_grid for the local atoms / home zone.
86 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t nbs,
87 const gmx_domdec_zones_t *zones,
91 nbnxn_atomdata_t *nbat);
93 /* Add simple grid type information to the local super/sub grid */
94 void nbnxn_grid_add_simple(nbnxn_search_t nbs,
95 nbnxn_atomdata_t *nbat);
97 /* Return the number of x and y cells in the local grid */
98 void nbnxn_get_ncells(nbnxn_search_t nbs,int *ncx,int *ncy);
100 /* Return the order indices *a of the atoms on the ns grid, size n */
101 void nbnxn_get_atomorder(nbnxn_search_t nbs,int **a,int *n);
103 /* Renumber the atom indices on the grid to consecutive order */
104 void nbnxn_set_atomorder(nbnxn_search_t nbs);
106 /* Initializes a set of pair lists stored in nbnxn_pairlist_set_t */
107 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
108 gmx_bool simple, gmx_bool combined,
109 gmx_nbat_alloc_t *alloc,
110 gmx_nbat_free_t *free);
112 /* Make a apir-list with radius rlist, store it in nbl.
113 * The parameter min_ci_balanced sets the minimum required
114 * number or roughly equally sized ci blocks in nbl.
115 * When set >0 ci lists will be chopped up when the estimate
116 * for the number of equally sized lists is below min_ci_balanced.
118 void nbnxn_make_pairlist(const nbnxn_search_t nbs,
119 const nbnxn_atomdata_t *nbat,
120 const t_blocka *excl,
123 nbnxn_pairlist_set_t *nbl_list,
128 /* Initialize the non-bonded atom data structure.
129 * The enum for nbatXFormat is in the file defining nbnxn_atomdata_t.
130 * Copy the ntypes*ntypes*2 sized nbfp non-bonded parameter list
131 * to the atom data structure.
133 void nbnxn_atomdata_init(FILE *fp,
134 nbnxn_atomdata_t *nbat,
136 int ntype,const real *nbfp,
139 gmx_nbat_alloc_t *alloc,
140 gmx_nbat_free_t *free);
142 /* Copy the atom data to the non-bonded atom data structure */
143 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
145 const nbnxn_search_t nbs,
146 const t_mdatoms *mdatoms,
149 /* Copy the shift vectors to nbat */
150 void nbnxn_atomdata_copy_shiftvec(gmx_bool dynamic_box,
152 nbnxn_atomdata_t *nbat);
154 /* Copy x to nbat->x.
155 * FillLocal tells if the local filler particle coordinates should be zeroed.
157 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
161 nbnxn_atomdata_t *nbat);
163 /* Add the forces stored in nbat to f, zeros the forces in nbat */
164 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
166 const nbnxn_atomdata_t *nbat,
169 /* Add the fshift force stored in nbat to fshift */
170 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,