Merge branch release-2021 into merge-2021-into-master
[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/options/basicoptions.h"
63 #include "gromacs/options/ioptionscontainerwithsections.h"
64 #include "gromacs/options/optionsection.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/keyvaluetreebuilder.h"
67 #include "gromacs/utility/keyvaluetreetransform.h"
68 #include "gromacs/utility/pleasecite.h"
69 #include "gromacs/utility/strconvert.h"
70
71 namespace gmx
72 {
73
74 namespace
75 {
76
77 /*! \internal
78  * \brief Describes an applied electric field in a coordinate
79  * dimension.
80  *
81  * Can compute the applied field strength at a time, and supports
82  * operations to read and form corresponding .mdp contents.
83  */
84 class ElectricFieldDimension
85 {
86 public:
87     /*! \brief
88      * Adds an option section to specify parameters for this field component.
89      */
90     void initMdpOptions(IOptionsContainerWithSections* options, const char* sectionName)
91     {
92         auto section = options->addSection(OptionSection(sectionName));
93         section.addOption(RealOption("E0").store(&a_));
94         section.addOption(RealOption("omega").store(&omega_));
95         section.addOption(RealOption("t0").store(&t0_));
96         section.addOption(RealOption("sigma").store(&sigma_));
97         if (sigma_ <= 0 && t0_ != 0.0)
98         {
99             GMX_THROW(
100                     InvalidInputError("Non-pulsed field (sigma = 0) ignores the value of t0. "
101                                       "Please, set t0 to 0 to avoid this error."));
102         }
103     }
104     /*! \brief
105      * Creates mdp parameters for this field component.
106      */
107     void buildMdpOutput(KeyValueTreeObjectBuilder* builder, const std::string& name) const
108     {
109         builder->addUniformArray<real>("electric-field-" + name, { a_, omega_, t0_, sigma_ });
110     }
111
112     /*! \brief Evaluates this field component at given time.
113      *
114      * \param[in] t The time to evualate at
115      * \return The electric field
116      */
117     real evaluate(real t) const
118     {
119         if (sigma_ > 0)
120         {
121             return a_ * (std::cos(omega_ * (t - t0_)) * std::exp(-square(t - t0_) / (2.0 * square(sigma_))));
122         }
123         else
124         {
125             return a_ * std::cos(omega_ * t);
126         }
127     }
128
129     //! Return the amplitude
130     real a() const { return a_; }
131
132 private:
133     //! Coeffient (V / nm)
134     real a_ = 0;
135     //! Frequency (1/ps)
136     real omega_ = 0;
137     //! Central time point (ps) for pulse
138     real t0_ = 0;
139     //! Width of pulse (ps, if zero there is no pulse)
140     real sigma_ = 0;
141 };
142
143 /*! \internal
144  * \brief Describe time dependent electric field
145  *
146  * Class that implements a force to be evaluated in mdrun.
147  * The electric field can be pulsed and oscillating, simply
148  * oscillating, or static, in each of X,Y,Z directions.
149  */
150 class ElectricField final : public IMDModule, public IMdpOptionProvider, public IMDOutputProvider, public IForceProvider
151 {
152 public:
153     ElectricField() : fpField_(nullptr) {}
154
155     // From IMDModule
156     IMdpOptionProvider* mdpOptionProvider() override { return this; }
157     IMDOutputProvider*  outputProvider() override { return this; }
158     void                initForceProviders(ForceProviders* forceProviders) override
159     {
160         if (isActive())
161         {
162             forceProviders->addForceProvider(this);
163         }
164     }
165
166     // From IMdpOptionProvider
167     void initMdpTransform(IKeyValueTreeTransformRules* transform) override;
168     void initMdpOptions(IOptionsContainerWithSections* options) override;
169     void buildMdpOutput(KeyValueTreeObjectBuilder* builder) const override;
170
171     // From IMDOutputProvider
172     void initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv) override;
173     void finishOutput() override;
174
175     // From IForceProvider
176     //! \copydoc IForceProvider::calculateForces()
177     void calculateForces(const ForceProviderInput& forceProviderInput,
178                          ForceProviderOutput*      forceProviderOutput) override;
179
180     void subscribeToSimulationSetupNotifications(MDModulesNotifiers* /* notifiers */) override {}
181     void subscribeToPreProcessingNotifications(MDModulesNotifiers* /* notifiers */) override {}
182
183 private:
184     //! Return whether or not to apply a field
185     bool isActive() const;
186
187     /*! \brief Return the field strength
188      *
189      * \param[in] dim The spatial direction
190      * \param[in] t   The time (ps)
191      * \return The field strength in V/nm units
192      */
193     real field(int dim, real t) const;
194
195     /*! \brief Print the field components to a file
196      *
197      * \param[in] t   The time
198      * Will throw and exit with fatal error if file is not open.
199      */
200     void printComponents(double t) const;
201
202     //! The components of the applied electric field in each coordinate dimension
203     ElectricFieldDimension efield_[DIM];
204     //! File pointer for electric field
205     FILE* fpField_;
206 };
207
208 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
209 void convertParameters(gmx::KeyValueTreeObjectBuilder* builder, const std::string& value)
210 {
211     const std::vector<std::string> sxt = splitString(value);
212     if (sxt.empty())
213     {
214         return;
215     }
216     if (sxt.size() != 4)
217     {
218         GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
219     }
220     builder->addValue<real>("E0", fromString<real>(sxt[0]));
221     builder->addValue<real>("omega", fromString<real>(sxt[1]));
222     builder->addValue<real>("t0", fromString<real>(sxt[2]));
223     builder->addValue<real>("sigma", fromString<real>(sxt[3]));
224 }
225
226 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules* rules)
227 {
228     rules->addRule().from<std::string>("/electric-field-x").toObject("/electric-field/x").transformWith(&convertParameters);
229     rules->addRule().from<std::string>("/electric-field-y").toObject("/electric-field/y").transformWith(&convertParameters);
230     rules->addRule().from<std::string>("/electric-field-z").toObject("/electric-field/z").transformWith(&convertParameters);
231 }
232
233 void ElectricField::initMdpOptions(IOptionsContainerWithSections* options)
234 {
235     auto section = options->addSection(OptionSection("electric-field"));
236     efield_[XX].initMdpOptions(&section, "x");
237     efield_[YY].initMdpOptions(&section, "y");
238     efield_[ZZ].initMdpOptions(&section, "z");
239 }
240
241 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
242 {
243     std::string comment =
244             R"(; Electric fields
245 ; Format for electric-field-x, etc. is: four real variables:
246 ; amplitude (V/nm), frequency omega (1/ps), time for the pulse peak (ps),
247 ; and sigma (ps) width of the pulse. Omega = 0 means static field,
248 ; sigma = 0 means no pulse, leaving the field to be a cosine function.)";
249     builder->addValue<std::string>("comment-electric-field", comment);
250     efield_[XX].buildMdpOutput(builder, "x");
251     efield_[YY].buildMdpOutput(builder, "y");
252     efield_[ZZ].buildMdpOutput(builder, "z");
253 }
254
255 void ElectricField::initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv)
256 {
257     if (isActive())
258     {
259         please_cite(fplog, "Caleman2008a");
260
261         // Optional output file showing the field, see manual.
262         if (opt2bSet("-field", nfile, fnm))
263         {
264             if (bAppendFiles)
265             {
266                 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
267             }
268             else
269             {
270                 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
271                                     "Applied electric field",
272                                     "Time (ps)",
273                                     "E (V/nm)",
274                                     oenv);
275             }
276         }
277     }
278 }
279
280 void ElectricField::finishOutput()
281 {
282     if (fpField_ != nullptr)
283     {
284         // This is opened sometimes with xvgropen, sometimes with
285         // gmx_fio_fopen, so we use the least common denominator for closing.
286         gmx_fio_fclose(fpField_);
287         fpField_ = nullptr;
288     }
289 }
290
291 real ElectricField::field(int dim, real t) const
292 {
293     return efield_[dim].evaluate(t);
294 }
295
296 bool ElectricField::isActive() const
297 {
298     return (efield_[XX].a() != 0 || efield_[YY].a() != 0 || efield_[ZZ].a() != 0);
299 }
300
301 void ElectricField::printComponents(double t) const
302 {
303     fprintf(fpField_, "%10g  %10g  %10g  %10g\n", t, field(XX, t), field(YY, t), field(ZZ, t));
304 }
305
306 void ElectricField::calculateForces(const ForceProviderInput& forceProviderInput,
307                                     ForceProviderOutput*      forceProviderOutput)
308 {
309     if (isActive())
310     {
311         const double     t  = forceProviderInput.t_;
312         const t_commrec& cr = forceProviderInput.cr_;
313
314         // NOTE: The non-conservative electric field does not have a virial
315         ArrayRef<RVec> f = forceProviderOutput->forceWithVirial_.force_;
316
317         auto chargeA = forceProviderInput.chargeA_;
318         for (int m = 0; (m < DIM); m++)
319         {
320             const real fieldStrength = gmx::c_fieldfac * field(m, t);
321
322             if (fieldStrength != 0)
323             {
324                 // TODO: Check parallellism
325                 for (int i = 0; i < forceProviderInput.homenr_; ++i)
326                 {
327                     // NOTE: Not correct with perturbed charges
328                     f[i][m] += chargeA[i] * fieldStrength;
329                 }
330             }
331         }
332         if (MASTER(&cr) && fpField_ != nullptr)
333         {
334             printComponents(t);
335         }
336     }
337 }
338
339 } // namespace
340
341 std::unique_ptr<IMDModule> createElectricFieldModule()
342 {
343     return std::make_unique<ElectricField>();
344 }
345
346 } // namespace gmx