1014a0adde1933fa711ff80d20f16746abd7d021
[alexxy/gromacs.git] / src / gromacs / mdrun / mdmodules.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,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 #include "gmxpre.h"
36
37 #include "mdmodules.h"
38
39 #include <memory>
40
41 #include "gromacs/applied_forces/electricfield.h"
42 #include "gromacs/applied_forces/densityfitting/densityfitting.h"
43 #include "gromacs/imd/imd.h"
44 #include "gromacs/mdtypes/iforceprovider.h"
45 #include "gromacs/mdtypes/imdmodule.h"
46 #include "gromacs/mdtypes/imdoutputprovider.h"
47 #include "gromacs/mdtypes/imdpoptionprovider.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/options/options.h"
50 #include "gromacs/options/optionsection.h"
51 #include "gromacs/options/treesupport.h"
52 #include "gromacs/swap/swapcoords.h"
53 #include "gromacs/utility/keyvaluetree.h"
54 #include "gromacs/utility/keyvaluetreebuilder.h"
55 #include "gromacs/utility/keyvaluetreetransform.h"
56 #include "gromacs/utility/mdmodulenotification.h"
57 #include "gromacs/utility/smalloc.h"
58
59 namespace gmx
60 {
61
62 class MDModules::Impl : public IMDOutputProvider
63 {
64 public:
65     Impl() :
66         densityFitting_(DensityFittingModuleInfo::create()),
67         field_(createElectricFieldModule()),
68         imd_(createInteractiveMolecularDynamicsModule()),
69         swapCoordinates_(createSwapCoordinatesModule())
70     {
71     }
72
73     void makeModuleOptions(Options* options)
74     {
75         // Create a section for applied-forces modules
76         auto appliedForcesOptions = options->addSection(OptionSection("applied-forces"));
77         field_->mdpOptionProvider()->initMdpOptions(&appliedForcesOptions);
78         densityFitting_->mdpOptionProvider()->initMdpOptions(&appliedForcesOptions);
79         // In future, other sections would also go here.
80     }
81
82     // From IMDOutputProvider
83     void initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv) override
84     {
85         field_->outputProvider()->initOutput(fplog, nfile, fnm, bAppendFiles, oenv);
86         densityFitting_->outputProvider()->initOutput(fplog, nfile, fnm, bAppendFiles, oenv);
87     }
88     void finishOutput() override
89     {
90         field_->outputProvider()->finishOutput();
91         densityFitting_->outputProvider()->finishOutput();
92     }
93
94     /*! \brief Manages callbacks and notifies the MD modules.
95      *
96      * \note The notifier must be constructed before the modules and shall
97      *       not be destructed before the modules are destructed.
98      */
99     MdModulesNotifier notifier_;
100
101     std::unique_ptr<IMDModule>      densityFitting_;
102     std::unique_ptr<IMDModule>      field_;
103     std::unique_ptr<ForceProviders> forceProviders_;
104     std::unique_ptr<IMDModule>      imd_;
105     std::unique_ptr<IMDModule>      swapCoordinates_;
106
107     /*! \brief List of registered MDModules
108      *
109      * Note that MDModules::Impl owns this container, but it is only used by
110      * the MDModules::initForceProviders() function. To be consistent with
111      * IMDModule's vision, as indicated by its docs, we should
112      * \todo update IMDModule docs to allow nullptr return values
113      * \todo check for nullptr returned by IMDModule methods.
114      * \todo include field_ in modules_
115      */
116     std::vector<std::shared_ptr<IMDModule>> modules_;
117 };
118
119 MDModules::MDModules() : impl_(new Impl) {}
120
121 MDModules::~MDModules() {}
122
123 void MDModules::initMdpTransform(IKeyValueTreeTransformRules* rules)
124 {
125     auto appliedForcesScope = rules->scopedTransform("/applied-forces");
126     impl_->field_->mdpOptionProvider()->initMdpTransform(appliedForcesScope.rules());
127     impl_->densityFitting_->mdpOptionProvider()->initMdpTransform(appliedForcesScope.rules());
128 }
129
130 void MDModules::buildMdpOutput(KeyValueTreeObjectBuilder* builder)
131 {
132     impl_->field_->mdpOptionProvider()->buildMdpOutput(builder);
133     impl_->densityFitting_->mdpOptionProvider()->buildMdpOutput(builder);
134 }
135
136 void MDModules::assignOptionsToModules(const KeyValueTreeObject& params, IKeyValueTreeErrorHandler* errorHandler)
137 {
138     Options moduleOptions;
139     impl_->makeModuleOptions(&moduleOptions);
140     // The actual output is in the data fields of the modules that
141     // were set up in the module options.
142     assignOptionsFromKeyValueTree(&moduleOptions, params, errorHandler);
143 }
144
145 void MDModules::adjustInputrecBasedOnModules(t_inputrec* ir)
146 {
147     Options moduleOptions;
148     impl_->makeModuleOptions(&moduleOptions);
149
150     checkForUnknownOptionsInKeyValueTree(*ir->params, moduleOptions);
151
152     std::unique_ptr<KeyValueTreeObject> params(
153             new KeyValueTreeObject(adjustKeyValueTreeFromOptions(*ir->params, moduleOptions)));
154     delete ir->params;
155     ir->params = params.release();
156 }
157
158 IMDOutputProvider* MDModules::outputProvider()
159 {
160     return impl_.get();
161 }
162
163 ForceProviders* MDModules::initForceProviders()
164 {
165     GMX_RELEASE_ASSERT(impl_->forceProviders_ == nullptr,
166                        "Force providers initialized multiple times");
167     impl_->forceProviders_ = std::make_unique<ForceProviders>();
168     impl_->field_->initForceProviders(impl_->forceProviders_.get());
169     impl_->densityFitting_->initForceProviders(impl_->forceProviders_.get());
170     for (auto&& module : impl_->modules_)
171     {
172         module->initForceProviders(impl_->forceProviders_.get());
173     }
174     return impl_->forceProviders_.get();
175 }
176
177 void MDModules::subscribeToPreProcessingNotifications()
178 {
179     impl_->densityFitting_->subscribeToPreProcessingNotifications(&impl_->notifier_);
180 }
181
182 void MDModules::subscribeToSimulationSetupNotifications()
183 {
184     impl_->densityFitting_->subscribeToSimulationSetupNotifications(&impl_->notifier_);
185 }
186
187 void MDModules::add(std::shared_ptr<gmx::IMDModule> module)
188 {
189     impl_->modules_.emplace_back(std::move(module));
190 }
191
192 const MdModulesNotifier& MDModules::notifier()
193 {
194     return impl_->notifier_;
195 }
196
197 } // namespace gmx