Remove no-inline-max-size and suppress remark
[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 #include "nblist.h"
40
41 #ifdef __cplusplus
42 extern "C" {
43 #endif
44
45 /* A buffer data structure of 64 bytes
46  * to be placed at the beginning and end of structs
47  * to avoid cache invalidation of the real contents
48  * of the struct by writes to neighboring memory.
49  */
50 typedef struct {
51     int dummy[16];
52 } gmx_cache_protect_t;
53
54 /* Abstract type for pair searching data */
55 typedef struct nbnxn_search * nbnxn_search_t;
56
57 /* Function that should return a pointer *ptr to memory
58  * of size nbytes.
59  * Error handling should be done within this function.
60  */
61 typedef void nbnxn_alloc_t (void **ptr, size_t nbytes);
62
63 /* Function that should free the memory pointed to by *ptr.
64  * NULL should not be passed to this function.
65  */
66 typedef void nbnxn_free_t (void *ptr);
67
68 /* This is the actual cluster-pair list j-entry.
69  * cj is the j-cluster.
70  * The interaction bits in excl are indexed i-major, j-minor.
71  * The cj entries are sorted such that ones with exclusions come first.
72  * This means that once a full mask (=NBNXN_INTERACTION_MASK_ALL)
73  * is found, all subsequent j-entries in the i-entry also have full masks.
74  */
75 typedef struct {
76     int          cj;    /* The j-cluster                    */
77     unsigned int excl;  /* The exclusion (interaction) bits */
78     /* Indices into the arrays of SIMD interaction masks. */
79     char         interaction_mask_indices[4];
80 } nbnxn_cj_t;
81
82 /* In nbnxn_ci_t the integer shift contains the shift in the lower 7 bits.
83  * The upper bits contain information for non-bonded kernel optimization.
84  * Simply calculating LJ and Coulomb for all pairs in a cluster pair is fine.
85  * But three flags can be used to skip interactions, currently only for subc=0
86  * !(shift & NBNXN_CI_DO_LJ(subc))   => we can skip LJ for all pairs
87  * shift & NBNXN_CI_HALF_LJ(subc)    => we can skip LJ for the second half of i
88  * !(shift & NBNXN_CI_DO_COUL(subc)) => we can skip Coulomb for all pairs
89  */
90 #define NBNXN_CI_SHIFT          127
91 #define NBNXN_CI_DO_LJ(subc)    (1<<(7+3*(subc)))
92 #define NBNXN_CI_HALF_LJ(subc)  (1<<(8+3*(subc)))
93 #define NBNXN_CI_DO_COUL(subc)  (1<<(9+3*(subc)))
94
95 /* Simple pair-list i-unit */
96 typedef struct {
97     int ci;             /* i-cluster             */
98     int shift;          /* Shift vector index plus possible flags, see above */
99     int cj_ind_start;   /* Start index into cj   */
100     int cj_ind_end;     /* End index into cj     */
101 } nbnxn_ci_t;
102
103 /* Grouped pair-list i-unit */
104 typedef struct {
105     int sci;            /* i-super-cluster       */
106     int shift;          /* Shift vector index plus possible flags */
107     int cj4_ind_start;  /* Start index into cj4  */
108     int cj4_ind_end;    /* End index into cj4    */
109 } nbnxn_sci_t;
110
111 typedef struct {
112     unsigned int imask;    /* The i-cluster interactions mask for 1 warp  */
113     int          excl_ind; /* Index into the exclusion array for 1 warp   */
114 } nbnxn_im_ei_t;
115
116 typedef struct {
117     int           cj[4];   /* The 4 j-clusters                            */
118     nbnxn_im_ei_t imei[2]; /* The i-cluster mask data       for 2 warps   */
119 } nbnxn_cj4_t;
120
121 typedef struct {
122     unsigned int pair[32]; /* Topology exclusion interaction bits for one warp,
123                             * each unsigned has bitS for 4*8 i clusters
124                             */
125 } nbnxn_excl_t;
126
127 typedef struct {
128     gmx_cache_protect_t cp0;
129
130     nbnxn_alloc_t      *alloc;
131     nbnxn_free_t       *free;
132
133     gmx_bool            bSimple;         /* Simple list has na_sc=na_s and uses cj   *
134                                           * Complex list uses cj4                    */
135
136     int                     na_ci;       /* The number of atoms per i-cluster        */
137     int                     na_cj;       /* The number of atoms per j-cluster        */
138     int                     na_sc;       /* The number of atoms per super cluster    */
139     real                    rlist;       /* The radius for constructing the list     */
140     int                     nci;         /* The number of i-clusters in the list     */
141     nbnxn_ci_t             *ci;          /* The i-cluster list, size nci             */
142     int                     ci_nalloc;   /* The allocation size of ci                */
143     int                     nsci;        /* The number of i-super-clusters in the list */
144     nbnxn_sci_t            *sci;         /* The i-super-cluster list                 */
145     int                     sci_nalloc;  /* The allocation size of sci               */
146
147     int                     ncj;         /* The number of j-clusters in the list     */
148     nbnxn_cj_t             *cj;          /* The j-cluster list, size ncj             */
149     int                     cj_nalloc;   /* The allocation size of cj                */
150
151     int                     ncj4;        /* The total number of 4*j clusters         */
152     nbnxn_cj4_t            *cj4;         /* The 4*j cluster list, size ncj4          */
153     int                     cj4_nalloc;  /* The allocation size of cj4               */
154     int                     nexcl;       /* The count for excl                       */
155     nbnxn_excl_t           *excl;        /* Atom interaction bits (non-exclusions)   */
156     int                     excl_nalloc; /* The allocation size for excl             */
157     int                     nci_tot;     /* The total number of i clusters           */
158
159     struct nbnxn_list_work *work;
160
161     gmx_cache_protect_t     cp1;
162 } nbnxn_pairlist_t;
163
164 typedef struct {
165     int                nnbl;        /* number of lists */
166     nbnxn_pairlist_t **nbl;         /* lists */
167     gmx_bool           bCombined;   /* TRUE if lists get combined into one (the 1st) */
168     gmx_bool           bSimple;     /* TRUE if the list of of type "simple"
169                                        (na_sc=na_s, no super-clusters used) */
170     int                natpair_ljq; /* Total number of atom pairs for LJ+Q kernel */
171     int                natpair_lj;  /* Total number of atom pairs for LJ kernel   */
172     int                natpair_q;   /* Total number of atom pairs for Q kernel    */
173     t_nblist         **nbl_fep;
174 } nbnxn_pairlist_set_t;
175
176 enum {
177     nbatXYZ, nbatXYZQ, nbatX4, nbatX8
178 };
179
180 typedef struct {
181     real *f;      /* f, size natoms*fstride                             */
182     real *fshift; /* Shift force array, size SHIFTS*DIM                 */
183     int   nV;     /* The size of *Vvdw and *Vc                          */
184     real *Vvdw;   /* Temporary Van der Waals group energy storage       */
185     real *Vc;     /* Temporary Coulomb group energy storage             */
186     int   nVS;    /* The size of *VSvdw and *VSc                        */
187     real *VSvdw;  /* Temporary SIMD Van der Waals group energy storage  */
188     real *VSc;    /* Temporary SIMD Coulomb group energy storage        */
189 } nbnxn_atomdata_output_t;
190
191 /* Block size in atoms for the non-bonded thread force-buffer reduction,
192  * should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8).
193  * Should be small to reduce the reduction and zeroing cost,
194  * but too small will result in overhead.
195  * Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.
196  */
197 #ifdef GMX_DOUBLE
198 #define NBNXN_BUFFERFLAG_SIZE   8
199 #else
200 #define NBNXN_BUFFERFLAG_SIZE  16
201 #endif
202
203 /* We currently store the reduction flags as bits in an unsigned int.
204  * In most cases this limits the number of flags to 32.
205  * The reduction will automatically disable the flagging and do a full
206  * reduction when the flags won't fit, but this will lead to very slow
207  * reduction. As we anyhow don't expect reasonable performance with
208  * more than 32 threads, we put in this hard limit.
209  * You can increase this number, but the reduction will be very slow.
210  */
211 #define NBNXN_BUFFERFLAG_MAX_THREADS  32
212
213 /* Flags for telling if threads write to force output buffers */
214 typedef struct {
215     int           nflag;       /* The number of flag blocks                         */
216     unsigned int *flag;        /* Bit i is set when thread i writes to a cell-block */
217     int           flag_nalloc; /* Allocation size of cxy_flag                       */
218 } nbnxn_buffer_flags_t;
219
220 /* LJ combination rules: geometric, Lorentz-Berthelot, none */
221 enum {
222     ljcrGEOM, ljcrLB, ljcrNONE, ljcrNR
223 };
224
225 /* TODO: Remove need for forward declare */
226 struct tMPI_Atomic;
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 int            *simd_exclusion_filter1;
264     unsigned int            *simd_exclusion_filter2;
265     real                    *simd_interaction_array; /* Array of masks needed for exclusions on QPX */
266     int                      nout;                   /* The number of force arrays                         */
267     nbnxn_atomdata_output_t *out;                    /* Output data structures               */
268     int                      nalloc;                 /* Allocation size of all arrays (for x/f *x/fstride) */
269     gmx_bool                 bUseBufferFlags;        /* Use the flags or operate on all atoms     */
270     nbnxn_buffer_flags_t     buffer_flags;           /* Flags for buffer zeroing+reduc.  */
271     gmx_bool                 bUseTreeReduce;         /* Use tree for force reduction */
272     struct tMPI_Atomic      *syncStep;               /* Synchronization step for tree reduce */
273 } nbnxn_atomdata_t;
274
275 #ifdef __cplusplus
276 }
277 #endif
278
279 #endif