2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 "statepropagatordata.h"
54 #include "topologyholder.h"
58 DomDecHelper::DomDecHelper(
60 int verbosePrintInterval
,
61 StatePropagatorData
*statePropagatorData
,
62 TopologyHolder
*topologyHolder
,
63 CheckBondedInteractionsCallbackPtr checkBondedInteractionsCallback
,
67 const MDLogger
&mdlog
,
72 gmx_wallcycle
*wcycle
,
75 ImdSession
*imdSession
,
78 isVerbose_(isVerbose
),
79 verbosePrintInterval_(verbosePrintInterval
),
80 nstglobalcomm_(nstglobalcomm
),
81 statePropagatorData_(statePropagatorData
),
82 topologyHolder_(topologyHolder
),
83 checkBondedInteractionsCallback_(std::move(checkBondedInteractionsCallback
)),
94 imdSession_(imdSession
),
97 GMX_ASSERT(DOMAINDECOMP(cr
), "Domain decomposition Helper constructed in non-DD simulation");
100 void DomDecHelper::setup()
102 std::unique_ptr
<t_state
> localState
= statePropagatorData_
->localState();
103 t_state
*globalState
= statePropagatorData_
->globalState();
104 PaddedHostVector
<RVec
> *forcePointer
= statePropagatorData_
->forcePointer();
106 // constant choices for this call to dd_partition_system
107 const bool verbose
= false;
108 const bool isMasterState
= true;
109 const int nstglobalcomm
= 1;
110 gmx_wallcycle
*wcycle
= nullptr;
112 // Distribute the charge groups over the nodes from the master node
114 fplog_
, mdlog_
, inputrec_
->init_step
, cr_
, isMasterState
,
115 nstglobalcomm
, globalState
, topologyHolder_
->globalTopology(),
116 inputrec_
, imdSession_
, pull_work_
, localState
.get(), forcePointer
,
117 mdAtoms_
, topologyHolder_
->localTopology_
.get(), fr_
,
118 vsite_
, constr_
, nrnb_
, wcycle
, verbose
);
119 topologyHolder_
->updateLocalTopology();
120 (*checkBondedInteractionsCallback_
)();
121 statePropagatorData_
->setLocalState(std::move(localState
));
124 void DomDecHelper::run(Step step
, Time gmx_unused time
)
126 if (step
!= nextNSStep_
|| (step
== inputrec_
->init_step
&& inputrec_
->bContinuation
))
130 std::unique_ptr
<t_state
> localState
= statePropagatorData_
->localState();
131 t_state
*globalState
= statePropagatorData_
->globalState();
132 PaddedHostVector
<RVec
> *forcePointer
= statePropagatorData_
->forcePointer();
134 // constant choices for this call to dd_partition_system
137 (step
% verbosePrintInterval_
== 0 || step
== inputrec_
->init_step
);
138 bool isMasterState
= false;
140 // Correct the new box if it is too skewed
141 if (inputrecDynamicBox(inputrec_
))
143 t_graph
*graph
= nullptr;
144 if (correct_box(fplog_
, step
, localState
->box
, graph
))
146 isMasterState
= true;
151 dd_collect_state(cr_
->dd
, localState
.get(), globalState
);
154 // Distribute the charge groups over the nodes from the master node
156 fplog_
, mdlog_
, step
, cr_
, isMasterState
,
157 nstglobalcomm_
, globalState
, topologyHolder_
->globalTopology(),
158 inputrec_
, imdSession_
, pull_work_
,
159 localState
.get(), forcePointer
, mdAtoms_
,
160 topologyHolder_
->localTopology_
.get(), fr_
,
161 vsite_
, constr_
, nrnb_
, wcycle_
, verbose
);
162 topologyHolder_
->updateLocalTopology();
163 (*checkBondedInteractionsCallback_
)();
164 statePropagatorData_
->setLocalState(std::move(localState
));
167 SignallerCallbackPtr
DomDecHelper::registerNSCallback()
169 return std::make_unique
<SignallerCallback
>(
170 [this](Step step
, Time gmx_unused time
)
171 {this->nextNSStep_
= step
; });