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