5bc278142d330498555d45d9a2f2e26e8992a80d
[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, 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  * Defines data structure and utilities for electric fields
38  *
39  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40  * \ingroup module_applied_forces
41  */
42 #include "gmxpre.h"
43
44 #include "electricfield.h"
45
46 #include <cmath>
47
48 #include <memory>
49 #include <string>
50
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/fileio/xvgr.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/mdtypes/commrec.h"
56 #include "gromacs/mdtypes/forceoutput.h"
57 #include "gromacs/mdtypes/iforceprovider.h"
58 #include "gromacs/mdtypes/imdmodule.h"
59 #include "gromacs/mdtypes/imdoutputprovider.h"
60 #include "gromacs/mdtypes/imdpoptionprovider.h"
61 #include "gromacs/mdtypes/mdatom.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     }
98     /*! \brief
99      * Creates mdp parameters for this field component.
100      */
101     void buildMdpOutput(KeyValueTreeObjectBuilder* builder, const std::string& name) const
102     {
103         builder->addUniformArray<real>("electric-field-" + name, { a_, omega_, t0_, sigma_ });
104     }
105
106     /*! \brief Evaluates this field component at given time.
107      *
108      * \param[in] t The time to evualate at
109      * \return The electric field
110      */
111     real evaluate(real t) const
112     {
113         if (sigma_ > 0)
114         {
115             return a_ * (std::cos(omega_ * (t - t0_)) * std::exp(-square(t - t0_) / (2.0 * square(sigma_))));
116         }
117         else
118         {
119             return a_ * std::cos(omega_ * t);
120         }
121     }
122
123     //! Return the amplitude
124     real a() const { return a_; }
125
126 private:
127     //! Coeffient (V / nm)
128     real a_ = 0;
129     //! Frequency (1/ps)
130     real omega_ = 0;
131     //! Central time point (ps) for pulse
132     real t0_ = 0;
133     //! Width of pulse (ps, if zero there is no pulse)
134     real sigma_ = 0;
135 };
136
137 /*! \internal
138  * \brief Describe time dependent electric field
139  *
140  * Class that implements a force to be evaluated in mdrun.
141  * The electric field can be pulsed and oscillating, simply
142  * oscillating, or static, in each of X,Y,Z directions.
143  */
144 class ElectricField final : public IMDModule, public IMdpOptionProvider, public IMDOutputProvider, public IForceProvider
145 {
146 public:
147     ElectricField() : fpField_(nullptr) {}
148
149     // From IMDModule
150     IMdpOptionProvider* mdpOptionProvider() override { return this; }
151     IMDOutputProvider*  outputProvider() override { return this; }
152     void                initForceProviders(ForceProviders* forceProviders) override
153     {
154         if (isActive())
155         {
156             forceProviders->addForceProvider(this);
157         }
158     }
159
160     // From IMdpOptionProvider
161     void initMdpTransform(IKeyValueTreeTransformRules* transform) override;
162     void initMdpOptions(IOptionsContainerWithSections* options) override;
163     void buildMdpOutput(KeyValueTreeObjectBuilder* builder) const override;
164
165     // From IMDOutputProvider
166     void initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv) override;
167     void finishOutput() override;
168
169     // From IForceProvider
170     //! \copydoc IForceProvider::calculateForces()
171     void calculateForces(const ForceProviderInput& forceProviderInput,
172                          ForceProviderOutput*      forceProviderOutput) override;
173
174 private:
175     //! Return whether or not to apply a field
176     bool isActive() const;
177
178     /*! \brief Return the field strength
179      *
180      * \param[in] dim The spatial direction
181      * \param[in] t   The time (ps)
182      * \return The field strength in V/nm units
183      */
184     real field(int dim, real t) const;
185
186     /*! \brief Print the field components to a file
187      *
188      * \param[in] t   The time
189      * Will throw and exit with fatal error if file is not open.
190      */
191     void printComponents(double t) const;
192
193     //! The components of the applied electric field in each coordinate dimension
194     ElectricFieldDimension efield_[DIM];
195     //! File pointer for electric field
196     FILE* fpField_;
197 };
198
199 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
200 void convertParameters(gmx::KeyValueTreeObjectBuilder* builder, const std::string& value)
201 {
202     const std::vector<std::string> sxt = splitString(value);
203     if (sxt.empty())
204     {
205         return;
206     }
207     if (sxt.size() != 4)
208     {
209         GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
210     }
211     builder->addValue<real>("E0", fromString<real>(sxt[0]));
212     builder->addValue<real>("omega", fromString<real>(sxt[1]));
213     builder->addValue<real>("t0", fromString<real>(sxt[2]));
214     builder->addValue<real>("sigma", fromString<real>(sxt[3]));
215 }
216
217 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules* rules)
218 {
219     rules->addRule().from<std::string>("/electric-field-x").toObject("/electric-field/x").transformWith(&convertParameters);
220     rules->addRule().from<std::string>("/electric-field-y").toObject("/electric-field/y").transformWith(&convertParameters);
221     rules->addRule().from<std::string>("/electric-field-z").toObject("/electric-field/z").transformWith(&convertParameters);
222 }
223
224 void ElectricField::initMdpOptions(IOptionsContainerWithSections* options)
225 {
226     auto section = options->addSection(OptionSection("electric-field"));
227     efield_[XX].initMdpOptions(&section, "x");
228     efield_[YY].initMdpOptions(&section, "y");
229     efield_[ZZ].initMdpOptions(&section, "z");
230 }
231
232 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
233 {
234     std::string comment =
235             R"(; Electric fields
236 ; Format for electric-field-x, etc. is: four real variables:
237 ; amplitude (V/nm), frequency omega (1/ps), time for the pulse peak (ps),
238 ; and sigma (ps) width of the pulse. Omega = 0 means static field,
239 ; sigma = 0 means no pulse, leaving the field to be a cosine function.)";
240     builder->addValue<std::string>("comment-electric-field", comment);
241     efield_[XX].buildMdpOutput(builder, "x");
242     efield_[YY].buildMdpOutput(builder, "y");
243     efield_[ZZ].buildMdpOutput(builder, "z");
244 }
245
246 void ElectricField::initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv)
247 {
248     if (isActive())
249     {
250         please_cite(fplog, "Caleman2008a");
251
252         // Optional output file showing the field, see manual.
253         if (opt2bSet("-field", nfile, fnm))
254         {
255             if (bAppendFiles)
256             {
257                 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
258             }
259             else
260             {
261                 fpField_ = xvgropen(opt2fn("-field", nfile, fnm), "Applied electric field",
262                                     "Time (ps)", "E (V/nm)", oenv);
263             }
264         }
265     }
266 }
267
268 void ElectricField::finishOutput()
269 {
270     if (fpField_ != nullptr)
271     {
272         // This is opened sometimes with xvgropen, sometimes with
273         // gmx_fio_fopen, so we use the least common denominator for closing.
274         gmx_fio_fclose(fpField_);
275         fpField_ = nullptr;
276     }
277 }
278
279 real ElectricField::field(int dim, real t) const
280 {
281     return efield_[dim].evaluate(t);
282 }
283
284 bool ElectricField::isActive() const
285 {
286     return (efield_[XX].a() != 0 || efield_[YY].a() != 0 || efield_[ZZ].a() != 0);
287 }
288
289 void ElectricField::printComponents(double t) const
290 {
291     fprintf(fpField_, "%10g  %10g  %10g  %10g\n", t, field(XX, t), field(YY, t), field(ZZ, t));
292 }
293
294 void ElectricField::calculateForces(const ForceProviderInput& forceProviderInput,
295                                     ForceProviderOutput*      forceProviderOutput)
296 {
297     if (isActive())
298     {
299         const t_mdatoms& mdatoms = forceProviderInput.mdatoms_;
300         const double     t       = forceProviderInput.t_;
301         const t_commrec& cr      = forceProviderInput.cr_;
302
303         // NOTE: The non-conservative electric field does not have a virial
304         ArrayRef<RVec> f = forceProviderOutput->forceWithVirial_.force_;
305
306         for (int m = 0; (m < DIM); m++)
307         {
308             const real fieldStrength = FIELDFAC * field(m, t);
309
310             if (fieldStrength != 0)
311             {
312                 // TODO: Check parallellism
313                 for (index i = 0; i != ssize(f); ++i)
314                 {
315                     // NOTE: Not correct with perturbed charges
316                     f[i][m] += mdatoms.chargeA[i] * fieldStrength;
317                 }
318             }
319         }
320         if (MASTER(&cr) && fpField_ != nullptr)
321         {
322             printComponents(t);
323         }
324     }
325 }
326
327 } // namespace
328
329 std::unique_ptr<IMDModule> createElectricFieldModule()
330 {
331     return std::make_unique<ElectricField>();
332 }
333
334 } // namespace gmx