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 Defines the domain decomposition helper for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "domdechelper.h"
46 #include "gromacs/domdec/collect.h"
47 #include "gromacs/domdec/partition.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/state.h"
51 #include "gromacs/pbcutil/pbc.h"
53 #include "freeenergyperturbationdata.h"
54 #include "statepropagatordata.h"
55 #include "topologyholder.h"
59 DomDecHelper::DomDecHelper(bool isVerbose,
60 int verbosePrintInterval,
61 StatePropagatorData* statePropagatorData,
62 TopologyHolder* topologyHolder,
66 const MDLogger& mdlog,
68 const t_inputrec* inputrec,
71 gmx_wallcycle* wcycle,
73 VirtualSitesHandler* vsite,
74 ImdSession* imdSession,
76 std::vector<DomDecCallback>&& domdecCallbacks) :
78 isVerbose_(isVerbose),
79 verbosePrintInterval_(verbosePrintInterval),
80 nstglobalcomm_(nstglobalcomm),
81 domdecCallbacks_(std::move(domdecCallbacks)),
82 statePropagatorData_(statePropagatorData),
83 topologyHolder_(topologyHolder),
94 imdSession_(imdSession),
97 GMX_ASSERT(haveDDAtomOrdering(*cr),
98 "Domain decomposition Helper constructed in non-DD simulation");
101 void DomDecHelper::setup()
103 // constant choices for this call to dd_partition_system
104 const bool verbose = false;
105 const bool isMasterState = true;
106 const int nstglobalcomm = 1;
107 gmx_wallcycle* wcycle = nullptr;
109 // Distribute the charge groups over the nodes from the master node
110 partitionSystem(verbose,
114 statePropagatorData_->localState(),
115 statePropagatorData_->globalState());
118 void DomDecHelper::run(Step step, Time gmx_unused time)
120 if (step != nextNSStep_ || (step == inputrec_->init_step && inputrec_->bContinuation))
124 t_state* localState = statePropagatorData_->localState();
125 t_state* globalState = statePropagatorData_->globalState();
127 // constant choices for this call to dd_partition_system
128 const bool verbose = isVerbose_ && (step % verbosePrintInterval_ == 0 || step == inputrec_->init_step);
129 bool isMasterState = false;
131 // Correct the new box if it is too skewed
132 if (inputrecDynamicBox(inputrec_))
134 // TODO: Correcting the box is done here (if using DD) or in ForceElement (non-DD simulations).
135 // Think about unifying this responsibility, could this be done in one place?
136 if (correct_box(fplog_, step, localState->box))
138 isMasterState = true;
143 dd_collect_state(cr_->dd, localState, globalState);
146 // Distribute the charge groups over the nodes from the master node
147 partitionSystem(verbose, isMasterState, nstglobalcomm_, wcycle_, localState, globalState);
150 void DomDecHelper::partitionSystem(bool verbose,
153 gmx_wallcycle* wcycle,
155 t_state* globalState)
157 ForceBuffers* forcePointer = statePropagatorData_->forcePointer();
159 // Work-around to keep dd_partition_system from failing -
160 // we're not actually using the information related to Nose-Hoover chains
161 localState->nhchainlength = inputrec_->opts.nhchainlength;
162 // Distribute the charge groups over the nodes from the master node
163 dd_partition_system(fplog_,
165 inputrec_->init_step,
170 topologyHolder_->globalTopology(),
177 topologyHolder_->localTopology_,
184 statePropagatorData_->setLocalState(localState);
185 for (const auto& callback : domdecCallbacks_)
191 std::optional<SignallerCallback> DomDecHelper::registerNSCallback()
193 return [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; };
196 void DomDecHelperBuilder::registerClient(IDomDecHelperClient* client)
202 if (state_ == ModularSimulatorBuilderState::NotAcceptingClientRegistrations)
204 GMX_THROW(SimulationAlgorithmSetupError(
205 "Tried to register to DomDecHelper after it was built."));
207 clients_.emplace_back(client);