Add index group option to density fitting
[alexxy/gromacs.git] / src / gromacs / applied_forces / densityfitting.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
36  * \brief
37  * Declares data structure and utilities for density fitting
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \ingroup module_applied_forces
41  */
42 #include "gmxpre.h"
43
44 #include "densityfitting.h"
45
46 #include "gromacs/mdrunutility/mdmodulenotification.h"
47 #include "gromacs/mdtypes/imdmodule.h"
48 #include "gromacs/utility/keyvaluetreebuilder.h"
49
50 #include "densityfittingforceprovider.h"
51 #include "densityfittingoptions.h"
52 #include "densityfittingoutputprovider.h"
53
54
55 namespace gmx
56 {
57
58 class IMdpOptionProvider;
59 class DensityFittingForceProvider;
60
61 namespace
62 {
63
64 /*! \internal
65  * \brief Density fitting
66  *
67  * Class that implements the density fitting forces and potential
68  * \note the virial calculation is not yet implemented
69  */
70 class DensityFitting final : public IMDModule
71 {
72     public:
73         /*! \brief Construct the density fitting module.
74          *
75          * \param[in] notifier allows the module to subscribe to notifications from MdModules.
76          *
77          * The density fitting code subscribes to these notifications:
78          *   - setting atom group indices in the densityFittingOptions_ by
79          *     taking a parmeter const IndexGroupsAndNames &
80          *   - storing its internal parameters in a tpr file by writing to a
81          *     key-value-tree during pre-processing by a function taking a
82          *     KeyValueTreeObjectBuilder as parameter
83          *   - reading its internal parameters from a key-value-tree during
84          *     simulation setup by taking a const KeyValueTreeObject & parameter
85          */
86         explicit DensityFitting(MdModulesNotifier *notifier)
87         {
88             // Callbacks for several kinds of MdModuleNotification are created
89             // and subscribed, and will be dispatched correctly at run time
90             // based on the type of the parameter required by the lambda.
91
92             // Setting atom group indices
93             const auto setFitGroupIndicesFunction = [this](const IndexGroupsAndNames &indexGroupsAndNames) {
94                     densityFittingOptions_.setFitGroupIndices(indexGroupsAndNames);
95                 };
96             notifier->notifier_.subscribe(setFitGroupIndicesFunction);
97
98             // Writing internal parameters during pre-processing
99             const auto writeInternalParametersFunction = [this](KeyValueTreeObjectBuilder treeBuilder) {
100                     densityFittingOptions_.writeInternalParametersToKvt(treeBuilder);
101                 };
102             notifier->notifier_.subscribe(writeInternalParametersFunction);
103
104             // Reading internal parameters during simulation setup
105             const auto readInternalParametersFunction = [this](const KeyValueTreeObject &tree) {
106                     densityFittingOptions_.readInternalParametersFromKvt(tree);
107                 };
108             notifier->notifier_.subscribe(readInternalParametersFunction);
109         }
110
111         //! From IMDModule; this class provides the mdpOptions itself
112         IMdpOptionProvider *mdpOptionProvider() override { return &densityFittingOptions_; }
113
114         //! Add this module to the force providers if active
115         void initForceProviders(ForceProviders *forceProviders) override
116         {
117             if (densityFittingOptions_.active())
118             {
119                 const auto &parameters = densityFittingOptions_.buildParameters();
120                 forceProvider_ = std::make_unique<DensityFittingForceProvider>(parameters);
121                 forceProviders->addForceProvider(forceProvider_.get());
122             }
123         }
124
125         //! This MDModule provides its own output
126         IMDOutputProvider *outputProvider() override { return &densityFittingOutputProvider_; }
127
128     private:
129         //! The output provider
130         DensityFittingOutputProvider                 densityFittingOutputProvider_;
131         //! The options provided for density fitting
132         DensityFittingOptions                        densityFittingOptions_;
133         //! Object that evaluates the forces
134         std::unique_ptr<DensityFittingForceProvider> forceProvider_;
135 };
136
137 }   // namespace
138
139 std::unique_ptr<IMDModule> DensityFittingModuleInfo::create(MdModulesNotifier * notifier)
140 {
141     return std::make_unique<DensityFitting>(notifier);
142 }
143
144 const std::string DensityFittingModuleInfo::name_ = "density-guided-simulation";
145
146 } // namespace gmx