Change nbnxn grid reallocation
[gromacs.git] / src / gromacs / mdlib / nbnxn_atomdata.cpp
blobfa5a923cb26a435714d4338fbc48f4766c00a5df
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
36 #include "gmxpre.h"
38 #include "nbnxn_atomdata.h"
40 #include <assert.h>
41 #include <stdlib.h>
42 #include <string.h>
44 #include <cmath>
46 #include <algorithm>
48 #include "thread_mpi/atomic.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdlib/nb_verlet.h"
55 #include "gromacs/mdlib/nbnxn_consts.h"
56 #include "gromacs/mdlib/nbnxn_internal.h"
57 #include "gromacs/mdlib/nbnxn_search.h"
58 #include "gromacs/mdlib/nbnxn_util.h"
59 #include "gromacs/mdtypes/forcerec.h" // only for GET_CGINFO_*
60 #include "gromacs/mdtypes/mdatom.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/simd/simd.h"
63 #include "gromacs/timing/wallcycle.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/gmxomp.h"
67 #include "gromacs/utility/logger.h"
68 #include "gromacs/utility/smalloc.h"
69 #include "gromacs/utility/stringutil.h"
71 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
73 /* Default nbnxn allocation routine, allocates NBNXN_MEM_ALIGN byte aligned */
74 void nbnxn_alloc_aligned(void **ptr, size_t nbytes)
76 *ptr = save_malloc_aligned("ptr", __FILE__, __LINE__, nbytes, 1, NBNXN_MEM_ALIGN);
79 /* Free function for memory allocated with nbnxn_alloc_aligned */
80 void nbnxn_free_aligned(void *ptr)
82 sfree_aligned(ptr);
85 /* Reallocation wrapper function for nbnxn data structures */
86 void nbnxn_realloc_void(void **ptr,
87 int nbytes_copy, int nbytes_new,
88 nbnxn_alloc_t *ma,
89 nbnxn_free_t *mf)
91 void *ptr_new;
93 ma(&ptr_new, nbytes_new);
95 if (nbytes_new > 0 && ptr_new == nullptr)
97 gmx_fatal(FARGS, "Allocation of %d bytes failed", nbytes_new);
100 if (nbytes_copy > 0)
102 if (nbytes_new < nbytes_copy)
104 gmx_incons("In nbnxn_realloc_void: new size less than copy size");
106 memcpy(ptr_new, *ptr, nbytes_copy);
108 if (*ptr != nullptr)
110 mf(*ptr);
112 *ptr = ptr_new;
115 /* Reallocate the nbnxn_atomdata_t for a size of n atoms */
116 void nbnxn_atomdata_realloc(nbnxn_atomdata_t *nbat, int n)
118 GMX_ASSERT(nbat->nalloc >= nbat->natoms, "We should have at least as many elelements allocated as there are set");
120 int t;
122 nbnxn_realloc_void((void **)&nbat->type,
123 nbat->natoms*sizeof(*nbat->type),
124 n*sizeof(*nbat->type),
125 nbat->alloc, nbat->free);
126 nbnxn_realloc_void((void **)&nbat->lj_comb,
127 nbat->natoms*2*sizeof(*nbat->lj_comb),
128 n*2*sizeof(*nbat->lj_comb),
129 nbat->alloc, nbat->free);
130 if (nbat->XFormat != nbatXYZQ)
132 nbnxn_realloc_void((void **)&nbat->q,
133 nbat->natoms*sizeof(*nbat->q),
134 n*sizeof(*nbat->q),
135 nbat->alloc, nbat->free);
137 if (nbat->nenergrp > 1)
139 nbnxn_realloc_void((void **)&nbat->energrp,
140 nbat->natoms/nbat->na_c*sizeof(*nbat->energrp),
141 n/nbat->na_c*sizeof(*nbat->energrp),
142 nbat->alloc, nbat->free);
144 nbnxn_realloc_void((void **)&nbat->x,
145 nbat->natoms*nbat->xstride*sizeof(*nbat->x),
146 n*nbat->xstride*sizeof(*nbat->x),
147 nbat->alloc, nbat->free);
148 for (t = 0; t < nbat->nout; t++)
150 /* Allocate one element extra for possible signaling with GPUs */
151 nbnxn_realloc_void((void **)&nbat->out[t].f,
152 nbat->natoms*nbat->fstride*sizeof(*nbat->out[t].f),
153 n*nbat->fstride*sizeof(*nbat->out[t].f),
154 nbat->alloc, nbat->free);
156 nbat->nalloc = n;
159 /* Initializes an nbnxn_atomdata_output_t data structure */
160 static void nbnxn_atomdata_output_init(nbnxn_atomdata_output_t *out,
161 int nb_kernel_type,
162 int nenergrp, int stride,
163 nbnxn_alloc_t *ma)
165 out->f = nullptr;
166 ma((void **)&out->fshift, SHIFTS*DIM*sizeof(*out->fshift));
167 out->nV = nenergrp*nenergrp;
168 ma((void **)&out->Vvdw, out->nV*sizeof(*out->Vvdw));
169 ma((void **)&out->Vc, out->nV*sizeof(*out->Vc ));
171 if (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
172 nb_kernel_type == nbnxnk4xN_SIMD_2xNN)
174 int cj_size = nbnxn_kernel_to_cluster_j_size(nb_kernel_type);
175 out->nVS = nenergrp*nenergrp*stride*(cj_size>>1)*cj_size;
176 ma((void **)&out->VSvdw, out->nVS*sizeof(*out->VSvdw));
177 ma((void **)&out->VSc, out->nVS*sizeof(*out->VSc ));
179 else
181 out->nVS = 0;
185 static void copy_int_to_nbat_int(const int *a, int na, int na_round,
186 const int *in, int fill, int *innb)
188 int i, j;
190 j = 0;
191 for (i = 0; i < na; i++)
193 innb[j++] = in[a[i]];
195 /* Complete the partially filled last cell with fill */
196 for (; i < na_round; i++)
198 innb[j++] = fill;
202 void copy_rvec_to_nbat_real(const int *a, int na, int na_round,
203 const rvec *x, int nbatFormat,
204 real *xnb, int a0)
206 /* We complete partially filled cells, can only be the last one in each
207 * column, with coordinates farAway. The actual coordinate value does
208 * not influence the results, since these filler particles do not interact.
209 * Clusters with normal atoms + fillers have a bounding box based only
210 * on the coordinates of the atoms. Clusters with only fillers have as
211 * the bounding box the coordinates of the first filler. Such clusters
212 * are not considered as i-entries, but they are considered as j-entries.
213 * So for performance it is better to have their bounding boxes far away,
214 * such that filler only clusters don't end up in the pair list.
216 const real farAway = -1000000;
218 int i, j, c;
220 switch (nbatFormat)
222 case nbatXYZ:
223 j = a0*STRIDE_XYZ;
224 for (i = 0; i < na; i++)
226 xnb[j++] = x[a[i]][XX];
227 xnb[j++] = x[a[i]][YY];
228 xnb[j++] = x[a[i]][ZZ];
230 /* Complete the partially filled last cell with farAway elements */
231 for (; i < na_round; i++)
233 xnb[j++] = farAway;
234 xnb[j++] = farAway;
235 xnb[j++] = farAway;
237 break;
238 case nbatXYZQ:
239 j = a0*STRIDE_XYZQ;
240 for (i = 0; i < na; i++)
242 xnb[j++] = x[a[i]][XX];
243 xnb[j++] = x[a[i]][YY];
244 xnb[j++] = x[a[i]][ZZ];
245 j++;
247 /* Complete the partially filled last cell with zeros */
248 for (; i < na_round; i++)
250 xnb[j++] = farAway;
251 xnb[j++] = farAway;
252 xnb[j++] = farAway;
253 j++;
255 break;
256 case nbatX4:
257 j = atom_to_x_index<c_packX4>(a0);
258 c = a0 & (c_packX4-1);
259 for (i = 0; i < na; i++)
261 xnb[j+XX*c_packX4] = x[a[i]][XX];
262 xnb[j+YY*c_packX4] = x[a[i]][YY];
263 xnb[j+ZZ*c_packX4] = x[a[i]][ZZ];
264 j++;
265 c++;
266 if (c == c_packX4)
268 j += (DIM-1)*c_packX4;
269 c = 0;
272 /* Complete the partially filled last cell with zeros */
273 for (; i < na_round; i++)
275 xnb[j+XX*c_packX4] = farAway;
276 xnb[j+YY*c_packX4] = farAway;
277 xnb[j+ZZ*c_packX4] = farAway;
278 j++;
279 c++;
280 if (c == c_packX4)
282 j += (DIM-1)*c_packX4;
283 c = 0;
286 break;
287 case nbatX8:
288 j = atom_to_x_index<c_packX8>(a0);
289 c = a0 & (c_packX8 - 1);
290 for (i = 0; i < na; i++)
292 xnb[j+XX*c_packX8] = x[a[i]][XX];
293 xnb[j+YY*c_packX8] = x[a[i]][YY];
294 xnb[j+ZZ*c_packX8] = x[a[i]][ZZ];
295 j++;
296 c++;
297 if (c == c_packX8)
299 j += (DIM-1)*c_packX8;
300 c = 0;
303 /* Complete the partially filled last cell with zeros */
304 for (; i < na_round; i++)
306 xnb[j+XX*c_packX8] = farAway;
307 xnb[j+YY*c_packX8] = farAway;
308 xnb[j+ZZ*c_packX8] = farAway;
309 j++;
310 c++;
311 if (c == c_packX8)
313 j += (DIM-1)*c_packX8;
314 c = 0;
317 break;
318 default:
319 gmx_incons("Unsupported nbnxn_atomdata_t format");
323 /* Stores the LJ parameter data in a format convenient for different kernels */
324 static void set_lj_parameter_data(nbnxn_atomdata_t *nbat, gmx_bool bSIMD)
326 real c6, c12;
328 int nt = nbat->ntype;
330 if (bSIMD)
332 #if GMX_SIMD
333 /* nbfp_aligned stores two parameters using the stride most suitable
334 * for the present SIMD architecture, as specified by the constant
335 * c_simdBestPairAlignment from the SIMD header.
336 * There's a slight inefficiency in allocating and initializing nbfp_aligned
337 * when it might not be used, but introducing the conditional code is not
338 * really worth it.
340 nbat->alloc((void **)&nbat->nbfp_aligned,
341 nt*nt*c_simdBestPairAlignment*sizeof(*nbat->nbfp_aligned));
342 for (int i = 0; i < nt; i++)
344 for (int j = 0; j < nt; j++)
346 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+0] = nbat->nbfp[(i*nt+j)*2+0];
347 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+1] = nbat->nbfp[(i*nt+j)*2+1];
348 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
349 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
352 #endif
355 /* We use combination rule data for SIMD combination rule kernels
356 * and with LJ-PME kernels. We then only need parameters per atom type,
357 * not per pair of atom types.
359 switch (nbat->comb_rule)
361 case ljcrGEOM:
362 nbat->comb_rule = ljcrGEOM;
364 for (int i = 0; i < nt; i++)
366 /* Store the sqrt of the diagonal from the nbfp matrix */
367 nbat->nbfp_comb[i*2 ] = std::sqrt(nbat->nbfp[(i*nt+i)*2 ]);
368 nbat->nbfp_comb[i*2+1] = std::sqrt(nbat->nbfp[(i*nt+i)*2+1]);
370 break;
371 case ljcrLB:
372 for (int i = 0; i < nt; i++)
374 /* Get 6*C6 and 12*C12 from the diagonal of the nbfp matrix */
375 c6 = nbat->nbfp[(i*nt+i)*2 ];
376 c12 = nbat->nbfp[(i*nt+i)*2+1];
377 if (c6 > 0 && c12 > 0)
379 /* We store 0.5*2^1/6*sigma and sqrt(4*3*eps),
380 * so we get 6*C6 and 12*C12 after combining.
382 nbat->nbfp_comb[i*2 ] = 0.5*gmx::sixthroot(c12/c6);
383 nbat->nbfp_comb[i*2+1] = std::sqrt(c6*c6/c12);
385 else
387 nbat->nbfp_comb[i*2 ] = 0;
388 nbat->nbfp_comb[i*2+1] = 0;
391 break;
392 case ljcrNONE:
393 /* We always store the full matrix (see code above) */
394 break;
395 default:
396 gmx_incons("Unknown combination rule");
397 break;
401 #if GMX_SIMD
402 static void
403 nbnxn_atomdata_init_simple_exclusion_masks(nbnxn_atomdata_t *nbat)
405 const int simd_width = GMX_SIMD_REAL_WIDTH;
406 int simd_excl_size;
407 /* Set the diagonal cluster pair exclusion mask setup data.
408 * In the kernel we check 0 < j - i to generate the masks.
409 * Here we store j - i for generating the mask for the first i,
410 * we substract 0.5 to avoid rounding issues.
411 * In the kernel we can subtract 1 to generate the subsequent mask.
413 int simd_4xn_diag_size;
415 simd_4xn_diag_size = std::max(NBNXN_CPU_CLUSTER_I_SIZE, simd_width);
416 snew_aligned(nbat->simd_4xn_diagonal_j_minus_i, simd_4xn_diag_size, NBNXN_MEM_ALIGN);
417 for (int j = 0; j < simd_4xn_diag_size; j++)
419 nbat->simd_4xn_diagonal_j_minus_i[j] = j - 0.5;
422 snew_aligned(nbat->simd_2xnn_diagonal_j_minus_i, simd_width, NBNXN_MEM_ALIGN);
423 for (int j = 0; j < simd_width/2; j++)
425 /* The j-cluster size is half the SIMD width */
426 nbat->simd_2xnn_diagonal_j_minus_i[j] = j - 0.5;
427 /* The next half of the SIMD width is for i + 1 */
428 nbat->simd_2xnn_diagonal_j_minus_i[simd_width/2+j] = j - 1 - 0.5;
431 /* We use up to 32 bits for exclusion masking.
432 * The same masks are used for the 4xN and 2x(N+N) kernels.
433 * The masks are read either into integer SIMD registers or into
434 * real SIMD registers (together with a cast).
435 * In single precision this means the real and integer SIMD registers
436 * are of equal size.
438 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*simd_width;
439 #if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
440 snew_aligned(nbat->simd_exclusion_filter64, simd_excl_size, NBNXN_MEM_ALIGN);
441 #else
442 snew_aligned(nbat->simd_exclusion_filter, simd_excl_size, NBNXN_MEM_ALIGN);
443 #endif
445 for (int j = 0; j < simd_excl_size; j++)
447 /* Set the consecutive bits for masking pair exclusions */
448 #if GMX_DOUBLE && !GMX_SIMD_HAVE_INT32_LOGICAL
449 nbat->simd_exclusion_filter64[j] = (1U << j);
450 #else
451 nbat->simd_exclusion_filter[j] = (1U << j);
452 #endif
455 #if !GMX_SIMD_HAVE_LOGICAL && !GMX_SIMD_HAVE_INT32_LOGICAL
456 // If the SIMD implementation has no bitwise logical operation support
457 // whatsoever we cannot use the normal masking. Instead,
458 // we generate a vector of all 2^4 possible ways an i atom
459 // interacts with its 4 j atoms. Each array entry contains
460 // GMX_SIMD_REAL_WIDTH values that are read with a single aligned SIMD load.
461 // Since there is no logical value representation in this case, we use
462 // any nonzero value to indicate 'true', while zero mean 'false'.
463 // This can then be converted to a SIMD boolean internally in the SIMD
464 // module by comparing to zero.
465 // Each array entry encodes how this i atom will interact with the 4 j atoms.
466 // Matching code exists in set_ci_top_excls() to generate indices into this array.
467 // Those indices are used in the kernels.
469 simd_excl_size = NBNXN_CPU_CLUSTER_I_SIZE*NBNXN_CPU_CLUSTER_I_SIZE;
470 const real simdFalse = 0.0;
471 const real simdTrue = 1.0;
472 real *simd_interaction_array;
474 snew_aligned(simd_interaction_array, simd_excl_size * GMX_SIMD_REAL_WIDTH, NBNXN_MEM_ALIGN);
475 for (int j = 0; j < simd_excl_size; j++)
477 int index = j * GMX_SIMD_REAL_WIDTH;
478 for (int i = 0; i < GMX_SIMD_REAL_WIDTH; i++)
480 simd_interaction_array[index + i] = (j & (1 << i)) ? simdTrue : simdFalse;
483 nbat->simd_interaction_array = simd_interaction_array;
484 #endif
486 #endif
488 /* Initializes an nbnxn_atomdata_t data structure */
489 void nbnxn_atomdata_init(const gmx::MDLogger &mdlog,
490 nbnxn_atomdata_t *nbat,
491 int nb_kernel_type,
492 int enbnxninitcombrule,
493 int ntype, const real *nbfp,
494 int n_energygroups,
495 int nout,
496 nbnxn_alloc_t *alloc,
497 nbnxn_free_t *free)
499 int nth;
500 real c6, c12, tol;
501 char *ptr;
502 gmx_bool simple, bCombGeom, bCombLB, bSIMD;
504 if (alloc == nullptr)
506 nbat->alloc = nbnxn_alloc_aligned;
508 else
510 nbat->alloc = alloc;
512 if (free == nullptr)
514 nbat->free = nbnxn_free_aligned;
516 else
518 nbat->free = free;
521 if (debug)
523 fprintf(debug, "There are %d atom types in the system, adding one for nbnxn_atomdata_t\n", ntype);
525 nbat->ntype = ntype + 1;
526 nbat->alloc((void **)&nbat->nbfp,
527 nbat->ntype*nbat->ntype*2*sizeof(*nbat->nbfp));
528 nbat->alloc((void **)&nbat->nbfp_comb, nbat->ntype*2*sizeof(*nbat->nbfp_comb));
530 /* A tolerance of 1e-5 seems reasonable for (possibly hand-typed)
531 * force-field floating point parameters.
533 tol = 1e-5;
534 ptr = getenv("GMX_LJCOMB_TOL");
535 if (ptr != nullptr)
537 double dbl;
539 sscanf(ptr, "%lf", &dbl);
540 tol = dbl;
542 bCombGeom = TRUE;
543 bCombLB = TRUE;
545 /* Temporarily fill nbat->nbfp_comb with sigma and epsilon
546 * to check for the LB rule.
548 for (int i = 0; i < ntype; i++)
550 c6 = nbfp[(i*ntype+i)*2 ]/6.0;
551 c12 = nbfp[(i*ntype+i)*2+1]/12.0;
552 if (c6 > 0 && c12 > 0)
554 nbat->nbfp_comb[i*2 ] = gmx::sixthroot(c12/c6);
555 nbat->nbfp_comb[i*2+1] = 0.25*c6*c6/c12;
557 else if (c6 == 0 && c12 == 0)
559 nbat->nbfp_comb[i*2 ] = 0;
560 nbat->nbfp_comb[i*2+1] = 0;
562 else
564 /* Can not use LB rule with only dispersion or repulsion */
565 bCombLB = FALSE;
569 for (int i = 0; i < nbat->ntype; i++)
571 for (int j = 0; j < nbat->ntype; j++)
573 if (i < ntype && j < ntype)
575 /* fr->nbfp has been updated, so that array too now stores c6/c12 including
576 * the 6.0/12.0 prefactors to save 2 flops in the most common case (force-only).
578 c6 = nbfp[(i*ntype+j)*2 ];
579 c12 = nbfp[(i*ntype+j)*2+1];
580 nbat->nbfp[(i*nbat->ntype+j)*2 ] = c6;
581 nbat->nbfp[(i*nbat->ntype+j)*2+1] = c12;
583 /* Compare 6*C6 and 12*C12 for geometric cobination rule */
584 bCombGeom = bCombGeom &&
585 gmx_within_tol(c6*c6, nbfp[(i*ntype+i)*2 ]*nbfp[(j*ntype+j)*2 ], tol) &&
586 gmx_within_tol(c12*c12, nbfp[(i*ntype+i)*2+1]*nbfp[(j*ntype+j)*2+1], tol);
588 /* Compare C6 and C12 for Lorentz-Berthelot combination rule */
589 c6 /= 6.0;
590 c12 /= 12.0;
591 bCombLB = bCombLB &&
592 ((c6 == 0 && c12 == 0 &&
593 (nbat->nbfp_comb[i*2+1] == 0 || nbat->nbfp_comb[j*2+1] == 0)) ||
594 (c6 > 0 && c12 > 0 &&
595 gmx_within_tol(gmx::sixthroot(c12/c6),
596 0.5*(nbat->nbfp_comb[i*2]+nbat->nbfp_comb[j*2]), tol) &&
597 gmx_within_tol(0.25*c6*c6/c12, std::sqrt(nbat->nbfp_comb[i*2+1]*nbat->nbfp_comb[j*2+1]), tol)));
599 else
601 /* Add zero parameters for the additional dummy atom type */
602 nbat->nbfp[(i*nbat->ntype+j)*2 ] = 0;
603 nbat->nbfp[(i*nbat->ntype+j)*2+1] = 0;
607 if (debug)
609 fprintf(debug, "Combination rules: geometric %d Lorentz-Berthelot %d\n",
610 bCombGeom, bCombLB);
613 simple = nbnxn_kernel_pairlist_simple(nb_kernel_type);
615 switch (enbnxninitcombrule)
617 case enbnxninitcombruleDETECT:
618 /* We prefer the geometic combination rule,
619 * as that gives a slightly faster kernel than the LB rule.
621 if (bCombGeom)
623 nbat->comb_rule = ljcrGEOM;
625 else if (bCombLB)
627 nbat->comb_rule = ljcrLB;
629 else
631 nbat->comb_rule = ljcrNONE;
633 nbat->free(nbat->nbfp_comb);
637 std::string mesg;
638 if (nbat->comb_rule == ljcrNONE)
640 mesg = "Using full Lennard-Jones parameter combination matrix";
642 else
644 mesg = gmx::formatString("Using %s Lennard-Jones combination rule",
645 nbat->comb_rule == ljcrGEOM ? "geometric" : "Lorentz-Berthelot");
647 GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
649 break;
650 case enbnxninitcombruleGEOM:
651 nbat->comb_rule = ljcrGEOM;
652 break;
653 case enbnxninitcombruleLB:
654 nbat->comb_rule = ljcrLB;
655 break;
656 case enbnxninitcombruleNONE:
657 nbat->comb_rule = ljcrNONE;
659 nbat->free(nbat->nbfp_comb);
660 break;
661 default:
662 gmx_incons("Unknown enbnxninitcombrule");
665 bSIMD = (nb_kernel_type == nbnxnk4xN_SIMD_4xN ||
666 nb_kernel_type == nbnxnk4xN_SIMD_2xNN);
668 set_lj_parameter_data(nbat, bSIMD);
670 nbat->natoms = 0;
671 nbat->type = nullptr;
672 nbat->lj_comb = nullptr;
673 if (simple)
675 int pack_x;
677 if (bSIMD)
679 pack_x = std::max(NBNXN_CPU_CLUSTER_I_SIZE,
680 nbnxn_kernel_to_cluster_j_size(nb_kernel_type));
681 switch (pack_x)
683 case 4:
684 nbat->XFormat = nbatX4;
685 break;
686 case 8:
687 nbat->XFormat = nbatX8;
688 break;
689 default:
690 gmx_incons("Unsupported packing width");
693 else
695 nbat->XFormat = nbatXYZ;
698 nbat->FFormat = nbat->XFormat;
700 else
702 nbat->XFormat = nbatXYZQ;
703 nbat->FFormat = nbatXYZ;
705 nbat->q = nullptr;
706 nbat->nenergrp = n_energygroups;
707 if (!simple)
709 // We now check for energy groups already when starting mdrun
710 GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
712 /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
713 if (nbat->nenergrp > 64)
715 gmx_fatal(FARGS, "With NxN kernels not more than 64 energy groups are supported\n");
717 nbat->neg_2log = 1;
718 while (nbat->nenergrp > (1<<nbat->neg_2log))
720 nbat->neg_2log++;
722 nbat->energrp = nullptr;
723 nbat->alloc((void **)&nbat->shift_vec, SHIFTS*sizeof(*nbat->shift_vec));
724 nbat->xstride = (nbat->XFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
725 nbat->fstride = (nbat->FFormat == nbatXYZQ ? STRIDE_XYZQ : DIM);
726 nbat->x = nullptr;
728 #if GMX_SIMD
729 if (simple)
731 nbnxn_atomdata_init_simple_exclusion_masks(nbat);
733 #endif
735 /* Initialize the output data structures */
736 nbat->nout = nout;
737 snew(nbat->out, nbat->nout);
738 nbat->nalloc = 0;
739 for (int i = 0; i < nbat->nout; i++)
741 nbnxn_atomdata_output_init(&nbat->out[i],
742 nb_kernel_type,
743 nbat->nenergrp, 1<<nbat->neg_2log,
744 nbat->alloc);
746 nbat->buffer_flags.flag = nullptr;
747 nbat->buffer_flags.flag_nalloc = 0;
749 nth = gmx_omp_nthreads_get(emntNonbonded);
751 ptr = getenv("GMX_USE_TREEREDUCE");
752 if (ptr != nullptr)
754 nbat->bUseTreeReduce = strtol(ptr, nullptr, 10);
756 #if defined __MIC__
757 else if (nth > 8) /*on the CPU we currently don't benefit even at 32*/
759 nbat->bUseTreeReduce = 1;
761 #endif
762 else
764 nbat->bUseTreeReduce = 0;
766 if (nbat->bUseTreeReduce)
768 GMX_LOG(mdlog.info).asParagraph().appendText("Using tree force reduction");
770 snew(nbat->syncStep, nth);
774 template<int packSize>
775 static void copy_lj_to_nbat_lj_comb(const real *ljparam_type,
776 const int *type, int na,
777 real *ljparam_at)
779 /* The LJ params follow the combination rule:
780 * copy the params for the type array to the atom array.
782 for (int is = 0; is < na; is += packSize)
784 for (int k = 0; k < packSize; k++)
786 int i = is + k;
787 ljparam_at[is*2 + k] = ljparam_type[type[i]*2 ];
788 ljparam_at[is*2 + packSize + k] = ljparam_type[type[i]*2 + 1];
793 /* Sets the atom type in nbnxn_atomdata_t */
794 static void nbnxn_atomdata_set_atomtypes(nbnxn_atomdata_t *nbat,
795 const nbnxn_search_t nbs,
796 const int *type)
798 for (const nbnxn_grid_t &grid : nbs->grid)
800 /* Loop over all columns and copy and fill */
801 for (int i = 0; i < grid.ncx*grid.ncy; i++)
803 int ncz = grid.cxy_ind[i+1] - grid.cxy_ind[i];
804 int ash = (grid.cell0 + grid.cxy_ind[i])*grid.na_sc;
806 copy_int_to_nbat_int(nbs->a.data() + ash, grid.cxy_na[i], ncz*grid.na_sc,
807 type, nbat->ntype-1, nbat->type+ash);
812 /* Sets the LJ combination rule parameters in nbnxn_atomdata_t */
813 static void nbnxn_atomdata_set_ljcombparams(nbnxn_atomdata_t *nbat,
814 const nbnxn_search_t nbs)
816 if (nbat->comb_rule != ljcrNONE)
818 for (const nbnxn_grid_t &grid : nbs->grid)
820 /* Loop over all columns and copy and fill */
821 for (int i = 0; i < grid.ncx*grid.ncy; i++)
823 int ncz = grid.cxy_ind[i+1] - grid.cxy_ind[i];
824 int ash = (grid.cell0 + grid.cxy_ind[i])*grid.na_sc;
826 if (nbat->XFormat == nbatX4)
828 copy_lj_to_nbat_lj_comb<c_packX4>(nbat->nbfp_comb,
829 nbat->type + ash,
830 ncz*grid.na_sc,
831 nbat->lj_comb + ash*2);
833 else if (nbat->XFormat == nbatX8)
835 copy_lj_to_nbat_lj_comb<c_packX8>(nbat->nbfp_comb,
836 nbat->type + ash,
837 ncz*grid.na_sc,
838 nbat->lj_comb + ash*2);
840 else if (nbat->XFormat == nbatXYZQ)
842 copy_lj_to_nbat_lj_comb<1>(nbat->nbfp_comb,
843 nbat->type + ash,
844 ncz*grid.na_sc,
845 nbat->lj_comb + ash*2);
852 /* Sets the charges in nbnxn_atomdata_t *nbat */
853 static void nbnxn_atomdata_set_charges(nbnxn_atomdata_t *nbat,
854 const nbnxn_search_t nbs,
855 const real *charge)
857 for (const nbnxn_grid_t &grid : nbs->grid)
859 /* Loop over all columns and copy and fill */
860 for (int cxy = 0; cxy < grid.ncx*grid.ncy; cxy++)
862 int ash = (grid.cell0 + grid.cxy_ind[cxy])*grid.na_sc;
863 int na = grid.cxy_na[cxy];
864 int na_round = (grid.cxy_ind[cxy+1] - grid.cxy_ind[cxy])*grid.na_sc;
866 if (nbat->XFormat == nbatXYZQ)
868 real *q = nbat->x + ash*STRIDE_XYZQ + ZZ + 1;
869 int i;
870 for (i = 0; i < na; i++)
872 *q = charge[nbs->a[ash+i]];
873 q += STRIDE_XYZQ;
875 /* Complete the partially filled last cell with zeros */
876 for (; i < na_round; i++)
878 *q = 0;
879 q += STRIDE_XYZQ;
882 else
884 real *q = nbat->q + ash;
885 int i;
886 for (i = 0; i < na; i++)
888 *q = charge[nbs->a[ash+i]];
889 q++;
891 /* Complete the partially filled last cell with zeros */
892 for (; i < na_round; i++)
894 *q = 0;
895 q++;
902 /* Set the charges of perturbed atoms in nbnxn_atomdata_t to 0.
903 * This is to automatically remove the RF/PME self term in the nbnxn kernels.
904 * Part of the zero interactions are still calculated in the normal kernels.
905 * All perturbed interactions are calculated in the free energy kernel,
906 * using the original charge and LJ data, not nbnxn_atomdata_t.
908 static void nbnxn_atomdata_mask_fep(nbnxn_atomdata_t *nbat,
909 const nbnxn_search_t nbs)
911 real *q;
912 int stride_q, nsubc;
914 if (nbat->XFormat == nbatXYZQ)
916 q = nbat->x + ZZ + 1;
917 stride_q = STRIDE_XYZQ;
919 else
921 q = nbat->q;
922 stride_q = 1;
925 for (const nbnxn_grid_t &grid : nbs->grid)
927 if (grid.bSimple)
929 nsubc = 1;
931 else
933 nsubc = c_gpuNumClusterPerCell;
936 int c_offset = grid.cell0*grid.na_sc;
938 /* Loop over all columns and copy and fill */
939 for (int c = 0; c < grid.nc*nsubc; c++)
941 /* Does this cluster contain perturbed particles? */
942 if (grid.fep[c] != 0)
944 for (int i = 0; i < grid.na_c; i++)
946 /* Is this a perturbed particle? */
947 if (grid.fep[c] & (1 << i))
949 int ind = c_offset + c*grid.na_c + i;
950 /* Set atom type and charge to non-interacting */
951 nbat->type[ind] = nbat->ntype - 1;
952 q[ind*stride_q] = 0;
960 /* Copies the energy group indices to a reordered and packed array */
961 static void copy_egp_to_nbat_egps(const int *a, int na, int na_round,
962 int na_c, int bit_shift,
963 const int *in, int *innb)
965 int i;
966 int comb;
968 int j = 0;
969 for (i = 0; i < na; i += na_c)
971 /* Store na_c energy group numbers into one int */
972 comb = 0;
973 for (int sa = 0; sa < na_c; sa++)
975 int at = a[i+sa];
976 if (at >= 0)
978 comb |= (GET_CGINFO_GID(in[at]) << (sa*bit_shift));
981 innb[j++] = comb;
983 /* Complete the partially filled last cell with fill */
984 for (; i < na_round; i += na_c)
986 innb[j++] = 0;
990 /* Set the energy group indices for atoms in nbnxn_atomdata_t */
991 static void nbnxn_atomdata_set_energygroups(nbnxn_atomdata_t *nbat,
992 const nbnxn_search_t nbs,
993 const int *atinfo)
995 if (nbat->nenergrp == 1)
997 return;
1000 for (const nbnxn_grid_t &grid : nbs->grid)
1002 /* Loop over all columns and copy and fill */
1003 for (int i = 0; i < grid.ncx*grid.ncy; i++)
1005 int ncz = grid.cxy_ind[i+1] - grid.cxy_ind[i];
1006 int ash = (grid.cell0 + grid.cxy_ind[i])*grid.na_sc;
1008 copy_egp_to_nbat_egps(nbs->a.data() + ash, grid.cxy_na[i], ncz*grid.na_sc,
1009 nbat->na_c, nbat->neg_2log,
1010 atinfo, nbat->energrp+(ash>>grid.na_c_2log));
1015 /* Sets all required atom parameter data in nbnxn_atomdata_t */
1016 void nbnxn_atomdata_set(nbnxn_atomdata_t *nbat,
1017 const nbnxn_search_t nbs,
1018 const t_mdatoms *mdatoms,
1019 const int *atinfo)
1021 nbnxn_atomdata_set_atomtypes(nbat, nbs, mdatoms->typeA);
1023 nbnxn_atomdata_set_charges(nbat, nbs, mdatoms->chargeA);
1025 if (nbs->bFEP)
1027 nbnxn_atomdata_mask_fep(nbat, nbs);
1030 /* This must be done after masking types for FEP */
1031 nbnxn_atomdata_set_ljcombparams(nbat, nbs);
1033 nbnxn_atomdata_set_energygroups(nbat, nbs, atinfo);
1036 /* Copies the shift vector array to nbnxn_atomdata_t */
1037 void nbnxn_atomdata_copy_shiftvec(gmx_bool bDynamicBox,
1038 rvec *shift_vec,
1039 nbnxn_atomdata_t *nbat)
1041 int i;
1043 nbat->bDynamicBox = bDynamicBox;
1044 for (i = 0; i < SHIFTS; i++)
1046 copy_rvec(shift_vec[i], nbat->shift_vec[i]);
1050 /* Copies (and reorders) the coordinates to nbnxn_atomdata_t */
1051 void nbnxn_atomdata_copy_x_to_nbat_x(const nbnxn_search_t nbs,
1052 int locality,
1053 gmx_bool FillLocal,
1054 rvec *x,
1055 nbnxn_atomdata_t *nbat,
1056 gmx_wallcycle *wcycle)
1058 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1059 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1061 int g0 = 0, g1 = 0;
1062 int nth, th;
1064 switch (locality)
1066 case eatAll:
1067 g0 = 0;
1068 g1 = nbs->grid.size();
1069 break;
1070 case eatLocal:
1071 g0 = 0;
1072 g1 = 1;
1073 break;
1074 case eatNonlocal:
1075 g0 = 1;
1076 g1 = nbs->grid.size();
1077 break;
1080 if (FillLocal)
1082 nbat->natoms_local = nbs->grid[0].nc*nbs->grid[0].na_sc;
1085 nth = gmx_omp_nthreads_get(emntPairsearch);
1087 #pragma omp parallel for num_threads(nth) schedule(static)
1088 for (th = 0; th < nth; th++)
1092 for (int g = g0; g < g1; g++)
1094 const nbnxn_grid_t *grid;
1095 int cxy0, cxy1;
1097 grid = &nbs->grid[g];
1099 cxy0 = (grid->ncx*grid->ncy* th +nth-1)/nth;
1100 cxy1 = (grid->ncx*grid->ncy*(th+1)+nth-1)/nth;
1102 for (int cxy = cxy0; cxy < cxy1; cxy++)
1104 int na, ash, na_fill;
1106 na = grid->cxy_na[cxy];
1107 ash = (grid->cell0 + grid->cxy_ind[cxy])*grid->na_sc;
1109 if (g == 0 && FillLocal)
1111 na_fill =
1112 (grid->cxy_ind[cxy+1] - grid->cxy_ind[cxy])*grid->na_sc;
1114 else
1116 /* We fill only the real particle locations.
1117 * We assume the filling entries at the end have been
1118 * properly set before during pair-list generation.
1120 na_fill = na;
1122 copy_rvec_to_nbat_real(nbs->a.data() + ash, na, na_fill, x,
1123 nbat->XFormat, nbat->x, ash);
1127 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1130 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1131 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1134 static void
1135 nbnxn_atomdata_clear_reals(real * gmx_restrict dest,
1136 int i0, int i1)
1138 for (int i = i0; i < i1; i++)
1140 dest[i] = 0;
1144 gmx_unused static void
1145 nbnxn_atomdata_reduce_reals(real * gmx_restrict dest,
1146 gmx_bool bDestSet,
1147 real ** gmx_restrict src,
1148 int nsrc,
1149 int i0, int i1)
1151 if (bDestSet)
1153 /* The destination buffer contains data, add to it */
1154 for (int i = i0; i < i1; i++)
1156 for (int s = 0; s < nsrc; s++)
1158 dest[i] += src[s][i];
1162 else
1164 /* The destination buffer is unitialized, set it first */
1165 for (int i = i0; i < i1; i++)
1167 dest[i] = src[0][i];
1168 for (int s = 1; s < nsrc; s++)
1170 dest[i] += src[s][i];
1176 gmx_unused static void
1177 nbnxn_atomdata_reduce_reals_simd(real gmx_unused * gmx_restrict dest,
1178 gmx_bool gmx_unused bDestSet,
1179 real gmx_unused ** gmx_restrict src,
1180 int gmx_unused nsrc,
1181 int gmx_unused i0, int gmx_unused i1)
1183 #if GMX_SIMD
1184 /* The SIMD width here is actually independent of that in the kernels,
1185 * but we use the same width for simplicity (usually optimal anyhow).
1187 SimdReal dest_SSE, src_SSE;
1189 if (bDestSet)
1191 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1193 dest_SSE = load<SimdReal>(dest+i);
1194 for (int s = 0; s < nsrc; s++)
1196 src_SSE = load<SimdReal>(src[s]+i);
1197 dest_SSE = dest_SSE + src_SSE;
1199 store(dest+i, dest_SSE);
1202 else
1204 for (int i = i0; i < i1; i += GMX_SIMD_REAL_WIDTH)
1206 dest_SSE = load<SimdReal>(src[0]+i);
1207 for (int s = 1; s < nsrc; s++)
1209 src_SSE = load<SimdReal>(src[s]+i);
1210 dest_SSE = dest_SSE + src_SSE;
1212 store(dest+i, dest_SSE);
1215 #endif
1218 /* Add part of the force array(s) from nbnxn_atomdata_t to f */
1219 static void
1220 nbnxn_atomdata_add_nbat_f_to_f_part(const nbnxn_search_t nbs,
1221 const nbnxn_atomdata_t *nbat,
1222 nbnxn_atomdata_output_t *out,
1223 int nfa,
1224 int a0, int a1,
1225 rvec *f)
1227 gmx::ArrayRef<const int> cell = nbs->cell;
1229 /* Loop over all columns and copy and fill */
1230 switch (nbat->FFormat)
1232 case nbatXYZ:
1233 case nbatXYZQ:
1234 if (nfa == 1)
1236 const real *fnb = out[0].f;
1238 for (int a = a0; a < a1; a++)
1240 int i = cell[a]*nbat->fstride;
1242 f[a][XX] += fnb[i];
1243 f[a][YY] += fnb[i+1];
1244 f[a][ZZ] += fnb[i+2];
1247 else
1249 for (int a = a0; a < a1; a++)
1251 int i = cell[a]*nbat->fstride;
1253 for (int fa = 0; fa < nfa; fa++)
1255 f[a][XX] += out[fa].f[i];
1256 f[a][YY] += out[fa].f[i+1];
1257 f[a][ZZ] += out[fa].f[i+2];
1261 break;
1262 case nbatX4:
1263 if (nfa == 1)
1265 const real *fnb = out[0].f;
1267 for (int a = a0; a < a1; a++)
1269 int i = atom_to_x_index<c_packX4>(cell[a]);
1271 f[a][XX] += fnb[i+XX*c_packX4];
1272 f[a][YY] += fnb[i+YY*c_packX4];
1273 f[a][ZZ] += fnb[i+ZZ*c_packX4];
1276 else
1278 for (int a = a0; a < a1; a++)
1280 int i = atom_to_x_index<c_packX4>(cell[a]);
1282 for (int fa = 0; fa < nfa; fa++)
1284 f[a][XX] += out[fa].f[i+XX*c_packX4];
1285 f[a][YY] += out[fa].f[i+YY*c_packX4];
1286 f[a][ZZ] += out[fa].f[i+ZZ*c_packX4];
1290 break;
1291 case nbatX8:
1292 if (nfa == 1)
1294 const real *fnb = out[0].f;
1296 for (int a = a0; a < a1; a++)
1298 int i = atom_to_x_index<c_packX8>(cell[a]);
1300 f[a][XX] += fnb[i+XX*c_packX8];
1301 f[a][YY] += fnb[i+YY*c_packX8];
1302 f[a][ZZ] += fnb[i+ZZ*c_packX8];
1305 else
1307 for (int a = a0; a < a1; a++)
1309 int i = atom_to_x_index<c_packX8>(cell[a]);
1311 for (int fa = 0; fa < nfa; fa++)
1313 f[a][XX] += out[fa].f[i+XX*c_packX8];
1314 f[a][YY] += out[fa].f[i+YY*c_packX8];
1315 f[a][ZZ] += out[fa].f[i+ZZ*c_packX8];
1319 break;
1320 default:
1321 gmx_incons("Unsupported nbnxn_atomdata_t format");
1325 static inline unsigned char reverse_bits(unsigned char b)
1327 /* http://graphics.stanford.edu/~seander/bithacks.html#ReverseByteWith64BitsDiv */
1328 return (b * 0x0202020202ULL & 0x010884422010ULL) % 1023;
1331 static void nbnxn_atomdata_add_nbat_f_to_f_treereduce(const nbnxn_atomdata_t *nbat,
1332 int nth)
1334 const nbnxn_buffer_flags_t *flags = &nbat->buffer_flags;
1336 int next_pow2 = 1<<(gmx::log2I(nth-1)+1);
1338 assert(nbat->nout == nth); /* tree-reduce currently only works for nout==nth */
1340 memset(nbat->syncStep, 0, sizeof(*(nbat->syncStep))*nth);
1342 #pragma omp parallel num_threads(nth)
1346 int b0, b1, b;
1347 int i0, i1;
1348 int group_size, th;
1350 th = gmx_omp_get_thread_num();
1352 for (group_size = 2; group_size < 2*next_pow2; group_size *= 2)
1354 int index[2], group_pos, partner_pos, wu;
1355 int partner_th = th ^ (group_size/2);
1357 if (group_size > 2)
1359 #ifdef TMPI_ATOMICS
1360 /* wait on partner thread - replaces full barrier */
1361 int sync_th, sync_group_size;
1363 tMPI_Atomic_memory_barrier(); /* gurantee data is saved before marking work as done */
1364 tMPI_Atomic_set(&(nbat->syncStep[th]), group_size/2); /* mark previous step as completed */
1366 /* find thread to sync with. Equal to partner_th unless nth is not a power of two. */
1367 for (sync_th = partner_th, sync_group_size = group_size; sync_th >= nth && sync_group_size > 2; sync_group_size /= 2)
1369 sync_th &= ~(sync_group_size/4);
1371 if (sync_th < nth) /* otherwise nothing to sync index[1] will be >=nout */
1373 /* wait on the thread which computed input data in previous step */
1374 while (tMPI_Atomic_get((volatile tMPI_Atomic_t*)&(nbat->syncStep[sync_th])) < group_size/2)
1376 gmx_pause();
1378 /* guarantee that no later load happens before wait loop is finisehd */
1379 tMPI_Atomic_memory_barrier();
1381 #else /* TMPI_ATOMICS */
1382 #pragma omp barrier
1383 #endif
1386 /* Calculate buffers to sum (result goes into first buffer) */
1387 group_pos = th % group_size;
1388 index[0] = th - group_pos;
1389 index[1] = index[0] + group_size/2;
1391 /* If no second buffer, nothing to do */
1392 if (index[1] >= nbat->nout && group_size > 2)
1394 continue;
1397 #if NBNXN_BUFFERFLAG_MAX_THREADS > 256
1398 #error reverse_bits assumes max 256 threads
1399 #endif
1400 /* Position is permuted so that one of the 2 vectors being added was computed on the same thread in the previous step.
1401 This improves locality and enables to sync with just a single thread between steps (=the levels in the btree).
1402 The permutation which allows this corresponds to reversing the bits of the group position.
1404 group_pos = reverse_bits(group_pos)/(256/group_size);
1406 partner_pos = group_pos ^ 1;
1408 /* loop over two work-units (own and partner) */
1409 for (wu = 0; wu < 2; wu++)
1411 if (wu == 1)
1413 if (partner_th < nth)
1415 break; /* partner exists we don't have to do his work */
1417 else
1419 group_pos = partner_pos;
1423 /* Calculate the cell-block range for our thread */
1424 b0 = (flags->nflag* group_pos )/group_size;
1425 b1 = (flags->nflag*(group_pos+1))/group_size;
1427 for (b = b0; b < b1; b++)
1429 i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1430 i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1432 if (bitmask_is_set(flags->flag[b], index[1]) || group_size > 2)
1434 #if GMX_SIMD
1435 nbnxn_atomdata_reduce_reals_simd
1436 #else
1437 nbnxn_atomdata_reduce_reals
1438 #endif
1439 (nbat->out[index[0]].f,
1440 bitmask_is_set(flags->flag[b], index[0]) || group_size > 2,
1441 &(nbat->out[index[1]].f), 1, i0, i1);
1444 else if (!bitmask_is_set(flags->flag[b], index[0]))
1446 nbnxn_atomdata_clear_reals(nbat->out[index[0]].f,
1447 i0, i1);
1453 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1458 static void nbnxn_atomdata_add_nbat_f_to_f_stdreduce(const nbnxn_atomdata_t *nbat,
1459 int nth)
1461 #pragma omp parallel for num_threads(nth) schedule(static)
1462 for (int th = 0; th < nth; th++)
1466 const nbnxn_buffer_flags_t *flags;
1467 int nfptr;
1468 real *fptr[NBNXN_BUFFERFLAG_MAX_THREADS];
1470 flags = &nbat->buffer_flags;
1472 /* Calculate the cell-block range for our thread */
1473 int b0 = (flags->nflag* th )/nth;
1474 int b1 = (flags->nflag*(th+1))/nth;
1476 for (int b = b0; b < b1; b++)
1478 int i0 = b *NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1479 int i1 = (b+1)*NBNXN_BUFFERFLAG_SIZE*nbat->fstride;
1481 nfptr = 0;
1482 for (int out = 1; out < nbat->nout; out++)
1484 if (bitmask_is_set(flags->flag[b], out))
1486 fptr[nfptr++] = nbat->out[out].f;
1489 if (nfptr > 0)
1491 #if GMX_SIMD
1492 nbnxn_atomdata_reduce_reals_simd
1493 #else
1494 nbnxn_atomdata_reduce_reals
1495 #endif
1496 (nbat->out[0].f,
1497 bitmask_is_set(flags->flag[b], 0),
1498 fptr, nfptr,
1499 i0, i1);
1501 else if (!bitmask_is_set(flags->flag[b], 0))
1503 nbnxn_atomdata_clear_reals(nbat->out[0].f,
1504 i0, i1);
1508 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1512 /* Add the force array(s) from nbnxn_atomdata_t to f */
1513 void nbnxn_atomdata_add_nbat_f_to_f(const nbnxn_search_t nbs,
1514 int locality,
1515 const nbnxn_atomdata_t *nbat,
1516 rvec *f,
1517 gmx_wallcycle *wcycle)
1519 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1520 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1522 int a0 = 0, na = 0;
1524 nbs_cycle_start(&nbs->cc[enbsCCreducef]);
1526 switch (locality)
1528 case eatAll:
1529 a0 = 0;
1530 na = nbs->natoms_nonlocal;
1531 break;
1532 case eatLocal:
1533 a0 = 0;
1534 na = nbs->natoms_local;
1535 break;
1536 case eatNonlocal:
1537 a0 = nbs->natoms_local;
1538 na = nbs->natoms_nonlocal - nbs->natoms_local;
1539 break;
1542 int nth = gmx_omp_nthreads_get(emntNonbonded);
1544 if (nbat->nout > 1)
1546 if (locality != eatAll)
1548 gmx_incons("add_f_to_f called with nout>1 and locality!=eatAll");
1551 /* Reduce the force thread output buffers into buffer 0, before adding
1552 * them to the, differently ordered, "real" force buffer.
1554 if (nbat->bUseTreeReduce)
1556 nbnxn_atomdata_add_nbat_f_to_f_treereduce(nbat, nth);
1558 else
1560 nbnxn_atomdata_add_nbat_f_to_f_stdreduce(nbat, nth);
1563 #pragma omp parallel for num_threads(nth) schedule(static)
1564 for (int th = 0; th < nth; th++)
1568 nbnxn_atomdata_add_nbat_f_to_f_part(nbs, nbat,
1569 nbat->out,
1571 a0+((th+0)*na)/nth,
1572 a0+((th+1)*na)/nth,
1575 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1578 nbs_cycle_stop(&nbs->cc[enbsCCreducef]);
1580 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1581 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1584 /* Adds the shift forces from nbnxn_atomdata_t to fshift */
1585 void nbnxn_atomdata_add_nbat_fshift_to_fshift(const nbnxn_atomdata_t *nbat,
1586 rvec *fshift)
1588 const nbnxn_atomdata_output_t * out = nbat->out;
1590 for (int s = 0; s < SHIFTS; s++)
1592 rvec sum;
1593 clear_rvec(sum);
1594 for (int th = 0; th < nbat->nout; th++)
1596 sum[XX] += out[th].fshift[s*DIM+XX];
1597 sum[YY] += out[th].fshift[s*DIM+YY];
1598 sum[ZZ] += out[th].fshift[s*DIM+ZZ];
1600 rvec_inc(fshift[s], sum);