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