1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustr
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
36 #ifndef _nbnxn_internal_h
37 #define _nbnxn_internal_h
41 #include "gmx_cyclecounter.h"
44 /* The include below sets the SIMD instruction type (precision+width)
45 * for all nbnxn SIMD search and non-bonded kernel code.
47 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
48 #define GMX_USE_HALF_WIDTH_SIMD_HERE
50 #include "gmx_simd_macros.h"
59 /* Use 4-way SIMD for, always, single precision bounding box calculations */
60 #define NBNXN_SEARCH_BB_SSE
65 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
66 #define NBNXN_MEM_ALIGN (GMX_SIMD_WIDTH_HERE*sizeof(real))
68 /* No alignment required, but set it so we can call the same routines */
69 #define NBNXN_MEM_ALIGN 32
73 /* A pair-search grid struct for one domain decomposition zone */
75 rvec c0; /* The lower corner of the (local) grid */
76 rvec c1; /* The upper corner of the (local) grid */
77 real atom_density; /* The atom number density for the local grid */
79 gmx_bool bSimple; /* Is this grid simple or super/sub */
80 int na_c; /* Number of atoms per cluster */
81 int na_cj; /* Number of atoms for list j-clusters */
82 int na_sc; /* Number of atoms per super-cluster */
83 int na_c_2log; /* 2log of na_c */
85 int ncx; /* Number of (super-)cells along x */
86 int ncy; /* Number of (super-)cells along y */
87 int nc; /* Total number of (super-)cells */
89 real sx; /* x-size of a (super-)cell */
90 real sy; /* y-size of a (super-)cell */
91 real inv_sx; /* 1/sx */
92 real inv_sy; /* 1/sy */
94 int cell0; /* Index in nbs->cell corresponding to cell 0 */
96 int *cxy_na; /* The number of atoms for each column in x,y */
97 int *cxy_ind; /* Grid (super)cell index, offset from cell0 */
98 int cxy_nalloc; /* Allocation size for cxy_na and cxy_ind */
100 int *nsubc; /* The number of sub cells for each super cell */
101 float *bbcz; /* Bounding boxes in z for the super cells */
102 float *bb; /* 3D bounding boxes for the sub cells */
103 float *bbj; /* 3D j-b.boxes for SSE-double or AVX-single */
104 int *flags; /* Flag for the super cells */
105 int nc_nalloc; /* Allocation size for the pointers above */
107 float *bbcz_simple; /* bbcz for simple grid converted from super */
108 float *bb_simple; /* bb for simple grid converted from super */
109 int *flags_simple; /* flags for simple grid converted from super */
110 int nc_nalloc_simple; /* Allocation size for the pointers above */
112 int nsubc_tot; /* Total number of subcell, used for printing */
115 #ifdef GMX_NBNXN_SIMD
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;
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;
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 */
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;
143 int cj_ind; /* The current cj_ind index for the current list */
144 int cj4_init; /* The first unitialized cj4 block */
146 float *d2; /* Bounding box distance work array */
148 nbnxn_cj_t *cj; /* The j-cell list */
149 int cj_nalloc; /* Allocation size of cj */
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 */
154 int *sort; /* Sort index */
155 int sort_nalloc; /* Allocation size of sort */
157 nbnxn_sci_t *sci_sort; /* Second sci array, for sorting */
158 int sci_sort_nalloc; /* Allocation size of sci_sort */
160 gmx_cache_protect_t cp1; /* Protect cache between threads */
163 /* Function type for setting the i-atom coordinate working data */
165 gmx_icell_set_x_t (int ci,
166 real shx, real shy, real shz,
168 int stride, const real *x,
169 nbnxn_list_work_t *work);
171 static gmx_icell_set_x_t icell_set_x_simple;
172 #ifdef GMX_NBNXN_SIMD
173 static gmx_icell_set_x_t icell_set_x_simple_simd_4xn;
174 static gmx_icell_set_x_t icell_set_x_simple_simd_2xnn;
176 static gmx_icell_set_x_t icell_set_x_supersub;
177 #ifdef NBNXN_SEARCH_SSE
178 static gmx_icell_set_x_t icell_set_x_supersub_sse8;
181 /* Local cycle count struct for profiling */
188 /* Local cycle count enum for profiling */
190 enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr
193 /* Thread-local work struct, contains part of nbnxn_grid_t */
195 gmx_cache_protect_t cp0;
201 int sort_work_nalloc;
203 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
205 int ndistc; /* Number of distance checks for flop counting */
207 nbnxn_cycle_t cc[enbsCCnr];
209 gmx_cache_protect_t cp1;
210 } nbnxn_search_work_t;
212 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
213 typedef struct nbnxn_search {
214 int ePBC; /* PBC type enum */
215 matrix box; /* The periodic unit-cell */
217 gmx_bool DomDec; /* Are we doing domain decomposition? */
218 ivec dd_dim; /* Are we doing DD in x,y,z? */
219 gmx_domdec_zones_t *zones; /* The domain decomposition zones */
221 int ngrid; /* The number of grids, equal to #DD-zones */
222 nbnxn_grid_t *grid; /* Array of grids, size ngrid */
223 int *cell; /* Actual allocated cell array for all grids */
224 int cell_nalloc; /* Allocation size of cell */
225 int *a; /* Atom index for grid, the inverse of cell */
226 int a_nalloc; /* Allocation size of a */
228 int natoms_local; /* The local atoms run from 0 to natoms_local */
229 int natoms_nonlocal; /* The non-local atoms run from natoms_local
230 * to natoms_nonlocal */
232 gmx_bool print_cycles;
234 nbnxn_cycle_t cc[enbsCCnr];
236 gmx_icell_set_x_t *icell_set_x; /* Function for setting i-coords */
238 int nthread_max; /* Maximum number of threads for pair-search */
239 nbnxn_search_work_t *work; /* Work array, size nthread_max */
243 static void nbs_cycle_start(nbnxn_cycle_t *cc)
245 cc->start = gmx_cycles_read();
248 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
250 cc->c += gmx_cycles_read() - cc->start;