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"
49 /* Use 4-way SIMD for, always, single precision bounding box calculations */
50 #define NBNXN_SEARCH_BB_SSE
55 /* Memory alignment in bytes as required by SIMD aligned loads/stores */
56 #define NBNXN_MEM_ALIGN (GMX_NBNXN_SIMD_BITWIDTH/8)
58 /* No alignment required, but set it so we can call the same routines */
59 #define NBNXN_MEM_ALIGN 32
63 /* A pair-search grid struct for one domain decomposition zone */
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 */
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 */
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 */
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 */
84 int cell0; /* Index in nbs->cell corresponding to cell 0 */
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 */
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 */
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 */
102 int nsubc_tot; /* Total number of subcell, used for printing */
105 #ifdef GMX_NBNXN_SIMD
106 #if GMX_NBNXN_SIMD_BITWIDTH == 128
107 #define GMX_MM128_HERE
109 #if GMX_NBNXN_SIMD_BITWIDTH == 256
110 #define GMX_MM256_HERE
112 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
115 #include "gmx_simd_macros.h"
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 gmx_cache_protect_t cp1; /* Protect cache between threads */
157 /* Function type for setting the i-atom coordinate working data */
159 gmx_icell_set_x_t(int ci,
160 real shx,real shy,real shz,
162 int stride,const real *x,
163 nbnxn_list_work_t *work);
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;
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;
175 #undef GMX_MM128_HERE
176 #undef GMX_MM256_HERE
178 /* Local cycle count struct for profiling */
185 /* Local cycle count enum for profiling */
186 enum { enbsCCgrid, enbsCCsearch, enbsCCcombine, enbsCCreducef, enbsCCnr };
188 /* Thread-local work struct, contains part of nbnxn_grid_t */
190 gmx_cache_protect_t cp0;
196 int sort_work_nalloc;
198 nbnxn_buffer_flags_t buffer_flags; /* Flags for force buffer access */
200 int ndistc; /* Number of distance checks for flop counting */
202 nbnxn_cycle_t cc[enbsCCnr];
204 gmx_cache_protect_t cp1;
205 } nbnxn_search_work_t;
207 /* Main pair-search struct, contains the grid(s), not the pair-list(s) */
208 typedef struct nbnxn_search {
209 int ePBC; /* PBC type enum */
210 matrix box; /* The periodic unit-cell */
212 gmx_bool DomDec; /* Are we doing domain decomposition? */
213 ivec dd_dim; /* Are we doing DD in x,y,z? */
214 gmx_domdec_zones_t *zones; /* The domain decomposition zones */
216 int ngrid; /* The number of grids, equal to #DD-zones */
217 nbnxn_grid_t *grid; /* Array of grids, size ngrid */
218 int *cell; /* Actual allocated cell array for all grids */
219 int cell_nalloc; /* Allocation size of cell */
220 int *a; /* Atom index for grid, the inverse of cell */
221 int a_nalloc; /* Allocation size of a */
223 int natoms_local; /* The local atoms run from 0 to natoms_local */
224 int natoms_nonlocal; /* The non-local atoms run from natoms_local
225 * to natoms_nonlocal */
227 gmx_bool print_cycles;
229 nbnxn_cycle_t cc[enbsCCnr];
231 gmx_icell_set_x_t *icell_set_x; /* Function for setting i-coords */
233 int nthread_max; /* Maximum number of threads for pair-search */
234 nbnxn_search_work_t *work; /* Work array, size nthread_max */
238 static void nbs_cycle_start(nbnxn_cycle_t *cc)
240 cc->start = gmx_cycles_read();
243 static void nbs_cycle_stop(nbnxn_cycle_t *cc)
245 cc->c += gmx_cycles_read() - cc->start;