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