413af358019c353bd972b3573bb6aad13a163ad5
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_pairlist.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, 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 _nbnxn_pairlist_h
37 #define _nbnxn_pairlist_h
38
39 #include <stddef.h>
40
41 #include "thread_mpi/atomic.h"
42
43 #include "gromacs/legacyheaders/types/nblist.h"
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
47
48 #ifdef __cplusplus
49 extern "C" {
50 #endif
51
52 /* A buffer data structure of 64 bytes
53  * to be placed at the beginning and end of structs
54  * to avoid cache invalidation of the real contents
55  * of the struct by writes to neighboring memory.
56  */
57 typedef struct {
58     int dummy[16];
59 } gmx_cache_protect_t;
60
61 /* Abstract type for pair searching data */
62 typedef struct nbnxn_search * nbnxn_search_t;
63
64 /* Function that should return a pointer *ptr to memory
65  * of size nbytes.
66  * Error handling should be done within this function.
67  */
68 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
69
70 /* Function that should free the memory pointed to by *ptr.
71  * NULL should not be passed to this function.
72  */
73 typedef void nbnxn_free_t (void *ptr);
74
75 /* This is the actual cluster-pair list j-entry.
76  * cj is the j-cluster.
77  * The interaction bits in excl are indexed i-major, j-minor.
78  * The cj entries are sorted such that ones with exclusions come first.
79  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
80  * is found, all subsequent j-entries in the i-entry also have full masks.
81  */
82 typedef struct {
83     int          cj;    /* The j-cluster                    */
84     unsigned int excl;  /* The exclusion (interaction) bits */
85     /* Indices into the arrays of SIMD interaction masks. */
86     char         interaction_mask_indices[4];
87 } nbnxn_cj_t;
88
89 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
90  * The upper bits contain information for non-bonded kernel optimization.
91  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
92  * But three flags can be used to skip interactions, currently only for subc=0
93  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
94  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
95  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
96  */
97 #define NBNXN_CI_SHIFT          127
98 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
99 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
100 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
101
102 /* Simple pair-list i-unit */
103 typedef struct {
104     int ci;             /* i-cluster             */
105     int shift;          /* Shift vector index plus possible flags, see above */
106     int cj_ind_start;   /* Start index into cj   */
107     int cj_ind_end;     /* End index into cj     */
108 } nbnxn_ci_t;
109
110 /* Grouped pair-list i-unit */
111 typedef struct {
112     int sci;            /* i-super-cluster       */
113     int shift;          /* Shift vector index plus possible flags */
114     int cj4_ind_start;  /* Start index into cj4  */
115     int cj4_ind_end;    /* End index into cj4    */
116 } nbnxn_sci_t;
117
118 typedef struct {
119     unsigned int imask;    /* The i-cluster interactions mask for 1 warp  */
120     int          excl_ind; /* Index into the exclusion array for 1 warp   */
121 } nbnxn_im_ei_t;
122
123 typedef struct {
124     int           cj[4];   /* The 4 j-clusters                            */
125     nbnxn_im_ei_t imei[2]; /* The i-cluster mask data       for 2 warps   */
126 } nbnxn_cj4_t;
127
128 typedef struct {
129     unsigned int pair[32]; /* Topology exclusion interaction bits for one warp,
130                             * each unsigned has bitS for 4*8 i clusters
131                             */
132 } nbnxn_excl_t;
133
134 typedef struct nbnxn_pairlist_t {
135     gmx_cache_protect_t cp0;
136
137     nbnxn_alloc_t      *alloc;
138     nbnxn_free_t       *free;
139
140     gmx_bool            bSimple;         /* Simple list has na_sc=na_s and uses cj   *
141                                           * Complex list uses cj4                    */
142
143     int                     na_ci;       /* The number of atoms per i-cluster        */
144     int                     na_cj;       /* The number of atoms per j-cluster        */
145     int                     na_sc;       /* The number of atoms per super cluster    */
146     real                    rlist;       /* The radius for constructing the list     */
147     int                     nci;         /* The number of i-clusters in the list     */
148     nbnxn_ci_t             *ci;          /* The i-cluster list, size nci             */
149     int                     ci_nalloc;   /* The allocation size of ci                */
150     int                     nsci;        /* The number of i-super-clusters in the list */
151     nbnxn_sci_t            *sci;         /* The i-super-cluster list                 */
152     int                     sci_nalloc;  /* The allocation size of sci               */
153
154     int                     ncj;         /* The number of j-clusters in the list     */
155     nbnxn_cj_t             *cj;          /* The j-cluster list, size ncj             */
156     int                     cj_nalloc;   /* The allocation size of cj                */
157
158     int                     ncj4;        /* The total number of 4*j clusters         */
159     nbnxn_cj4_t            *cj4;         /* The 4*j cluster list, size ncj4          */
160     int                     cj4_nalloc;  /* The allocation size of cj4               */
161     int                     nexcl;       /* The count for excl                       */
162     nbnxn_excl_t           *excl;        /* Atom interaction bits (non-exclusions)   */
163     int                     excl_nalloc; /* The allocation size for excl             */
164     int                     nci_tot;     /* The total number of i clusters           */
165
166     struct nbnxn_list_work *work;
167
168     gmx_cache_protect_t     cp1;
169 } nbnxn_pairlist_t;
170
171 typedef struct {
172     int                nnbl;        /* number of lists */
173     nbnxn_pairlist_t **nbl;         /* lists */
174     gmx_bool           bCombined;   /* TRUE if lists get combined into one (the 1st) */
175     gmx_bool           bSimple;     /* TRUE if the list of of type "simple"
176                                        (na_sc=na_s, no super-clusters used) */
177     int                natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
178     int                natpair_lj;  /* Total number of atom pairs for LJ kernel   */
179     int                natpair_q;   /* Total number of atom pairs for Q kernel    */
180     t_nblist         **nbl_fep;
181 } nbnxn_pairlist_set_t;
182
183 enum {
184     nbatXYZ, nbatXYZQ, nbatX4, nbatX8
185 };
186
187 typedef struct {
188     real *f;      /* f, size natoms*fstride                             */
189     real *fshift; /* Shift force array, size SHIFTS*DIM                 */
190     int   nV;     /* The size of *Vvdw and *Vc                          */
191     real *Vvdw;   /* Temporary Van der Waals group energy storage       */
192     real *Vc;     /* Temporary Coulomb group energy storage             */
193     int   nVS;    /* The size of *VSvdw and *VSc                        */
194     real *VSvdw;  /* Temporary SIMD Van der Waals group energy storage  */
195     real *VSc;    /* Temporary SIMD Coulomb group energy storage        */
196 } nbnxn_atomdata_output_t;
197
198 /* Block size in atoms for the non-bonded thread force-buffer reduction,
199  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
200  * Should be small to reduce the reduction and zeroing cost,
201  * but too small will result in overhead.
202  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
203  */
204 #ifdef GMX_DOUBLE
205 #define NBNXN_BUFFERFLAG_SIZE   8
206 #else
207 #define NBNXN_BUFFERFLAG_SIZE  16
208 #endif
209
210 /* We currently store the reduction flags as bits in an unsigned int.
211  * In most cases this limits the number of flags to 32.
212  * The reduction will automatically disable the flagging and do a full
213  * reduction when the flags won't fit, but this will lead to very slow
214  * reduction. As we anyhow don't expect reasonable performance with
215  * more than 32 threads, we put in this hard limit.
216  * You can increase this number, but the reduction will be very slow.
217  */
218 #define NBNXN_BUFFERFLAG_MAX_THREADS  32
219
220 /* Flags for telling if threads write to force output buffers */
221 typedef struct {
222     int           nflag;       /* The number of flag blocks                         */
223     unsigned int *flag;        /* Bit i is set when thread i writes to a cell-block */
224     int           flag_nalloc; /* Allocation size of cxy_flag                       */
225 } nbnxn_buffer_flags_t;
226
227 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
228 enum {
229     ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
230 };
231
232 typedef struct nbnxn_atomdata_t {
233     nbnxn_alloc_t           *alloc;
234     nbnxn_free_t            *free;
235     int                      ntype;           /* The number of different atom types                 */
236     real                    *nbfp;            /* Lennard-Jones 6*C6 and 12*C12 params, size ntype^2*2 */
237     int                      comb_rule;       /* Combination rule, see enum above                   */
238     real                    *nbfp_comb;       /* LJ parameter per atom type, size ntype*2           */
239     real                    *nbfp_s4;         /* As nbfp, but with stride 4, size ntype^2*4. This
240                                                * might suit 4-wide SIMD loads of two values (e.g.
241                                                * two floats in single precision on x86).            */
242     int                      natoms;          /* Number of atoms                                    */
243     int                      natoms_local;    /* Number of local atoms                           */
244     int                     *type;            /* Atom types                                         */
245     real                    *lj_comb;         /* LJ parameters per atom for combining for pairs     */
246     int                      XFormat;         /* The format of x (and q), enum                      */
247     int                      FFormat;         /* The format of f, enum                              */
248     real                    *q;               /* Charges, can be NULL if incorporated in x          */
249     int                      na_c;            /* The number of atoms per cluster                    */
250     int                      nenergrp;        /* The number of energy groups                        */
251     int                      neg_2log;        /* Log2 of nenergrp                                   */
252     int                     *energrp;         /* The energy groups per cluster, can be NULL         */
253     gmx_bool                 bDynamicBox;     /* Do we need to update shift_vec every step?    */
254     rvec                    *shift_vec;       /* Shift vectors, copied from t_forcerec              */
255     int                      xstride;         /* stride for a coordinate in x (usually 3 or 4)      */
256     int                      fstride;         /* stride for a coordinate in f (usually 3 or 4)      */
257     real                    *x;               /* x and possibly q, size natoms*xstride              */
258
259     /* j-atom minus i-atom index for generating self and Newton exclusions
260      * cluster-cluster pairs of the diagonal, for 4xn and 2xnn kernels.
261      */
262     real                    *simd_4xn_diagonal_j_minus_i;
263     real                    *simd_2xnn_diagonal_j_minus_i;
264     /* Filters for topology exclusion masks for the SIMD kernels.
265      * filter2 is the same as filter1, but with each element duplicated.
266      */
267     unsigned int            *simd_exclusion_filter1;
268     unsigned int            *simd_exclusion_filter2;
269     real                    *simd_interaction_array; /* Array of masks needed for exclusions on QPX */
270     int                      nout;                   /* The number of force arrays                         */
271     nbnxn_atomdata_output_t *out;                    /* Output data structures               */
272     int                      nalloc;                 /* Allocation size of all arrays (for x/f *x/fstride) */
273     gmx_bool                 bUseBufferFlags;        /* Use the flags or operate on all atoms     */
274     nbnxn_buffer_flags_t     buffer_flags;           /* Flags for buffer zeroing+reduc.  */
275     gmx_bool                 bUseTreeReduce;         /* Use tree for force reduction */
276     tMPI_Atomic_t           *syncStep;               /* Synchronization step for tree reduce */
277 } nbnxn_atomdata_t;
278
279 #ifdef __cplusplus
280 }
281 #endif
282
283 #endif