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/partition.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/state.h"
51 #include "statepropagatordata.h"
52 #include "topologyholder.h"
56 DomDecHelper::DomDecHelper(
58 int verbosePrintInterval
,
59 StatePropagatorData
*statePropagatorData
,
60 TopologyHolder
*topologyHolder
,
61 CheckBondedInteractionsCallbackPtr checkBondedInteractionsCallback
,
65 const MDLogger
&mdlog
,
70 gmx_wallcycle
*wcycle
,
73 ImdSession
*imdSession
,
76 isVerbose_(isVerbose
),
77 verbosePrintInterval_(verbosePrintInterval
),
78 nstglobalcomm_(nstglobalcomm
),
79 statePropagatorData_(statePropagatorData
),
80 topologyHolder_(topologyHolder
),
81 checkBondedInteractionsCallback_(std::move(checkBondedInteractionsCallback
)),
92 imdSession_(imdSession
),
95 GMX_ASSERT(DOMAINDECOMP(cr
), "Domain decomposition Helper constructed in non-DD simulation");
98 void DomDecHelper::setup()
100 std::unique_ptr
<t_state
> localState
= statePropagatorData_
->localState();
101 t_state
*globalState
= statePropagatorData_
->globalState();
102 PaddedVector
<RVec
> *forcePointer
= statePropagatorData_
->forcePointer();
104 // constant choices for this call to dd_partition_system
105 const bool verbose
= false;
106 const bool isMasterState
= true;
107 const int nstglobalcomm
= 1;
108 gmx_wallcycle
*wcycle
= nullptr;
110 // Distribute the charge groups over the nodes from the master node
112 fplog_
, mdlog_
, inputrec_
->init_step
, cr_
, isMasterState
,
113 nstglobalcomm
, globalState
, topologyHolder_
->globalTopology(),
114 inputrec_
, imdSession_
, pull_work_
, localState
.get(), forcePointer
,
115 mdAtoms_
, topologyHolder_
->localTopology_
.get(), fr_
,
116 vsite_
, constr_
, nrnb_
, wcycle
, verbose
);
117 topologyHolder_
->updateLocalTopology();
118 (*checkBondedInteractionsCallback_
)();
119 statePropagatorData_
->setLocalState(std::move(localState
));
122 void DomDecHelper::run(Step step
, Time gmx_unused time
)
124 if (step
!= nextNSStep_
|| (step
== inputrec_
->init_step
&& inputrec_
->bContinuation
))
128 std::unique_ptr
<t_state
> localState
= statePropagatorData_
->localState();
129 t_state
*globalState
= statePropagatorData_
->globalState();
130 PaddedVector
<RVec
> *forcePointer
= statePropagatorData_
->forcePointer();
132 // constant choices for this call to dd_partition_system
135 (step
% verbosePrintInterval_
== 0 || step
== inputrec_
->init_step
);
136 const bool isMasterState
= false;
138 // Distribute the charge groups over the nodes from the master node
140 fplog_
, mdlog_
, inputrec_
->init_step
, cr_
, isMasterState
,
141 nstglobalcomm_
, globalState
, topologyHolder_
->globalTopology(),
142 inputrec_
, imdSession_
, pull_work_
,
143 localState
.get(), forcePointer
, mdAtoms_
,
144 topologyHolder_
->localTopology_
.get(), fr_
,
145 vsite_
, constr_
, nrnb_
, wcycle_
, verbose
);
146 topologyHolder_
->updateLocalTopology();
147 (*checkBondedInteractionsCallback_
)();
148 statePropagatorData_
->setLocalState(std::move(localState
));
151 SignallerCallbackPtr
DomDecHelper::registerNSCallback()
153 return std::make_unique
<SignallerCallback
>(
154 [this](Step step
, Time gmx_unused time
)
155 {this->nextNSStep_
= step
; });