Use more forward declarations to removed header dependencies
[alexxy/gromacs.git] / src / gromacs / modularsimulator / pmeloadbalancehelper.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 Declares the PME load balancing 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 "pmeloadbalancehelper.h"
45
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"
54
55 #include "statepropagatordata.h"
56
57 namespace gmx
58 {
59 bool PmeLoadBalanceHelper::doPmeLoadBalancing(const MdrunOptions& mdrunOptions,
60                                               const t_inputrec*   inputrec,
61                                               const t_forcerec*   fr)
62 {
63     return (mdrunOptions.tunePme && EEL_PME(fr->ic->eeltype) && !mdrunOptions.reproducible
64             && inputrec->cutoff_scheme != ecutsGROUP);
65 }
66
67 PmeLoadBalanceHelper::PmeLoadBalanceHelper(bool                 isVerbose,
68                                            StatePropagatorData* statePropagatorData,
69                                            FILE*                fplog,
70                                            t_commrec*           cr,
71                                            const MDLogger&      mdlog,
72                                            const t_inputrec*    inputrec,
73                                            gmx_wallcycle*       wcycle,
74                                            t_forcerec*          fr) :
75     pme_loadbal_(nullptr),
76     nextNSStep_(-1),
77     isVerbose_(isVerbose),
78     bPMETunePrinting_(false),
79     statePropagatorData_(statePropagatorData),
80     fplog_(fplog),
81     cr_(cr),
82     mdlog_(mdlog),
83     inputrec_(inputrec),
84     wcycle_(wcycle),
85     fr_(fr)
86 {
87 }
88
89 void PmeLoadBalanceHelper::setup()
90 {
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.");
94     pme_loadbal_init(&pme_loadbal_, cr_, mdlog_, *inputrec_, box, *fr_->ic, *fr_->nbv, fr_->pmedata,
95                      fr_->nbv->useGpu());
96 }
97
98 void PmeLoadBalanceHelper::run(gmx::Step step, gmx::Time gmx_unused time)
99 {
100     if (step != nextNSStep_ || step == inputrec_->init_step)
101     {
102         return;
103     }
104
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_, cr_, (isVerbose_ && MASTER(cr_)) ? stderr : nullptr, fplog_,
109                    mdlog_, *inputrec_, fr_, statePropagatorData_->constBox(),
110                    statePropagatorData_->constPositionsView().paddedArrayRef(), wcycle_, step,
111                    step - inputrec_->init_step, &bPMETunePrinting_, false);
112 }
113
114 void PmeLoadBalanceHelper::teardown()
115 {
116     pme_loadbal_done(pme_loadbal_, fplog_, mdlog_, fr_->nbv->useGpu());
117 }
118
119 bool PmeLoadBalanceHelper::pmePrinting()
120 {
121     return bPMETunePrinting_;
122 }
123
124 const pme_load_balancing_t* PmeLoadBalanceHelper::loadBalancingObject()
125 {
126     return pme_loadbal_;
127 }
128
129 SignallerCallbackPtr PmeLoadBalanceHelper::registerNSCallback()
130 {
131     return std::make_unique<SignallerCallback>(
132             [this](Step step, Time gmx_unused time) { nextNSStep_ = step; });
133 }
134
135 } // namespace gmx