2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 * \brief Declares the PME load balancing helper for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "pmeloadbalancehelper.h"
46 #include "gromacs/ewald/pme_load_balancing.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/forcerec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/interaction_const.h"
51 #include "gromacs/mdtypes/mdrunoptions.h"
52 #include "gromacs/mdtypes/state.h"
53 #include "gromacs/nbnxm/nbnxm.h"
55 #include "statepropagatordata.h"
59 bool PmeLoadBalanceHelper::doPmeLoadBalancing(const MdrunOptions& mdrunOptions,
60 const t_inputrec* inputrec,
63 return (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
64 && inputrec->cutoff_scheme != CutoffScheme::Group);
67 PmeLoadBalanceHelper::PmeLoadBalanceHelper(bool isVerbose,
68 StatePropagatorData* statePropagatorData,
71 const MDLogger& mdlog,
72 const t_inputrec* inputrec,
73 gmx_wallcycle* wcycle,
75 pme_loadbal_(nullptr),
77 isVerbose_(isVerbose),
78 bPMETunePrinting_(false),
79 statePropagatorData_(statePropagatorData),
89 void PmeLoadBalanceHelper::setup()
91 auto box = statePropagatorData_->constBox();
92 GMX_RELEASE_ASSERT(box[0][0] != 0 && box[1][1] != 0 && box[2][2] != 0,
93 "PmeLoadBalanceHelper cannot be initialized with zero box.");
95 &pme_loadbal_, cr_, mdlog_, *inputrec_, box, *fr_->ic, *fr_->nbv, fr_->pmedata, fr_->nbv->useGpu());
98 void PmeLoadBalanceHelper::run(gmx::Step step, gmx::Time gmx_unused time)
100 if (step != nextNSStep_ || step == inputrec_->init_step)
105 // PME grid + cut-off optimization with GPUs or PME nodes
106 // TODO pass SimulationWork object into this function, such that last argument can be set as
107 // simulationWork.useGpuPmePpCommunication as is done in main MD loop.
108 pme_loadbal_do(pme_loadbal_,
110 (isVerbose_ && MASTER(cr_)) ? stderr : nullptr,
115 statePropagatorData_->constBox(),
116 statePropagatorData_->constPositionsView().paddedArrayRef(),
119 step - inputrec_->init_step,
124 void PmeLoadBalanceHelper::teardown()
126 pme_loadbal_done(pme_loadbal_, fplog_, mdlog_, fr_->nbv->useGpu());
129 bool PmeLoadBalanceHelper::pmePrinting()
131 return bPMETunePrinting_;
134 const pme_load_balancing_t* PmeLoadBalanceHelper::loadBalancingObject()
139 std::optional<SignallerCallback> PmeLoadBalanceHelper::registerNSCallback()
141 return [this](Step step, Time gmx_unused time) { nextNSStep_ = step; };