2f1119c06eed5f746b685f10ada7aaa331fbdbcd
[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         std::unique_ptr<IMDModule>      densityFitting_;
97         std::unique_ptr<IMDModule>      field_;
98         std::unique_ptr<ForceProviders> forceProviders_;
99         std::unique_ptr<IMDModule>      imd_;
100         std::unique_ptr<IMDModule>      swapCoordinates_;
101
102         /*! \brief List of registered MDModules
103          *
104          * Note that MDModules::Impl owns this container, but it is only used by
105          * the MDModules::initForceProviders() function. To be consistent with
106          * IMDModule's vision, as indicated by its docs, we should
107          * \todo update IMDModule docs to allow nullptr return values
108          * \todo check for nullptr returned by IMDModule methods.
109          * \todo include field_ in modules_
110          */
111         std::vector< std::shared_ptr<IMDModule> > modules_;
112
113         //! Manages resources and notifies the MD modules when available
114         MdModulesNotifier notifier_;
115 };
116
117 MDModules::MDModules() : impl_(new Impl)
118 {
119 }
120
121 MDModules::~MDModules()
122 {
123 }
124
125 void MDModules::initMdpTransform(IKeyValueTreeTransformRules *rules)
126 {
127     auto appliedForcesScope = rules->scopedTransform("/applied-forces");
128     impl_->field_->mdpOptionProvider()->initMdpTransform(appliedForcesScope.rules());
129     impl_->densityFitting_->mdpOptionProvider()->initMdpTransform(appliedForcesScope.rules());
130 }
131
132 void MDModules::buildMdpOutput(KeyValueTreeObjectBuilder *builder)
133 {
134     impl_->field_->mdpOptionProvider()->buildMdpOutput(builder);
135     impl_->densityFitting_->mdpOptionProvider()->buildMdpOutput(builder);
136 }
137
138 void MDModules::assignOptionsToModules(const KeyValueTreeObject  &params,
139                                        IKeyValueTreeErrorHandler *errorHandler)
140 {
141     Options moduleOptions;
142     impl_->makeModuleOptions(&moduleOptions);
143     // The actual output is in the data fields of the modules that
144     // were set up in the module options.
145     assignOptionsFromKeyValueTree(&moduleOptions, params, errorHandler);
146 }
147
148 void MDModules::adjustInputrecBasedOnModules(t_inputrec *ir)
149 {
150     Options moduleOptions;
151     impl_->makeModuleOptions(&moduleOptions);
152
153     checkForUnknownOptionsInKeyValueTree(*ir->params, moduleOptions);
154
155     std::unique_ptr<KeyValueTreeObject> params(
156             new KeyValueTreeObject(
157                     adjustKeyValueTreeFromOptions(*ir->params, moduleOptions)));
158     delete ir->params;
159     ir->params = params.release();
160 }
161
162 IMDOutputProvider *MDModules::outputProvider()
163 {
164     return impl_.get();
165 }
166
167 ForceProviders *MDModules::initForceProviders()
168 {
169     GMX_RELEASE_ASSERT(impl_->forceProviders_ == nullptr,
170                        "Force providers initialized multiple times");
171     impl_->forceProviders_ = std::make_unique<ForceProviders>();
172     impl_->field_->initForceProviders(impl_->forceProviders_.get());
173     impl_->densityFitting_->initForceProviders(impl_->forceProviders_.get());
174     for (auto && module : impl_->modules_)
175     {
176         module->initForceProviders(impl_->forceProviders_.get());
177     }
178     return impl_->forceProviders_.get();
179 }
180
181 void MDModules::add(std::shared_ptr<gmx::IMDModule> module)
182 {
183     impl_->modules_.emplace_back(std::move(module));
184 }
185
186 const MdModulesNotifier &MDModules::notifier()
187 {
188     return impl_->notifier_;
189 }
190
191 } // namespace gmx