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