Code beautification with uncrustify
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_internal.h
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustr
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_internal_h
37 #define _nbnxn_internal_h
38
39 #include "typedefs.h"
40 #include "domdec.h"
41 #include "gmx_cyclecounter.h"
42
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
46
47
48 #ifdef GMX_X86_SSE2
49 /* Use 4-way SIMD for, always, single precision bounding box calculations */
50 #define NBNXN_SEARCH_BB_SSE
51 #endif
52
53
54 #ifdef GMX_NBNXN_SIMD
55 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
56 #define NBNXN_MEM_ALIGN  (GMX_NBNXN_SIMD_BITWIDTH/8)
57 #else
58 /* No alignment required, but set it so we can call the same routines */
59 #define NBNXN_MEM_ALIGN  32
60 #endif
61
62
63 /* A pair-search grid struct for one domain decomposition zone */
64 typedef struct {
65     rvec     c0;               /* The lower corner of the (local) grid        */
66     rvec     c1;               /* The upper corner of the (local) grid        */
67     real     atom_density;     /* The atom number density for the local grid  */
68
69     gmx_bool bSimple;          /* Is this grid simple or super/sub            */
70     int      na_c;             /* Number of atoms per cluster                 */
71     int      na_cj;            /* Number of atoms for list j-clusters         */
72     int      na_sc;            /* Number of atoms per super-cluster           */
73     int      na_c_2log;        /* 2log of na_c                                */
74
75     int      ncx;              /* Number of (super-)cells along x             */
76     int      ncy;              /* Number of (super-)cells along y             */
77     int      nc;               /* Total number of (super-)cells               */
78
79     real     sx;               /* x-size of a (super-)cell                    */
80     real     sy;               /* y-size of a (super-)cell                    */
81     real     inv_sx;           /* 1/sx                                        */
82     real     inv_sy;           /* 1/sy                                        */
83
84     int      cell0;            /* Index in nbs->cell corresponding to cell 0  */
85
86     int     *cxy_na;           /* The number of atoms for each column in x,y  */
87     int     *cxy_ind;          /* Grid (super)cell index, offset from cell0   */
88     int      cxy_nalloc;       /* Allocation size for cxy_na and cxy_ind      */
89
90     int     *nsubc;            /* The number of sub cells for each super cell */
91     float   *bbcz;             /* Bounding boxes in z for the super cells     */
92     float   *bb;               /* 3D bounding boxes for the sub cells         */
93     float   *bbj;              /* 3D j-b.boxes for SSE-double or AVX-single   */
94     int     *flags;            /* Flag for the super cells                    */
95     int      nc_nalloc;        /* Allocation size for the pointers above      */
96
97     float   *bbcz_simple;      /* bbcz for simple grid converted from super   */
98     float   *bb_simple;        /* bb for simple grid converted from super     */
99     int     *flags_simple;     /* flags for simple grid converted from super  */
100     int      nc_nalloc_simple; /* Allocation size for the pointers above   */
101
102     int      nsubc_tot;        /* Total number of subcell, used for printing  */
103 } nbnxn_grid_t;
104
105 #ifdef GMX_NBNXN_SIMD
106 #if GMX_NBNXN_SIMD_BITWIDTH == 128
107 #define GMX_MM128_HERE
108 #else
109 #if GMX_NBNXN_SIMD_BITWIDTH == 256
110 #define GMX_MM256_HERE
111 #else
112 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
113 #endif
114 #endif
115 #include "gmx_simd_macros.h"
116
117 typedef struct nbnxn_x_ci_simd_4xn {
118     /* The i-cluster coordinates for simple search */
119     gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
120     gmx_mm_pr ix_SSE1, iy_SSE1, iz_SSE1;
121     gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
122     gmx_mm_pr ix_SSE3, iy_SSE3, iz_SSE3;
123 } nbnxn_x_ci_simd_4xn_t;
124
125 typedef struct nbnxn_x_ci_simd_2xnn {
126     /* The i-cluster coordinates for simple search */
127     gmx_mm_pr ix_SSE0, iy_SSE0, iz_SSE0;
128     gmx_mm_pr ix_SSE2, iy_SSE2, iz_SSE2;
129 } nbnxn_x_ci_simd_2xnn_t;
130
131 #endif
132
133 /* Working data for the actual i-supercell during pair search */
134 typedef struct nbnxn_list_work {
135     gmx_cache_protect_t     cp0;   /* Protect cache between threads               */
136
137     float                  *bb_ci; /* The bounding boxes, pbc shifted, for each cluster */
138     real                   *x_ci;  /* The coordinates, pbc shifted, for each atom       */
139 #ifdef GMX_NBNXN_SIMD
140     nbnxn_x_ci_simd_4xn_t  *x_ci_simd_4xn;
141     nbnxn_x_ci_simd_2xnn_t *x_ci_simd_2xnn;
142 #endif
143     int                     cj_ind;    /* The current cj_ind index for the current list     */
144     int                     cj4_init;  /* The first unitialized cj4 block                   */
145
146     float                  *d2;        /* Bounding box distance work array                  */
147
148     nbnxn_cj_t             *cj;        /* The j-cell list                                   */
149     int                     cj_nalloc; /* Allocation size of cj                             */
150
151     int                     ncj_noq;   /* Nr. of cluster pairs without Coul for flop count  */
152     int                     ncj_hlj;   /* Nr. of cluster pairs with 1/2 LJ for flop count   */
153
154     gmx_cache_protect_t     cp1;       /* Protect cache between threads               */
155 } nbnxn_list_work_t;
156
157 /* Function type for setting the i-atom coordinate working data */
158 typedef void
159     gmx_icell_set_x_t (int ci,
160                        real shx, real shy, real shz,
161                        int na_c,
162                        int stride, const real *x,
163                        nbnxn_list_work_t *work);
164
165 static gmx_icell_set_x_t icell_set_x_simple;
166 #ifdef GMX_NBNXN_SIMD
167 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
168 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
169 #endif
170 static gmx_icell_set_x_t icell_set_x_supersub;
171 #ifdef NBNXN_SEARCH_SSE
172 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
173 #endif
174
175 #undef GMX_MM128_HERE
176 #undef GMX_MM256_HERE
177
178 /* Local cycle count struct for profiling */
179 typedef struct {
180     int          count;
181     gmx_cycles_t c;
182     gmx_cycles_t start;
183 } nbnxn_cycle_t;
184
185 /* Local cycle count enum for profiling */
186 enum {
187     enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
188 };
189
190 /* Thread-local work struct, contains part of nbnxn_grid_t */
191 typedef struct {
192     gmx_cache_protect_t  cp0;
193
194     int                 *cxy_na;
195     int                  cxy_na_nalloc;
196
197     int                 *sort_work;
198     int                  sort_work_nalloc;
199
200     nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
201
202     int                  ndistc;       /* Number of distance checks for flop counting */
203
204     nbnxn_cycle_t        cc[enbsCCnr];
205
206     gmx_cache_protect_t  cp1;
207 } nbnxn_search_work_t;
208
209 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
210 typedef struct nbnxn_search {
211     int                 ePBC;            /* PBC type enum                              */
212     matrix              box;             /* The periodic unit-cell                     */
213
214     gmx_bool            DomDec;          /* Are we doing domain decomposition?         */
215     ivec                dd_dim;          /* Are we doing DD in x,y,z?                  */
216     gmx_domdec_zones_t *zones;           /* The domain decomposition zones        */
217
218     int                 ngrid;           /* The number of grids, equal to #DD-zones    */
219     nbnxn_grid_t       *grid;            /* Array of grids, size ngrid                 */
220     int                *cell;            /* Actual allocated cell array for all grids  */
221     int                 cell_nalloc;     /* Allocation size of cell                    */
222     int                *a;               /* Atom index for grid, the inverse of cell   */
223     int                 a_nalloc;        /* Allocation size of a                       */
224
225     int                 natoms_local;    /* The local atoms run from 0 to natoms_local */
226     int                 natoms_nonlocal; /* The non-local atoms run from natoms_local
227                                           * to natoms_nonlocal */
228
229     gmx_bool             print_cycles;
230     int                  search_count;
231     nbnxn_cycle_t        cc[enbsCCnr];
232
233     gmx_icell_set_x_t   *icell_set_x; /* Function for setting i-coords    */
234
235     int                  nthread_max; /* Maximum number of threads for pair-search  */
236     nbnxn_search_work_t *work;        /* Work array, size nthread_max          */
237 } nbnxn_search_t_t;
238
239
240 static void nbs_cycle_start(nbnxn_cycle_t *cc)
241 {
242     cc->start = gmx_cycles_read();
243 }
244
245 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
246 {
247     cc->c += gmx_cycles_read() - cc->start;
248     cc->count++;
249 }
250
251
252 #ifdef __cplusplus
253 }
254 #endif
255
256 #endif