Use const refs in init_nb_verlet
[alexxy/gromacs.git] / src / gromacs / nbnxm / pairlist_tuning.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief Implements functions for tuning adjustable parameters for the nbnxn non-bonded search and interaction kernels
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_nbnxm
42  */
43
44 #include "gmxpre.h"
45
46 #include "pairlist_tuning.h"
47
48 #include <cassert>
49 #include <cmath>
50 #include <cstdlib>
51
52 #include <algorithm>
53 #include <string>
54
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/multipletimestepping.h"
63 #include "gromacs/mdtypes/state.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"
72
73 #include "nbnxm_geometry.h"
74 #include "pairlistsets.h"
75
76 /*! \brief Returns if we can (heuristically) change nstlist and rlist
77  *
78  * \param [in] ir  The input parameter record
79  */
80 static bool supportsDynamicPairlistGenerationInterval(const t_inputrec& ir)
81 {
82     return ir.cutoff_scheme == ecutsVERLET && EI_DYNAMICS(ir.eI)
83            && !(EI_MD(ir.eI) && ir.etc == TemperatureCoupling::No) && ir.verletbuf_tol > 0;
84 }
85
86 /*! \brief Cost of non-bonded kernels
87  *
88  * We determine the extra cost of the non-bonded kernels compared to
89  * a reference nstlist value of 10 (which is the default in grompp).
90  */
91 static const int nbnxnReferenceNstlist = 10;
92 //! The values to try when switching
93 const int nstlist_try[] = { 20, 25, 40, 50, 80, 100 };
94 //! Number of elements in the neighborsearch list trials.
95 #define NNSTL (sizeof(nstlist_try) / sizeof(nstlist_try[0]))
96 /* Increase nstlist until the size of the pair-list increased by
97  * \p c_nbnxnListSizeFactor??? or more, but never more than
98  * \p c_nbnxnListSizeFactor??? + \p c_nbnxnListSizeFactorMargin.
99  * Since we have dynamic pair list pruning, the force kernel cost depends
100  * only very weakly on nstlist. It depends strongly on nstlistPrune.
101  * Increasing nstlist mainly affects the cost of the pair search (down due
102  * to lower frequency, up due to larger list) and the list pruning kernel.
103  * We increase nstlist conservatively with regard to kernel performance.
104  * In serial the search cost is not high and thus we don't gain much by
105  * increasing nstlist a lot. In parallel the MPI and CPU-GPU communication
106  * volume as well as the communication buffer preparation and reduction time
107  * increase quickly with rlist and thus nslist. Therefore we should avoid
108  * large nstlist, even if that also reduces the domain decomposition cost.
109  * With GPUs we perform the dynamic pruning in a rolling fashion and this
110  * overlaps with the update on the CPU, which allows even larger nstlist.
111  */
112 // CPU: pair-search is a factor ~1.5 slower than the non-bonded kernel.
113 //! Target pair-list size increase ratio for CPU
114 static const float c_nbnxnListSizeFactorCpu = 1.25;
115 // Intel KNL: pair-search is a factor ~2-3 slower than the non-bonded kernel.
116 //! Target pair-list size increase ratio for Intel KNL
117 static const float c_nbnxnListSizeFactorIntelXeonPhi = 1.4;
118 // GPU: pair-search is a factor 1.5-3 slower than the non-bonded kernel.
119 //! Target pair-list size increase ratio for GPU
120 static const float c_nbnxnListSizeFactorGPU = 1.4;
121 //! Never increase the size of the pair-list more than the factor above plus this margin
122 static const float c_nbnxnListSizeFactorMargin = 0.1;
123
124 void increaseNstlist(FILE*               fp,
125                      t_commrec*          cr,
126                      t_inputrec*         ir,
127                      int                 nstlist_cmdline,
128                      const gmx_mtop_t*   mtop,
129                      const matrix        box,
130                      bool                useOrEmulateGpuForNonbondeds,
131                      const gmx::CpuInfo& cpuinfo)
132 {
133     if (!EI_DYNAMICS(ir->eI))
134     {
135         /* Can only increase nstlist with dynamics */
136         return;
137     }
138
139     size_t      nstlist_ind = 0;
140     const char* nstl_gpu =
141             "\nFor optimal performance with a GPU nstlist (now %d) should be larger.\nThe "
142             "optimum depends on your CPU and GPU resources.\nYou might want to try several "
143             "nstlist values.\n";
144     const char* nve_err = "Can not increase nstlist because an NVE ensemble is used";
145     const char* vbd_err =
146             "Can not increase nstlist because verlet-buffer-tolerance is not set or used";
147     const char* box_err = "Can not increase nstlist because the box is too small";
148     const char* dd_err  = "Can not increase nstlist because of domain decomposition limitations";
149     char        buf[STRLEN];
150
151     /* When most of the computation, and in particular the non-bondeds is only
152      * performed every ir->mtsFactor steps due to multiple time stepping,
153      * we scale all nstlist values by this factor.
154      */
155     const int mtsFactor = gmx::nonbondedMtsFactor(*ir);
156
157     if (nstlist_cmdline <= 0)
158     {
159         if (ir->nstlist <= mtsFactor)
160         {
161             /* The user probably set nstlist<=mtsFactor for a reason,
162              * don't mess with the settings, except when < mtsFactor.
163              */
164             ir->nstlist = mtsFactor;
165
166             return;
167         }
168
169         /* With a GPU and fixed nstlist suggest tuning nstlist */
170         if (fp != nullptr && useOrEmulateGpuForNonbondeds && ir->nstlist < nstlist_try[0] * mtsFactor
171             && !supportsDynamicPairlistGenerationInterval(*ir))
172         {
173             fprintf(fp, nstl_gpu, ir->nstlist);
174         }
175
176         nstlist_ind = 0;
177         while (nstlist_ind < NNSTL && ir->nstlist >= nstlist_try[nstlist_ind] * mtsFactor)
178         {
179             nstlist_ind++;
180         }
181         if (nstlist_ind == NNSTL)
182         {
183             /* There are no larger nstlist value to try */
184             return;
185         }
186     }
187
188     if (EI_MD(ir->eI) && ir->etc == TemperatureCoupling::No)
189     {
190         if (MASTER(cr))
191         {
192             fprintf(stderr, "%s\n", nve_err);
193         }
194         if (fp != nullptr)
195         {
196             fprintf(fp, "%s\n", nve_err);
197         }
198
199         return;
200     }
201
202     if (ir->verletbuf_tol == 0 && useOrEmulateGpuForNonbondeds)
203     {
204         gmx_fatal(FARGS,
205                   "You are using an old tpr file with a GPU, please generate a new tpr file with "
206                   "an up to date version of grompp");
207     }
208
209     if (ir->verletbuf_tol < 0)
210     {
211         if (MASTER(cr))
212         {
213             fprintf(stderr, "%s\n", vbd_err);
214         }
215         if (fp != nullptr)
216         {
217             fprintf(fp, "%s\n", vbd_err);
218         }
219
220         return;
221     }
222
223     GMX_RELEASE_ASSERT(supportsDynamicPairlistGenerationInterval(*ir),
224                        "In all cases that do not support dynamic nstlist, we should have returned "
225                        "with an appropriate message above");
226
227     const bool  runningOnXeonPhi = (cpuinfo.brandString().find("Xeon Phi") != std::string::npos);
228     const float listfac_ok       = useOrEmulateGpuForNonbondeds
229                                      ? c_nbnxnListSizeFactorGPU
230                                      : runningOnXeonPhi ? c_nbnxnListSizeFactorIntelXeonPhi
231                                                         : c_nbnxnListSizeFactorCpu;
232     float listfac_max = listfac_ok + c_nbnxnListSizeFactorMargin;
233
234     const int nstlist_orig = ir->nstlist;
235     if (nstlist_cmdline > 0)
236     {
237         if (fp)
238         {
239             sprintf(buf, "Getting nstlist=%d from command line option", nstlist_cmdline);
240         }
241         ir->nstlist = nstlist_cmdline;
242     }
243
244     ListSetupType listType =
245             (useOrEmulateGpuForNonbondeds ? ListSetupType::Gpu : ListSetupType::CpuSimdWhenSupported);
246     VerletbufListSetup listSetup = verletbufGetSafeListSetup(listType);
247
248     /* Allow rlist to make the list a given factor larger than the list
249      * would be with the reference value for nstlist (10*mtsFactor).
250      */
251     int nstlist_prev = ir->nstlist;
252     ir->nstlist      = nbnxnReferenceNstlist * mtsFactor;
253     const real rlistWithReferenceNstlist =
254             calcVerletBufferSize(*mtop, det(box), *ir, ir->nstlist, ir->nstlist - 1, -1, listSetup);
255     ir->nstlist = nstlist_prev;
256
257     /* Determine the pair list size increase due to zero interactions */
258     real rlist_inc = nbnxn_get_rlist_effective_inc(listSetup.cluster_size_j, mtop->natoms / det(box));
259     real rlist_ok  = (rlistWithReferenceNstlist + rlist_inc) * std::cbrt(listfac_ok) - rlist_inc;
260     real rlist_max = (rlistWithReferenceNstlist + rlist_inc) * std::cbrt(listfac_max) - rlist_inc;
261     if (debug)
262     {
263         fprintf(debug, "nstlist tuning: rlist_inc %.3f rlist_ok %.3f rlist_max %.3f\n", rlist_inc, rlist_ok, rlist_max);
264     }
265
266     nstlist_prev    = nstlist_orig;
267     real rlist_prev = ir->rlist;
268     real rlist_new  = 0;
269     bool bBox = false, bDD = false, bCont = false;
270     do
271     {
272         if (nstlist_cmdline <= 0)
273         {
274             ir->nstlist = nstlist_try[nstlist_ind] * mtsFactor;
275         }
276
277         /* Set the pair-list buffer size in ir */
278         rlist_new = calcVerletBufferSize(
279                 *mtop, det(box), *ir, ir->nstlist, ir->nstlist - mtsFactor, -1, listSetup);
280
281         /* Does rlist fit in the box? */
282         bBox = (gmx::square(rlist_new) < max_cutoff2(ir->pbcType, box));
283         bDD  = true;
284         if (bBox && DOMAINDECOMP(cr))
285         {
286             /* Currently (as of July 2020), the code in this if clause is never executed.
287              * increaseNstlist(...) is only called from prepare_verlet_scheme, which in turns
288              * gets called by the runner _before_ setting up DD. DOMAINDECOMP(cr) will therefore
289              * always be false here. See #3334.
290              */
291             /* Check if rlist fits in the domain decomposition */
292             if (inputrec2nboundeddim(ir) < DIM)
293             {
294                 gmx_incons(
295                         "Changing nstlist with domain decomposition and unbounded dimensions is "
296                         "not implemented yet");
297             }
298             bDD = change_dd_cutoff(cr, box, gmx::ArrayRef<const gmx::RVec>(), rlist_new);
299         }
300
301         if (debug)
302         {
303             fprintf(debug,
304                     "nstlist %d rlist %.3f bBox %s bDD %s\n",
305                     ir->nstlist,
306                     rlist_new,
307                     gmx::boolToString(bBox),
308                     gmx::boolToString(bDD));
309         }
310
311         bCont = false;
312
313         if (nstlist_cmdline <= 0)
314         {
315             if (bBox && bDD && rlist_new <= rlist_max)
316             {
317                 /* Increase nstlist */
318                 nstlist_prev = ir->nstlist;
319                 rlist_prev   = rlist_new;
320                 bCont        = (nstlist_ind + 1 < NNSTL && rlist_new < rlist_ok);
321             }
322             else
323             {
324                 /* Stick with the previous nstlist */
325                 ir->nstlist = nstlist_prev;
326                 rlist_new   = rlist_prev;
327                 bBox        = true;
328                 bDD         = true;
329             }
330         }
331
332         nstlist_ind++;
333     } while (bCont);
334
335     if (!bBox || !bDD)
336     {
337         gmx_warning("%s", !bBox ? box_err : dd_err);
338         if (fp != nullptr)
339         {
340             fprintf(fp, "\n%s\n", !bBox ? box_err : dd_err);
341         }
342         ir->nstlist = nstlist_orig;
343     }
344     else if (ir->nstlist != nstlist_orig || rlist_new != ir->rlist)
345     {
346         sprintf(buf,
347                 "Changing nstlist from %d to %d, rlist from %g to %g",
348                 nstlist_orig,
349                 ir->nstlist,
350                 ir->rlist,
351                 rlist_new);
352         if (MASTER(cr))
353         {
354             fprintf(stderr, "%s\n\n", buf);
355         }
356         if (fp != nullptr)
357         {
358             fprintf(fp, "%s\n\n", buf);
359         }
360         ir->rlist = rlist_new;
361     }
362 }
363
364 /*! \brief The interval in steps at which we perform dynamic, rolling pruning on a GPU.
365  *
366  * Ideally we should auto-tune this value.
367  * Not considering overheads, 1 would be the ideal value. But 2 seems
368  * a reasonable compromise that reduces GPU kernel launch overheads and
369  * also avoids inefficiency on large GPUs when pruning small lists.
370  * Because with domain decomposition we alternate local/non-local pruning
371  * at even/odd steps, which gives a period of 2, this value currenly needs
372  * to be 2, which is indirectly asserted when the GPU pruning is dispatched
373  * during the force evaluation.
374  */
375 static const int c_nbnxnGpuRollingListPruningInterval = 2;
376
377 /*! \brief The minimum nstlist for dynamic pair list pruning.
378  *
379  * In most cases going lower than 4 will lead to a too high pruning cost.
380  * This value should be a multiple of \p c_nbnxnGpuRollingListPruningInterval
381  */
382 static const int c_nbnxnDynamicListPruningMinLifetime = 4;
383
384 /*! \brief Set the dynamic pairlist pruning parameters in \p ic
385  *
386  * \param[in]     inputrec          The input parameter record
387  * \param[in]     mtop        The global topology
388  * \param[in]     box         The unit cell
389  * \param[in]     useGpuList  Tells if we are using a GPU type pairlist
390  * \param[in]     listSetup   The nbnxn pair list setup
391  * \param[in]     userSetNstlistPrune  The user set ic->nstlistPrune (using an env.var.)
392  * \param[in]     interactionConst              The nonbonded interactions constants
393  * \param[in,out] listParams  The list setup parameters
394  */
395 static void setDynamicPairlistPruningParameters(const t_inputrec&          inputrec,
396                                                 const gmx_mtop_t&          mtop,
397                                                 const matrix               box,
398                                                 const bool                 useGpuList,
399                                                 const VerletbufListSetup&  listSetup,
400                                                 const bool                 userSetNstlistPrune,
401                                                 const interaction_const_t& interactionConst,
402                                                 PairlistParams*            listParams)
403 {
404     /* When applying multiple time stepping to the non-bonded forces,
405      * we only compute them every mtsFactor steps, so all parameters here
406      * should be a multiple of mtsFactor.
407      */
408     listParams->mtsFactor = gmx::nonbondedMtsFactor(inputrec);
409
410     const int mtsFactor = listParams->mtsFactor;
411
412     GMX_RELEASE_ASSERT(inputrec.nstlist % mtsFactor == 0,
413                        "nstlist should be a multiple of mtsFactor");
414
415     listParams->lifetime = inputrec.nstlist - mtsFactor;
416
417     /* When nstlistPrune was set by the user, we need to execute one loop
418      * iteration to determine rlistInner.
419      * Otherwise we compute rlistInner and increase nstlist as long as
420      * we have a pairlist buffer of length 0 (i.e. rlistInner == cutoff).
421      */
422     const real interactionCutoff = std::max(interactionConst.rcoulomb, interactionConst.rvdw);
423     int        tunedNstlistPrune = listParams->nstlistPrune;
424     do
425     {
426         /* Dynamic pruning on the GPU is performed on the list for
427          * the next step on the coordinates of the current step,
428          * so the list lifetime is nstlistPrune (not the usual nstlist-mtsFactor).
429          */
430         int listLifetime         = tunedNstlistPrune - (useGpuList ? 0 : mtsFactor);
431         listParams->nstlistPrune = tunedNstlistPrune;
432         listParams->rlistInner   = calcVerletBufferSize(
433                 mtop, det(box), inputrec, tunedNstlistPrune, listLifetime, -1, listSetup);
434
435         /* On the GPU we apply the dynamic pruning in a rolling fashion
436          * every c_nbnxnGpuRollingListPruningInterval steps,
437          * so keep nstlistPrune a multiple of the interval.
438          */
439         tunedNstlistPrune += (useGpuList ? c_nbnxnGpuRollingListPruningInterval : 1) * mtsFactor;
440     } while (!userSetNstlistPrune && tunedNstlistPrune < inputrec.nstlist
441              && listParams->rlistInner == interactionCutoff);
442
443     if (userSetNstlistPrune)
444     {
445         listParams->useDynamicPruning = true;
446     }
447     else
448     {
449         /* Determine the pair list size increase due to zero interactions */
450         real rlistInc = nbnxn_get_rlist_effective_inc(listSetup.cluster_size_j, mtop.natoms / det(box));
451
452         /* Dynamic pruning is only useful when the inner list is smaller than
453          * the outer. The factor 0.99 ensures at least 3% list size reduction.
454          *
455          * With dynamic pruning on the CPU we prune after updating,
456          * so nstlistPrune=nstlist-1 would add useless extra work.
457          * With the GPU there will probably be more overhead than gain
458          * with nstlistPrune=nstlist-1, so we disable dynamic pruning.
459          * Note that in such cases the first sub-condition is likely also false.
460          */
461         listParams->useDynamicPruning =
462                 (listParams->rlistInner + rlistInc < 0.99 * (listParams->rlistOuter + rlistInc)
463                  && listParams->nstlistPrune < listParams->lifetime);
464     }
465
466     if (!listParams->useDynamicPruning)
467     {
468         /* These parameters should not be used, but set them to useful values */
469         listParams->nstlistPrune = -1;
470         listParams->rlistInner   = listParams->rlistOuter;
471     }
472 }
473
474 /*! \brief Returns a string describing the setup of a single pair-list
475  *
476  * \param[in] listName           Short name of the list, can be ""
477  * \param[in] nstList            The list update interval in steps
478  * \param[in] nstListForSpacing  Update interval for setting the number characters for printing \p nstList
479  * \param[in] rList              List cut-off radius
480  * \param[in] interactionCutoff  The interaction cut-off, use for printing the list buffer size
481  */
482 static std::string formatListSetup(const std::string& listName,
483                                    int                nstList,
484                                    int                nstListForSpacing,
485                                    real               rList,
486                                    real               interactionCutoff)
487 {
488     std::string listSetup = "  ";
489     if (!listName.empty())
490     {
491         listSetup += listName + " list: ";
492     }
493     listSetup += "updated every ";
494     // Make the shortest int format string that fits nstListForSpacing
495     std::string nstListFormat =
496             "%" + gmx::formatString("%zu", gmx::formatString("%d", nstListForSpacing).size()) + "d";
497     listSetup += gmx::formatString(nstListFormat.c_str(), nstList);
498     listSetup += gmx::formatString(
499             " steps, buffer %.3f nm, rlist %.3f nm\n", rList - interactionCutoff, rList);
500
501     return listSetup;
502 }
503
504 void setupDynamicPairlistPruning(const gmx::MDLogger&       mdlog,
505                                  const t_inputrec&          inputrec,
506                                  const gmx_mtop_t&          mtop,
507                                  matrix                     box,
508                                  const interaction_const_t& interactionConst,
509                                  PairlistParams*            listParams)
510 {
511     GMX_RELEASE_ASSERT(listParams->rlistOuter > 0, "With the nbnxn setup rlist should be > 0");
512
513     /* Initialize the parameters to no dynamic list pruning */
514     listParams->useDynamicPruning = false;
515
516     const VerletbufListSetup ls = { IClusterSizePerListType[listParams->pairlistType],
517                                     JClusterSizePerListType[listParams->pairlistType] };
518
519     /* Currently emulation mode does not support dual pair-lists */
520     const bool useGpuList = sc_isGpuPairListType[listParams->pairlistType];
521
522     if (supportsDynamicPairlistGenerationInterval(inputrec) && getenv("GMX_DISABLE_DYNAMICPRUNING") == nullptr)
523     {
524         /* Note that nstlistPrune can have any value independently of nstlist.
525          * Actually applying rolling pruning is only useful when
526          * nstlistPrune < nstlist -1
527          */
528         char* env                 = getenv("GMX_NSTLIST_DYNAMICPRUNING");
529         bool  userSetNstlistPrune = (env != nullptr);
530
531         if (userSetNstlistPrune)
532         {
533             char* end                = nullptr;
534             listParams->nstlistPrune = strtol(env, &end, 10);
535             if (!end || (*end != 0)
536                 || !(listParams->nstlistPrune > 0 && listParams->nstlistPrune < inputrec.nstlist))
537             {
538                 gmx_fatal(FARGS,
539                           "Invalid value passed in GMX_NSTLIST_DYNAMICPRUNING=%s, should be > 0 "
540                           "and < nstlist",
541                           env);
542             }
543         }
544         else
545         {
546             static_assert(c_nbnxnDynamicListPruningMinLifetime % c_nbnxnGpuRollingListPruningInterval == 0,
547                           "c_nbnxnDynamicListPruningMinLifetime sets the starting value for "
548                           "nstlistPrune, which should be divisible by the rolling pruning interval "
549                           "for efficiency reasons.");
550
551             // TODO: Use auto-tuning to determine nstlistPrune
552             listParams->nstlistPrune = c_nbnxnDynamicListPruningMinLifetime;
553         }
554
555         setDynamicPairlistPruningParameters(
556                 inputrec, mtop, box, useGpuList, ls, userSetNstlistPrune, interactionConst, listParams);
557
558         if (listParams->useDynamicPruning && useGpuList)
559         {
560             /* Note that we can round down here. This makes the effective
561              * rolling pruning interval slightly shorter than nstlistTune,
562              * thus giving correct results, but a slightly lower efficiency.
563              */
564             GMX_RELEASE_ASSERT(listParams->nstlistPrune >= c_nbnxnGpuRollingListPruningInterval,
565                                ("With dynamic list pruning on GPUs pruning frequency must be at "
566                                 "least as large as the rolling pruning interval ("
567                                 + std::to_string(c_nbnxnGpuRollingListPruningInterval) + ").")
568                                        .c_str());
569             listParams->numRollingPruningParts =
570                     listParams->nstlistPrune / c_nbnxnGpuRollingListPruningInterval;
571         }
572         else
573         {
574             listParams->numRollingPruningParts = 1;
575         }
576     }
577
578     std::string mesg;
579
580     const real interactionCutoff = std::max(interactionConst.rcoulomb, interactionConst.rvdw);
581     if (listParams->useDynamicPruning)
582     {
583         mesg += gmx::formatString(
584                 "Using a dual %dx%d pair-list setup updated with dynamic%s pruning:\n",
585                 ls.cluster_size_i,
586                 ls.cluster_size_j,
587                 listParams->numRollingPruningParts > 1 ? ", rolling" : "");
588         mesg += formatListSetup(
589                 "outer", inputrec.nstlist, inputrec.nstlist, listParams->rlistOuter, interactionCutoff);
590         mesg += formatListSetup(
591                 "inner", listParams->nstlistPrune, inputrec.nstlist, listParams->rlistInner, interactionCutoff);
592     }
593     else
594     {
595         mesg += gmx::formatString("Using a %dx%d pair-list setup:\n", ls.cluster_size_i, ls.cluster_size_j);
596         mesg += formatListSetup(
597                 "", inputrec.nstlist, inputrec.nstlist, listParams->rlistOuter, interactionCutoff);
598     }
599     if (supportsDynamicPairlistGenerationInterval(inputrec))
600     {
601         const VerletbufListSetup listSetup1x1 = { 1, 1 };
602         const real               rlistOuter   = calcVerletBufferSize(
603                 mtop, det(box), inputrec, inputrec.nstlist, inputrec.nstlist - 1, -1, listSetup1x1);
604         real rlistInner = rlistOuter;
605         if (listParams->useDynamicPruning)
606         {
607             int listLifeTime = listParams->nstlistPrune - (useGpuList ? 0 : 1);
608             rlistInner       = calcVerletBufferSize(
609                     mtop, det(box), inputrec, listParams->nstlistPrune, listLifeTime, -1, listSetup1x1);
610         }
611
612         mesg += gmx::formatString(
613                 "At tolerance %g kJ/mol/ps per atom, equivalent classical 1x1 list would be:\n",
614                 inputrec.verletbuf_tol);
615         if (listParams->useDynamicPruning)
616         {
617             mesg += formatListSetup(
618                     "outer", inputrec.nstlist, inputrec.nstlist, rlistOuter, interactionCutoff);
619             mesg += formatListSetup(
620                     "inner", listParams->nstlistPrune, inputrec.nstlist, rlistInner, interactionCutoff);
621         }
622         else
623         {
624             mesg += formatListSetup("", inputrec.nstlist, inputrec.nstlist, rlistOuter, interactionCutoff);
625         }
626     }
627
628     GMX_LOG(mdlog.info).asParagraph().appendText(mesg);
629 }