0e92eb3cde7c16a92a96bdd624a8d2da3a603c4d
[alexxy/gromacs.git] / src / gromacs / nbnxm / pairlistset.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 #ifndef GMX_NBNXM_PAIRLISTSET_H
37 #define GMX_NBNXM_PAIRLISTSET_H
38
39 #include "gromacs/math/vectypes.h"
40 #include "gromacs/utility/basedefinitions.h"
41 #include "gromacs/utility/real.h"
42
43 struct gmx_domdec_zones_t;
44 struct gmx_groups_t;
45 struct nbnxn_atomdata_t;
46 struct nbnxn_pairlist_set_t;
47 struct nbnxn_search;
48 struct t_blocka;
49 struct t_nrnb;
50
51 /* Function that should return a pointer *ptr to memory
52  * of size nbytes.
53  * Error handling should be done within this function.
54  */
55 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
56
57 /* Function that should free the memory pointed to by *ptr.
58  * NULL should not be passed to this function.
59  */
60 typedef void nbnxn_free_t (void *ptr);
61
62 /* Tells if the pair-list corresponding to nb_kernel_type is simple.
63  * Returns FALSE for super-sub type pair-list.
64  */
65 gmx_bool nbnxn_kernel_pairlist_simple(int nb_kernel_type);
66
67 /* Due to the cluster size the effective pair-list is longer than
68  * that of a simple atom pair-list. This function gives the extra distance.
69  */
70 real nbnxn_get_rlist_effective_inc(int cluster_size, real atom_density);
71
72 /* Allocates and initializes a pair search data structure */
73 nbnxn_search *nbnxn_init_search(const ivec                *n_dd_cells,
74                                 const gmx_domdec_zones_t  *zones,
75                                 gmx_bool                   bFEP,
76                                 int                        nthread_max);
77
78 /* Initializes a set of pair lists stored in nbnxn_pairlist_set_t */
79 void nbnxn_init_pairlist_set(nbnxn_pairlist_set_t *nbl_list,
80                              gmx_bool simple, gmx_bool combined);
81
82 /* Make a pair-list with radius rlist, store it in nbl.
83  * The parameter min_ci_balanced sets the minimum required
84  * number or roughly equally sized ci blocks in nbl.
85  * When set >0 ci lists will be chopped up when the estimate
86  * for the number of equally sized lists is below min_ci_balanced.
87  * With perturbed particles, also a group scheme style nbl_fep list is made.
88  */
89 void nbnxn_make_pairlist(nbnxn_search         *nbs,
90                          nbnxn_atomdata_t     *nbat,
91                          const t_blocka       *excl,
92                          real                  rlist,
93                          int                   min_ci_balanced,
94                          nbnxn_pairlist_set_t *nbl_list,
95                          int                   iloc,
96                          int                   nb_kernel_type,
97                          t_nrnb               *nrnb);
98
99 /*! \brief Prepare the list-set produced by the search for dynamic pruning
100  *
101  * \param[in,out] listSet  The list-set to prepare for dynamic pruning.
102  */
103 void nbnxnPrepareListForDynamicPruning(nbnxn_pairlist_set_t *listSet);
104
105 #endif