730aa55c432f63d8ff4a64d77175e3ced92d529c
[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, 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 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 <string>
49
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/vec.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 #include "gromacs/utility/stringutil.h"
71
72 namespace gmx
73 {
74
75 namespace
76 {
77
78 /*! \internal
79  * \brief Declaration of storage unit for fields
80  */
81 class ElectricFieldData
82 {
83     public:
84         ElectricFieldData() : a_(0), omega_(0), t0_(0), sigma_(0)
85         {
86         }
87
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_))
117                              * std::exp(-square(t-t0_)/(2.0*square(sigma_))));
118             }
119             else
120             {
121                 return a_ * std::cos(omega_*t);
122             }
123         }
124
125         //! Return the amplitude
126         real a()     const { return a_; }
127
128     private:
129         //! Coeffient (V / nm)
130         real a_;
131         //! Frequency (1/ps)
132         real omega_;
133         //! Central time point (ps) for pulse
134         real t0_;
135         //! Width of pulse (ps, if zero there is no pulse)
136         real sigma_;
137 };
138
139 /*! \internal
140  * \brief Describe time dependent electric field
141  *
142  * Class that implements a force to be evaluated in mdrun.
143  * The electric field can be pulsed and oscillating, simply
144  * oscillating, or static, in each of X,Y,Z directions.
145  */
146 class ElectricField final : public IMDModule,
147                             public IMdpOptionProvider, public IMDOutputProvider,
148                             public IForceProvider
149 {
150     public:
151         ElectricField() : fpField_(nullptr) {}
152
153         // From IMDModule
154         IMdpOptionProvider *mdpOptionProvider() override { return this; }
155         IMDOutputProvider *outputProvider() override { return this; }
156         void initForceProviders(ForceProviders *forceProviders) override
157         {
158             if (isActive())
159             {
160                 forceProviders->addForceProvider(this);
161             }
162         }
163
164         // From IMdpOptionProvider
165         void initMdpTransform(IKeyValueTreeTransformRules *transform) override;
166         void initMdpOptions(IOptionsContainerWithSections *options) override;
167         void buildMdpOutput(KeyValueTreeObjectBuilder *builder) const override;
168
169         // From IMDOutputProvider
170         void initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
171                         bool bAppendFiles, const gmx_output_env_t *oenv) override;
172         void finishOutput() override;
173
174         // From IForceProvider
175         //! \copydoc IForceProvider::calculateForces()
176         void calculateForces(const ForceProviderInput &forceProviderInput,
177                              ForceProviderOutput      *forceProviderOutput) override;
178
179     private:
180         //! Return whether or not to apply a field
181         bool isActive() const;
182
183         /*! \brief Return the field strength
184          *
185          * \param[in] dim The spatial direction
186          * \param[in] t   The time (ps)
187          * \return The field strength in V/nm units
188          */
189         real field(int dim, real t) const;
190
191         /*! \brief Print the field components to a file
192          *
193          * \param[in] t   The time
194          * Will throw and exit with fatal error if file is not open.
195          */
196         void printComponents(double t) const;
197
198         //! The field strength in each dimension
199         ElectricFieldData efield_[DIM];
200         //! File pointer for electric field
201         FILE             *fpField_;
202 };
203
204 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
205 void convertParameters(gmx::KeyValueTreeObjectBuilder *builder,
206                        const std::string              &value)
207 {
208     const std::vector<std::string> sxt = splitString(value);
209     if (sxt.empty())
210     {
211         return;
212     }
213     if (sxt.size() != 4)
214     {
215         GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
216     }
217     builder->addValue<real>("E0", fromString<real>(sxt[0]));
218     builder->addValue<real>("omega", fromString<real>(sxt[1]));
219     builder->addValue<real>("t0", fromString<real>(sxt[2]));
220     builder->addValue<real>("sigma", fromString<real>(sxt[3]));
221 }
222
223 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules *rules)
224 {
225     rules->addRule().from<std::string>("/electric-field-x").toObject("/electric-field/x")
226         .transformWith(&convertParameters);
227     rules->addRule().from<std::string>("/electric-field-y").toObject("/electric-field/y")
228         .transformWith(&convertParameters);
229     rules->addRule().from<std::string>("/electric-field-z").toObject("/electric-field/z")
230         .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     const char *const comment[] = {
244         "; 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     };
250     builder->addValue<std::string>("comment-electric-field", joinStrings(comment, "\n"));
251     efield_[XX].buildMdpOutput(builder, "x");
252     efield_[YY].buildMdpOutput(builder, "y");
253     efield_[ZZ].buildMdpOutput(builder, "z");
254 }
255
256 void ElectricField::initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
257                                bool bAppendFiles, const gmx_output_env_t *oenv)
258 {
259     if (isActive())
260     {
261         please_cite(fplog, "Caleman2008a");
262
263         // Optional outpuf file showing the field, see manual.
264         if (opt2bSet("-field", nfile, fnm))
265         {
266             if (bAppendFiles)
267             {
268                 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
269             }
270             else
271             {
272                 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
273                                     "Applied electric field", "Time (ps)",
274                                     "E (V/nm)", 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 ||
299             efield_[YY].a() != 0 ||
300             efield_[ZZ].a() != 0);
301 }
302
303 void ElectricField::printComponents(double t) const
304 {
305     fprintf(fpField_, "%10g  %10g  %10g  %10g\n", t,
306             field(XX, t), field(YY, t), field(ZZ, t));
307 }
308
309 void ElectricField::calculateForces(const ForceProviderInput &forceProviderInput,
310                                     ForceProviderOutput      *forceProviderOutput)
311 {
312     if (isActive())
313     {
314         const t_mdatoms &mdatoms = forceProviderInput.mdatoms_;
315         const double     t       = forceProviderInput.t_;
316         const t_commrec &cr      = forceProviderInput.cr_;
317
318         // NOTE: The non-conservative electric field does not have a virial
319         rvec *f = as_rvec_array(forceProviderOutput->forceWithVirial_.force_.data());
320
321         for (int m = 0; (m < DIM); m++)
322         {
323             real Ext = FIELDFAC*field(m, t);
324
325             if (Ext != 0)
326             {
327                 // TODO: Check parallellism
328                 for (int i = 0; i < mdatoms.homenr; ++i)
329                 {
330                     // NOTE: Not correct with perturbed charges
331                     f[i][m] += mdatoms.chargeA[i]*Ext;
332                 }
333             }
334         }
335         if (MASTER(&cr) && fpField_ != nullptr)
336         {
337             printComponents(t);
338         }
339     }
340 }
341
342 }   // namespace
343
344 std::unique_ptr<IMDModule> createElectricFieldModule()
345 {
346     return std::unique_ptr<IMDModule>(new ElectricField());
347 }
348
349 } // namespace gmx