2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,2018,2019, 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.
38 * \brief Implements functions for tuning adjustable parameters for the nbnxn non-bonded search and interaction kernels
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup __module_nb_verlet
46 #include "pairlist_tuning.h"
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/hardware/cpuinfo.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/mdlib/calc_verletbuf.h"
59 #include "gromacs/mdtypes/commrec.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/interaction_const.h"
62 #include "gromacs/mdtypes/state.h"
63 #include "gromacs/nbnxm/nbnxm_geometry.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/topology/topology.h"
66 #include "gromacs/utility/cstringutil.h"
67 #include "gromacs/utility/fatalerror.h"
68 #include "gromacs/utility/gmxassert.h"
69 #include "gromacs/utility/logger.h"
70 #include "gromacs/utility/strconvert.h"
71 #include "gromacs/utility/stringutil.h"
73 #include "pairlistsets.h"
75 /*! \brief Returns if we can (heuristically) change nstlist and rlist
77 * \param [in] ir The input parameter record
79 static bool supportsDynamicPairlistGenerationInterval(const t_inputrec
&ir
)
82 ir
.cutoff_scheme
== ecutsVERLET
&&
84 !(EI_MD(ir
.eI
) && ir
.etc
== etcNO
) &&
88 /*! \brief Cost of non-bonded kernels
90 * We determine the extra cost of the non-bonded kernels compared to
91 * a reference nstlist value of 10 (which is the default in grompp).
93 static const int nbnxnReferenceNstlist
= 10;
94 //! The values to try when switching
95 const int nstlist_try
[] = { 20, 25, 40, 50, 80, 100 };
96 //! Number of elements in the neighborsearch list trials.
97 #define NNSTL (sizeof(nstlist_try)/sizeof(nstlist_try[0]))
98 /* Increase nstlist until the size of the pair-list increased by
99 * \p c_nbnxnListSizeFactor??? or more, but never more than
100 * \p c_nbnxnListSizeFactor??? + \p c_nbnxnListSizeFactorMargin.
101 * Since we have dynamic pair list pruning, the force kernel cost depends
102 * only very weakly on nstlist. It depends strongly on nstlistPrune.
103 * Increasing nstlist mainly affects the cost of the pair search (down due
104 * to lower frequency, up due to larger list) and the list pruning kernel.
105 * We increase nstlist conservatively with regard to kernel performance.
106 * In serial the search cost is not high and thus we don't gain much by
107 * increasing nstlist a lot. In parallel the MPI and CPU-GPU communication
108 * volume as well as the communication buffer preparation and reduction time
109 * increase quickly with rlist and thus nslist. Therefore we should avoid
110 * large nstlist, even if that also reduces the domain decomposition cost.
111 * With GPUs we perform the dynamic pruning in a rolling fashion and this
112 * overlaps with the update on the CPU, which allows even larger nstlist.
114 // CPU: pair-search is a factor ~1.5 slower than the non-bonded kernel.
115 //! Target pair-list size increase ratio for CPU
116 static const float c_nbnxnListSizeFactorCpu
= 1.25;
117 // Intel KNL: pair-search is a factor ~2-3 slower than the non-bonded kernel.
118 //! Target pair-list size increase ratio for Intel KNL
119 static const float c_nbnxnListSizeFactorIntelXeonPhi
= 1.4;
120 // GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel.
121 //! Target pair-list size increase ratio for GPU
122 static const float c_nbnxnListSizeFactorGPU
= 1.4;
123 //! Never increase the size of the pair-list more than the factor above plus this margin
124 static const float c_nbnxnListSizeFactorMargin
= 0.1;
126 void increaseNstlist(FILE *fp
, t_commrec
*cr
,
127 t_inputrec
*ir
, int nstlist_cmdline
,
128 const gmx_mtop_t
*mtop
,
130 bool useOrEmulateGpuForNonbondeds
,
131 const gmx::CpuInfo
&cpuinfo
)
133 if (!EI_DYNAMICS(ir
->eI
))
135 /* Can only increase nstlist with dynamics */
139 float listfac_ok
, listfac_max
;
140 int nstlist_orig
, nstlist_prev
;
141 real rlist_inc
, rlist_ok
, rlist_max
;
142 real rlist_new
, rlist_prev
;
143 size_t nstlist_ind
= 0;
144 gmx_bool bBox
, bDD
, bCont
;
145 const char *nstl_gpu
= "\nFor optimal performance with a GPU nstlist (now %d) should be larger.\nThe optimum depends on your CPU and GPU resources.\nYou might want to try several nstlist values.\n";
146 const char *nve_err
= "Can not increase nstlist because an NVE ensemble is used";
147 const char *vbd_err
= "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
148 const char *box_err
= "Can not increase nstlist because the box is too small";
149 const char *dd_err
= "Can not increase nstlist because of domain decomposition limitations";
152 if (nstlist_cmdline
<= 0)
154 if (ir
->nstlist
== 1)
156 /* The user probably set nstlist=1 for a reason,
157 * don't mess with the settings.
162 /* With a GPU and fixed nstlist suggest tuning nstlist */
164 useOrEmulateGpuForNonbondeds
&&
165 ir
->nstlist
< nstlist_try
[0] &&
166 !supportsDynamicPairlistGenerationInterval(*ir
))
168 fprintf(fp
, nstl_gpu
, ir
->nstlist
);
172 while (nstlist_ind
< NNSTL
&& ir
->nstlist
>= nstlist_try
[nstlist_ind
])
176 if (nstlist_ind
== NNSTL
)
178 /* There are no larger nstlist value to try */
183 if (EI_MD(ir
->eI
) && ir
->etc
== etcNO
)
187 fprintf(stderr
, "%s\n", nve_err
);
191 fprintf(fp
, "%s\n", nve_err
);
197 if (ir
->verletbuf_tol
== 0 && useOrEmulateGpuForNonbondeds
)
199 gmx_fatal(FARGS
, "You are using an old tpr file with a GPU, please generate a new tpr file with an up to date version of grompp");
202 if (ir
->verletbuf_tol
< 0)
206 fprintf(stderr
, "%s\n", vbd_err
);
210 fprintf(fp
, "%s\n", vbd_err
);
216 GMX_RELEASE_ASSERT(supportsDynamicPairlistGenerationInterval(*ir
), "In all cases that do not support dynamic nstlist, we should have returned with an appropriate message above");
218 if (useOrEmulateGpuForNonbondeds
)
220 listfac_ok
= c_nbnxnListSizeFactorGPU
;
222 else if (cpuinfo
.brandString().find("Xeon Phi") != std::string::npos
)
224 listfac_ok
= c_nbnxnListSizeFactorIntelXeonPhi
;
228 listfac_ok
= c_nbnxnListSizeFactorCpu
;
230 listfac_max
= listfac_ok
+ c_nbnxnListSizeFactorMargin
;
232 nstlist_orig
= ir
->nstlist
;
233 if (nstlist_cmdline
> 0)
237 sprintf(buf
, "Getting nstlist=%d from command line option",
240 ir
->nstlist
= nstlist_cmdline
;
243 ListSetupType listType
= (useOrEmulateGpuForNonbondeds
? ListSetupType::Gpu
: ListSetupType::CpuSimdWhenSupported
);
244 VerletbufListSetup listSetup
= verletbufGetSafeListSetup(listType
);
246 /* Allow rlist to make the list a given factor larger than the list
247 * would be with the reference value for nstlist (10).
249 nstlist_prev
= ir
->nstlist
;
250 ir
->nstlist
= nbnxnReferenceNstlist
;
251 const real rlistWithReferenceNstlist
=
252 calcVerletBufferSize(*mtop
, det(box
), *ir
, ir
->nstlist
, ir
->nstlist
- 1,
254 ir
->nstlist
= nstlist_prev
;
256 /* Determine the pair list size increase due to zero interactions */
257 rlist_inc
= nbnxn_get_rlist_effective_inc(listSetup
.cluster_size_j
,
258 mtop
->natoms
/det(box
));
259 rlist_ok
= (rlistWithReferenceNstlist
+ rlist_inc
)*std::cbrt(listfac_ok
) - rlist_inc
;
260 rlist_max
= (rlistWithReferenceNstlist
+ rlist_inc
)*std::cbrt(listfac_max
) - rlist_inc
;
263 fprintf(debug
, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n",
264 rlist_inc
, rlist_ok
, rlist_max
);
267 nstlist_prev
= nstlist_orig
;
268 rlist_prev
= ir
->rlist
;
271 if (nstlist_cmdline
<= 0)
273 ir
->nstlist
= nstlist_try
[nstlist_ind
];
276 /* Set the pair-list buffer size in ir */
278 calcVerletBufferSize(*mtop
, det(box
), *ir
, ir
->nstlist
, ir
->nstlist
- 1, -1, listSetup
);
280 /* Does rlist fit in the box? */
281 bBox
= (gmx::square(rlist_new
) < max_cutoff2(ir
->ePBC
, box
));
283 if (bBox
&& DOMAINDECOMP(cr
))
285 /* Check if rlist fits in the domain decomposition */
286 if (inputrec2nboundeddim(ir
) < DIM
)
288 gmx_incons("Changing nstlist with domain decomposition and unbounded dimensions is not implemented yet");
291 copy_mat(box
, state_tmp
.box
);
292 bDD
= change_dd_cutoff(cr
, state_tmp
, rlist_new
);
297 fprintf(debug
, "nstlist %d rlist %.3f bBox %s bDD %s\n",
298 ir
->nstlist
, rlist_new
, gmx::boolToString(bBox
), gmx::boolToString(bDD
));
303 if (nstlist_cmdline
<= 0)
305 if (bBox
&& bDD
&& rlist_new
<= rlist_max
)
307 /* Increase nstlist */
308 nstlist_prev
= ir
->nstlist
;
309 rlist_prev
= rlist_new
;
310 bCont
= (nstlist_ind
+1 < NNSTL
&& rlist_new
< rlist_ok
);
314 /* Stick with the previous nstlist */
315 ir
->nstlist
= nstlist_prev
;
316 rlist_new
= rlist_prev
;
328 gmx_warning("%s", !bBox
? box_err
: dd_err
);
331 fprintf(fp
, "\n%s\n", !bBox
? box_err
: dd_err
);
333 ir
->nstlist
= nstlist_orig
;
335 else if (ir
->nstlist
!= nstlist_orig
|| rlist_new
!= ir
->rlist
)
337 sprintf(buf
, "Changing nstlist from %d to %d, rlist from %g to %g",
338 nstlist_orig
, ir
->nstlist
,
339 ir
->rlist
, rlist_new
);
342 fprintf(stderr
, "%s\n\n", buf
);
346 fprintf(fp
, "%s\n\n", buf
);
348 ir
->rlist
= rlist_new
;
352 /*! \brief The interval in steps at which we perform dynamic, rolling pruning on a GPU.
354 * Ideally we should auto-tune this value.
355 * Not considering overheads, 1 would be the ideal value. But 2 seems
356 * a reasonable compromise that reduces GPU kernel launch overheads and
357 * also avoids inefficiency on large GPUs when pruning small lists.
358 * Because with domain decomposition we alternate local/non-local pruning
359 * at even/odd steps, which gives a period of 2, this value currenly needs
360 * to be 2, which is indirectly asserted when the GPU pruning is dispatched
361 * during the force evaluation.
363 static const int c_nbnxnGpuRollingListPruningInterval
= 2;
365 /*! \brief The minimum nstlist for dynamic pair list pruning.
367 * In most cases going lower than 4 will lead to a too high pruning cost.
368 * This value should be a multiple of \p c_nbnxnGpuRollingListPruningInterval
370 static const int c_nbnxnDynamicListPruningMinLifetime
= 4;
372 /*! \brief Set the dynamic pairlist pruning parameters in \p ic
374 * \param[in] ir The input parameter record
375 * \param[in] mtop The global topology
376 * \param[in] box The unit cell
377 * \param[in] useGpuList Tells if we are using a GPU type pairlist
378 * \param[in] listSetup The nbnxn pair list setup
379 * \param[in] userSetNstlistPrune The user set ic->nstlistPrune (using an env.var.)
380 * \param[in] ic The nonbonded interactions constants
381 * \param[in,out] listParams The list setup parameters
384 setDynamicPairlistPruningParameters(const t_inputrec
*ir
,
385 const gmx_mtop_t
*mtop
,
387 const bool useGpuList
,
388 const VerletbufListSetup
&listSetup
,
389 const bool userSetNstlistPrune
,
390 const interaction_const_t
*ic
,
391 PairlistParams
*listParams
)
393 listParams
->lifetime
= ir
->nstlist
- 1;
395 /* When nstlistPrune was set by the user, we need to execute one loop
396 * iteration to determine rlistInner.
397 * Otherwise we compute rlistInner and increase nstlist as long as
398 * we have a pairlist buffer of length 0 (i.e. rlistInner == cutoff).
400 const real interactionCutoff
= std::max(ic
->rcoulomb
, ic
->rvdw
);
401 int tunedNstlistPrune
= listParams
->nstlistPrune
;
404 /* Dynamic pruning on the GPU is performed on the list for
405 * the next step on the coordinates of the current step,
406 * so the list lifetime is nstlistPrune (not the usual nstlist-1).
408 int listLifetime
= tunedNstlistPrune
- (useGpuList
? 0 : 1);
409 listParams
->nstlistPrune
= tunedNstlistPrune
;
410 listParams
->rlistInner
=
411 calcVerletBufferSize(*mtop
, det(box
), *ir
,
412 tunedNstlistPrune
, listLifetime
,
415 /* On the GPU we apply the dynamic pruning in a rolling fashion
416 * every c_nbnxnGpuRollingListPruningInterval steps,
417 * so keep nstlistPrune a multiple of the interval.
419 tunedNstlistPrune
+= useGpuList
? c_nbnxnGpuRollingListPruningInterval
: 1;
421 while (!userSetNstlistPrune
&&
422 tunedNstlistPrune
< ir
->nstlist
&&
423 listParams
->rlistInner
== interactionCutoff
);
425 if (userSetNstlistPrune
)
427 listParams
->useDynamicPruning
= true;
431 /* Determine the pair list size increase due to zero interactions */
432 real rlistInc
= nbnxn_get_rlist_effective_inc(listSetup
.cluster_size_j
,
433 mtop
->natoms
/det(box
));
435 /* Dynamic pruning is only useful when the inner list is smaller than
436 * the outer. The factor 0.99 ensures at least 3% list size reduction.
438 * With dynamic pruning on the CPU we prune after updating,
439 * so nstlistPrune=nstlist-1 would add useless extra work.
440 * With the GPU there will probably be more overhead than gain
441 * with nstlistPrune=nstlist-1, so we disable dynamic pruning.
442 * Note that in such cases the first sub-condition is likely also false.
444 listParams
->useDynamicPruning
=
445 (listParams
->rlistInner
+ rlistInc
< 0.99*(listParams
->rlistOuter
+ rlistInc
) &&
446 listParams
->nstlistPrune
< listParams
->lifetime
);
449 if (!listParams
->useDynamicPruning
)
451 /* These parameters should not be used, but set them to useful values */
452 listParams
->nstlistPrune
= -1;
453 listParams
->rlistInner
= listParams
->rlistOuter
;
457 /*! \brief Returns a string describing the setup of a single pair-list
459 * \param[in] listName Short name of the list, can be ""
460 * \param[in] nstList The list update interval in steps
461 * \param[in] nstListForSpacing Update interval for setting the number characters for printing \p nstList
462 * \param[in] rList List cut-off radius
463 * \param[in] interactionCutoff The interaction cut-off, use for printing the list buffer size
465 static std::string
formatListSetup(const std::string
&listName
,
467 int nstListForSpacing
,
469 real interactionCutoff
)
471 std::string listSetup
= " ";
472 if (!listName
.empty())
474 listSetup
+= listName
+ " list: ";
476 listSetup
+= "updated every ";
477 // Make the shortest int format string that fits nstListForSpacing
478 std::string nstListFormat
= "%" + gmx::formatString("%zu", gmx::formatString("%d", nstListForSpacing
).size()) + "d";
479 listSetup
+= gmx::formatString(nstListFormat
.c_str(), nstList
);
480 listSetup
+= gmx::formatString(" steps, buffer %.3f nm, rlist %.3f nm\n",
481 rList
- interactionCutoff
, rList
);
486 void setupDynamicPairlistPruning(const gmx::MDLogger
&mdlog
,
487 const t_inputrec
*ir
,
488 const gmx_mtop_t
*mtop
,
490 const interaction_const_t
*ic
,
491 PairlistParams
*listParams
)
493 GMX_RELEASE_ASSERT(listParams
->rlistOuter
> 0, "With the nbnxn setup rlist should be > 0");
495 /* Initialize the parameters to no dynamic list pruning */
496 listParams
->useDynamicPruning
= false;
498 const VerletbufListSetup ls
=
500 IClusterSizePerListType
[listParams
->pairlistType
],
501 JClusterSizePerListType
[listParams
->pairlistType
]
504 /* Currently emulation mode does not support dual pair-lists */
505 const bool useGpuList
= (listParams
->pairlistType
== PairlistType::HierarchicalNxN
);
507 if (supportsDynamicPairlistGenerationInterval(*ir
) &&
508 getenv("GMX_DISABLE_DYNAMICPRUNING") == nullptr)
510 /* Note that nstlistPrune can have any value independently of nstlist.
511 * Actually applying rolling pruning is only useful when
512 * nstlistPrune < nstlist -1
514 char *env
= getenv("GMX_NSTLIST_DYNAMICPRUNING");
515 bool userSetNstlistPrune
= (env
!= nullptr);
517 if (userSetNstlistPrune
)
520 listParams
->nstlistPrune
= strtol(env
, &end
, 10);
521 if (!end
|| (*end
!= 0) ||
522 !(listParams
->nstlistPrune
> 0 && listParams
->nstlistPrune
< ir
->nstlist
))
524 gmx_fatal(FARGS
, "Invalid value passed in GMX_NSTLIST_DYNAMICPRUNING=%s, should be > 0 and < nstlist", env
);
529 static_assert(c_nbnxnDynamicListPruningMinLifetime
% c_nbnxnGpuRollingListPruningInterval
== 0,
530 "c_nbnxnDynamicListPruningMinLifetime sets the starting value for nstlistPrune, which should be divisible by the rolling pruning interval for efficiency reasons.");
532 // TODO: Use auto-tuning to determine nstlistPrune
533 listParams
->nstlistPrune
= c_nbnxnDynamicListPruningMinLifetime
;
536 setDynamicPairlistPruningParameters(ir
, mtop
, box
, useGpuList
, ls
,
537 userSetNstlistPrune
, ic
,
540 if (listParams
->useDynamicPruning
&& useGpuList
)
542 /* Note that we can round down here. This makes the effective
543 * rolling pruning interval slightly shorter than nstlistTune,
544 * thus giving correct results, but a slightly lower efficiency.
546 GMX_RELEASE_ASSERT(listParams
->nstlistPrune
>= c_nbnxnGpuRollingListPruningInterval
,
547 ( "With dynamic list pruning on GPUs pruning frequency must be at least as large as the rolling pruning interval (" +
548 std::to_string(c_nbnxnGpuRollingListPruningInterval
) +
550 listParams
->numRollingPruningParts
= listParams
->nstlistPrune
/c_nbnxnGpuRollingListPruningInterval
;
554 listParams
->numRollingPruningParts
= 1;
560 const real interactionCutoff
= std::max(ic
->rcoulomb
, ic
->rvdw
);
561 if (listParams
->useDynamicPruning
)
563 mesg
+= gmx::formatString("Using a dual %dx%d pair-list setup updated with dynamic%s pruning:\n",
564 ls
.cluster_size_i
, ls
.cluster_size_j
,
565 listParams
->numRollingPruningParts
> 1 ? ", rolling" : "");
566 mesg
+= formatListSetup("outer", ir
->nstlist
, ir
->nstlist
, listParams
->rlistOuter
, interactionCutoff
);
567 mesg
+= formatListSetup("inner", listParams
->nstlistPrune
, ir
->nstlist
, listParams
->rlistInner
, interactionCutoff
);
571 mesg
+= gmx::formatString("Using a %dx%d pair-list setup:\n",
572 ls
.cluster_size_i
, ls
.cluster_size_j
);
573 mesg
+= formatListSetup("", ir
->nstlist
, ir
->nstlist
, listParams
->rlistOuter
, interactionCutoff
);
575 if (supportsDynamicPairlistGenerationInterval(*ir
))
577 const VerletbufListSetup listSetup1x1
= { 1, 1 };
578 const real rlistOuter
=
579 calcVerletBufferSize(*mtop
, det(box
), *ir
, ir
->nstlist
, ir
->nstlist
- 1,
581 real rlistInner
= rlistOuter
;
582 if (listParams
->useDynamicPruning
)
584 int listLifeTime
= listParams
->nstlistPrune
- (useGpuList
? 0 : 1);
586 calcVerletBufferSize(*mtop
, det(box
), *ir
, listParams
->nstlistPrune
, listLifeTime
,
590 mesg
+= gmx::formatString("At tolerance %g kJ/mol/ps per atom, equivalent classical 1x1 list would be:\n",
592 if (listParams
->useDynamicPruning
)
594 mesg
+= formatListSetup("outer", ir
->nstlist
, ir
->nstlist
, rlistOuter
, interactionCutoff
);
595 mesg
+= formatListSetup("inner", listParams
->nstlistPrune
, ir
->nstlist
, rlistInner
, interactionCutoff
);
599 mesg
+= formatListSetup("", ir
->nstlist
, ir
->nstlist
, rlistOuter
, interactionCutoff
);
603 GMX_LOG(mdlog
.info
).asParagraph().appendText(mesg
);