Move responsibility for checking bondeds in reverse topology
[alexxy/gromacs.git] / src / gromacs / modularsimulator / domdechelper.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  * \brief Defines the domain decomposition helper for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  */
41
42 #include "gmxpre.h"
43
44 #include "domdechelper.h"
45
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"
52
53 #include "freeenergyperturbationdata.h"
54 #include "statepropagatordata.h"
55 #include "topologyholder.h"
56
57 namespace gmx
58 {
59 DomDecHelper::DomDecHelper(bool                        isVerbose,
60                            int                         verbosePrintInterval,
61                            StatePropagatorData*        statePropagatorData,
62                            FreeEnergyPerturbationData* freeEnergyPerturbationData,
63                            TopologyHolder*             topologyHolder,
64                            int                         nstglobalcomm,
65                            FILE*                       fplog,
66                            t_commrec*                  cr,
67                            const MDLogger&             mdlog,
68                            Constraints*                constr,
69                            const t_inputrec*           inputrec,
70                            MDAtoms*                    mdAtoms,
71                            t_nrnb*                     nrnb,
72                            gmx_wallcycle*              wcycle,
73                            t_forcerec*                 fr,
74                            VirtualSitesHandler*        vsite,
75                            ImdSession*                 imdSession,
76                            pull_t*                     pull_work) :
77     nextNSStep_(-1),
78     isVerbose_(isVerbose),
79     verbosePrintInterval_(verbosePrintInterval),
80     nstglobalcomm_(nstglobalcomm),
81     statePropagatorData_(statePropagatorData),
82     freeEnergyPerturbationData_(freeEnergyPerturbationData),
83     topologyHolder_(topologyHolder),
84     fplog_(fplog),
85     cr_(cr),
86     mdlog_(mdlog),
87     constr_(constr),
88     inputrec_(inputrec),
89     mdAtoms_(mdAtoms),
90     nrnb_(nrnb),
91     wcycle_(wcycle),
92     fr_(fr),
93     vsite_(vsite),
94     imdSession_(imdSession),
95     pull_work_(pull_work)
96 {
97     GMX_ASSERT(DOMAINDECOMP(cr), "Domain decomposition Helper constructed in non-DD simulation");
98 }
99
100 void DomDecHelper::setup()
101 {
102     // constant choices for this call to dd_partition_system
103     const bool     verbose       = false;
104     const bool     isMasterState = true;
105     const int      nstglobalcomm = 1;
106     gmx_wallcycle* wcycle        = nullptr;
107
108     // Distribute the charge groups over the nodes from the master node
109     partitionSystem(verbose,
110                     isMasterState,
111                     nstglobalcomm,
112                     wcycle,
113                     statePropagatorData_->localState(),
114                     statePropagatorData_->globalState());
115 }
116
117 void DomDecHelper::run(Step step, Time gmx_unused time)
118 {
119     if (step != nextNSStep_ || (step == inputrec_->init_step && inputrec_->bContinuation))
120     {
121         return;
122     }
123     std::unique_ptr<t_state> localState  = statePropagatorData_->localState();
124     t_state*                 globalState = statePropagatorData_->globalState();
125
126     // constant choices for this call to dd_partition_system
127     const bool verbose = isVerbose_ && (step % verbosePrintInterval_ == 0 || step == inputrec_->init_step);
128     bool isMasterState = false;
129
130     // Correct the new box if it is too skewed
131     if (inputrecDynamicBox(inputrec_))
132     {
133         // TODO: Correcting the box is done here (if using DD) or in ForceElement (non-DD simulations).
134         //       Think about unifying this responsibility, could this be done in one place?
135         if (correct_box(fplog_, step, localState->box))
136         {
137             isMasterState = true;
138         }
139     }
140     if (isMasterState)
141     {
142         dd_collect_state(cr_->dd, localState.get(), globalState);
143     }
144
145     // Distribute the charge groups over the nodes from the master node
146     partitionSystem(verbose, isMasterState, nstglobalcomm_, wcycle_, std::move(localState), globalState);
147 }
148
149 void DomDecHelper::partitionSystem(bool                     verbose,
150                                    bool                     isMasterState,
151                                    int                      nstglobalcomm,
152                                    gmx_wallcycle*           wcycle,
153                                    std::unique_ptr<t_state> localState,
154                                    t_state*                 globalState)
155 {
156     ForceBuffers* forcePointer = statePropagatorData_->forcePointer();
157
158     // Distribute the charge groups over the nodes from the master node
159     dd_partition_system(fplog_,
160                         mdlog_,
161                         inputrec_->init_step,
162                         cr_,
163                         isMasterState,
164                         nstglobalcomm,
165                         globalState,
166                         topologyHolder_->globalTopology(),
167                         *inputrec_,
168                         imdSession_,
169                         pull_work_,
170                         localState.get(),
171                         forcePointer,
172                         mdAtoms_,
173                         topologyHolder_->localTopology_.get(),
174                         fr_,
175                         vsite_,
176                         constr_,
177                         nrnb_,
178                         wcycle,
179                         verbose);
180     topologyHolder_->updateLocalTopology();
181     statePropagatorData_->setLocalState(std::move(localState));
182     if (freeEnergyPerturbationData_)
183     {
184         freeEnergyPerturbationData_->updateMDAtoms();
185     }
186 }
187
188 std::optional<SignallerCallback> DomDecHelper::registerNSCallback()
189 {
190     return [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; };
191 }
192
193 } // namespace gmx