Merge branch release-2021
[alexxy/gromacs.git] / src / gromacs / applied_forces / electricfield.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,2019, The GROMACS development team.
5  * Copyright (c) 2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * Defines data structure and utilities for electric fields
39  *
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  * \ingroup module_applied_forces
42  */
43 #include "gmxpre.h"
44
45 #include "electricfield.h"
46
47 #include <cmath>
48
49 #include <memory>
50 #include <string>
51
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/forceoutput.h"
58 #include "gromacs/mdtypes/iforceprovider.h"
59 #include "gromacs/mdtypes/imdmodule.h"
60 #include "gromacs/mdtypes/imdoutputprovider.h"
61 #include "gromacs/mdtypes/imdpoptionprovider.h"
62 #include "gromacs/mdtypes/mdatom.h"
63 #include "gromacs/options/basicoptions.h"
64 #include "gromacs/options/ioptionscontainerwithsections.h"
65 #include "gromacs/options/optionsection.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/keyvaluetreebuilder.h"
68 #include "gromacs/utility/keyvaluetreetransform.h"
69 #include "gromacs/utility/pleasecite.h"
70 #include "gromacs/utility/strconvert.h"
71
72 namespace gmx
73 {
74
75 namespace
76 {
77
78 /*! \internal
79  * \brief Describes an applied electric field in a coordinate
80  * dimension.
81  *
82  * Can compute the applied field strength at a time, and supports
83  * operations to read and form corresponding .mdp contents.
84  */
85 class ElectricFieldDimension
86 {
87 public:
88     /*! \brief
89      * Adds an option section to specify parameters for this field component.
90      */
91     void initMdpOptions(IOptionsContainerWithSections* options, const char* sectionName)
92     {
93         auto section = options->addSection(OptionSection(sectionName));
94         section.addOption(RealOption("E0").store(&a_));
95         section.addOption(RealOption("omega").store(&omega_));
96         section.addOption(RealOption("t0").store(&t0_));
97         section.addOption(RealOption("sigma").store(&sigma_));
98     }
99     /*! \brief
100      * Creates mdp parameters for this field component.
101      */
102     void buildMdpOutput(KeyValueTreeObjectBuilder* builder, const std::string& name) const
103     {
104         builder->addUniformArray<real>("electric-field-" + name, { a_, omega_, t0_, sigma_ });
105     }
106
107     /*! \brief Evaluates this field component at given time.
108      *
109      * \param[in] t The time to evualate at
110      * \return The electric field
111      */
112     real evaluate(real t) const
113     {
114         if (sigma_ > 0)
115         {
116             return a_ * (std::cos(omega_ * (t - t0_)) * std::exp(-square(t - t0_) / (2.0 * square(sigma_))));
117         }
118         else
119         {
120             return a_ * std::cos(omega_ * t);
121         }
122     }
123
124     //! Return the amplitude
125     real a() const { return a_; }
126
127 private:
128     //! Coeffient (V / nm)
129     real a_ = 0;
130     //! Frequency (1/ps)
131     real omega_ = 0;
132     //! Central time point (ps) for pulse
133     real t0_ = 0;
134     //! Width of pulse (ps, if zero there is no pulse)
135     real sigma_ = 0;
136 };
137
138 /*! \internal
139  * \brief Describe time dependent electric field
140  *
141  * Class that implements a force to be evaluated in mdrun.
142  * The electric field can be pulsed and oscillating, simply
143  * oscillating, or static, in each of X,Y,Z directions.
144  */
145 class ElectricField final : public IMDModule, public IMdpOptionProvider, public IMDOutputProvider, public IForceProvider
146 {
147 public:
148     ElectricField() : fpField_(nullptr) {}
149
150     // From IMDModule
151     IMdpOptionProvider* mdpOptionProvider() override { return this; }
152     IMDOutputProvider*  outputProvider() override { return this; }
153     void                initForceProviders(ForceProviders* forceProviders) override
154     {
155         if (isActive())
156         {
157             forceProviders->addForceProvider(this);
158         }
159     }
160
161     // From IMdpOptionProvider
162     void initMdpTransform(IKeyValueTreeTransformRules* transform) override;
163     void initMdpOptions(IOptionsContainerWithSections* options) override;
164     void buildMdpOutput(KeyValueTreeObjectBuilder* builder) const override;
165
166     // From IMDOutputProvider
167     void initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv) override;
168     void finishOutput() override;
169
170     // From IForceProvider
171     //! \copydoc IForceProvider::calculateForces()
172     void calculateForces(const ForceProviderInput& forceProviderInput,
173                          ForceProviderOutput*      forceProviderOutput) override;
174
175     void subscribeToSimulationSetupNotifications(MdModulesNotifier* /* notifier */) override {}
176     void subscribeToPreProcessingNotifications(MdModulesNotifier* /* notifier */) override {}
177
178 private:
179     //! Return whether or not to apply a field
180     bool isActive() const;
181
182     /*! \brief Return the field strength
183      *
184      * \param[in] dim The spatial direction
185      * \param[in] t   The time (ps)
186      * \return The field strength in V/nm units
187      */
188     real field(int dim, real t) const;
189
190     /*! \brief Print the field components to a file
191      *
192      * \param[in] t   The time
193      * Will throw and exit with fatal error if file is not open.
194      */
195     void printComponents(double t) const;
196
197     //! The components of the applied electric field in each coordinate dimension
198     ElectricFieldDimension efield_[DIM];
199     //! File pointer for electric field
200     FILE* fpField_;
201 };
202
203 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
204 void convertParameters(gmx::KeyValueTreeObjectBuilder* builder, const std::string& value)
205 {
206     const std::vector<std::string> sxt = splitString(value);
207     if (sxt.empty())
208     {
209         return;
210     }
211     if (sxt.size() != 4)
212     {
213         GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
214     }
215     builder->addValue<real>("E0", fromString<real>(sxt[0]));
216     builder->addValue<real>("omega", fromString<real>(sxt[1]));
217     builder->addValue<real>("t0", fromString<real>(sxt[2]));
218     builder->addValue<real>("sigma", fromString<real>(sxt[3]));
219 }
220
221 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules* rules)
222 {
223     rules->addRule().from<std::string>("/electric-field-x").toObject("/electric-field/x").transformWith(&convertParameters);
224     rules->addRule().from<std::string>("/electric-field-y").toObject("/electric-field/y").transformWith(&convertParameters);
225     rules->addRule().from<std::string>("/electric-field-z").toObject("/electric-field/z").transformWith(&convertParameters);
226 }
227
228 void ElectricField::initMdpOptions(IOptionsContainerWithSections* options)
229 {
230     auto section = options->addSection(OptionSection("electric-field"));
231     efield_[XX].initMdpOptions(&section, "x");
232     efield_[YY].initMdpOptions(&section, "y");
233     efield_[ZZ].initMdpOptions(&section, "z");
234 }
235
236 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
237 {
238     std::string comment =
239             R"(; Electric fields
240 ; Format for electric-field-x, etc. is: four real variables:
241 ; amplitude (V/nm), frequency omega (1/ps), time for the pulse peak (ps),
242 ; and sigma (ps) width of the pulse. Omega = 0 means static field,
243 ; sigma = 0 means no pulse, leaving the field to be a cosine function.)";
244     builder->addValue<std::string>("comment-electric-field", comment);
245     efield_[XX].buildMdpOutput(builder, "x");
246     efield_[YY].buildMdpOutput(builder, "y");
247     efield_[ZZ].buildMdpOutput(builder, "z");
248 }
249
250 void ElectricField::initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv)
251 {
252     if (isActive())
253     {
254         please_cite(fplog, "Caleman2008a");
255
256         // Optional output file showing the field, see manual.
257         if (opt2bSet("-field", nfile, fnm))
258         {
259             if (bAppendFiles)
260             {
261                 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
262             }
263             else
264             {
265                 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
266                                     "Applied electric field",
267                                     "Time (ps)",
268                                     "E (V/nm)",
269                                     oenv);
270             }
271         }
272     }
273 }
274
275 void ElectricField::finishOutput()
276 {
277     if (fpField_ != nullptr)
278     {
279         // This is opened sometimes with xvgropen, sometimes with
280         // gmx_fio_fopen, so we use the least common denominator for closing.
281         gmx_fio_fclose(fpField_);
282         fpField_ = nullptr;
283     }
284 }
285
286 real ElectricField::field(int dim, real t) const
287 {
288     return efield_[dim].evaluate(t);
289 }
290
291 bool ElectricField::isActive() const
292 {
293     return (efield_[XX].a() != 0 || efield_[YY].a() != 0 || efield_[ZZ].a() != 0);
294 }
295
296 void ElectricField::printComponents(double t) const
297 {
298     fprintf(fpField_, "%10g  %10g  %10g  %10g\n", t, field(XX, t), field(YY, t), field(ZZ, t));
299 }
300
301 void ElectricField::calculateForces(const ForceProviderInput& forceProviderInput,
302                                     ForceProviderOutput*      forceProviderOutput)
303 {
304     if (isActive())
305     {
306         const t_mdatoms& mdatoms = forceProviderInput.mdatoms_;
307         const double     t       = forceProviderInput.t_;
308         const t_commrec& cr      = forceProviderInput.cr_;
309
310         // NOTE: The non-conservative electric field does not have a virial
311         ArrayRef<RVec> f = forceProviderOutput->forceWithVirial_.force_;
312
313         for (int m = 0; (m < DIM); m++)
314         {
315             const real fieldStrength = FIELDFAC * field(m, t);
316
317             if (fieldStrength != 0)
318             {
319                 // TODO: Check parallellism
320                 for (int i = 0; i < mdatoms.homenr; ++i)
321                 {
322                     // NOTE: Not correct with perturbed charges
323                     f[i][m] += mdatoms.chargeA[i] * fieldStrength;
324                 }
325             }
326         }
327         if (MASTER(&cr) && fpField_ != nullptr)
328         {
329             printComponents(t);
330         }
331     }
332 }
333
334 } // namespace
335
336 std::unique_ptr<IMDModule> createElectricFieldModule()
337 {
338     return std::make_unique<ElectricField>();
339 }
340
341 } // namespace gmx