Write copyright output files as utf8
[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,2020, 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     void subscribeToSimulationSetupNotifications(MdModulesNotifier* /* notifier */) override {}
175     void subscribeToPreProcessingNotifications(MdModulesNotifier* /* notifier */) override {}
176
177 private:
178     //! Return whether or not to apply a field
179     bool isActive() const;
180
181     /*! \brief Return the field strength
182      *
183      * \param[in] dim The spatial direction
184      * \param[in] t   The time (ps)
185      * \return The field strength in V/nm units
186      */
187     real field(int dim, real t) const;
188
189     /*! \brief Print the field components to a file
190      *
191      * \param[in] t   The time
192      * Will throw and exit with fatal error if file is not open.
193      */
194     void printComponents(double t) const;
195
196     //! The components of the applied electric field in each coordinate dimension
197     ElectricFieldDimension efield_[DIM];
198     //! File pointer for electric field
199     FILE* fpField_;
200 };
201
202 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
203 void convertParameters(gmx::KeyValueTreeObjectBuilder* builder, const std::string& value)
204 {
205     const std::vector<std::string> sxt = splitString(value);
206     if (sxt.empty())
207     {
208         return;
209     }
210     if (sxt.size() != 4)
211     {
212         GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
213     }
214     builder->addValue<real>("E0", fromString<real>(sxt[0]));
215     builder->addValue<real>("omega", fromString<real>(sxt[1]));
216     builder->addValue<real>("t0", fromString<real>(sxt[2]));
217     builder->addValue<real>("sigma", fromString<real>(sxt[3]));
218 }
219
220 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules* rules)
221 {
222     rules->addRule().from<std::string>("/electric-field-x").toObject("/electric-field/x").transformWith(&convertParameters);
223     rules->addRule().from<std::string>("/electric-field-y").toObject("/electric-field/y").transformWith(&convertParameters);
224     rules->addRule().from<std::string>("/electric-field-z").toObject("/electric-field/z").transformWith(&convertParameters);
225 }
226
227 void ElectricField::initMdpOptions(IOptionsContainerWithSections* options)
228 {
229     auto section = options->addSection(OptionSection("electric-field"));
230     efield_[XX].initMdpOptions(&section, "x");
231     efield_[YY].initMdpOptions(&section, "y");
232     efield_[ZZ].initMdpOptions(&section, "z");
233 }
234
235 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
236 {
237     std::string comment =
238             R"(; Electric fields
239 ; Format for electric-field-x, etc. is: four real variables:
240 ; amplitude (V/nm), frequency omega (1/ps), time for the pulse peak (ps),
241 ; and sigma (ps) width of the pulse. Omega = 0 means static field,
242 ; sigma = 0 means no pulse, leaving the field to be a cosine function.)";
243     builder->addValue<std::string>("comment-electric-field", comment);
244     efield_[XX].buildMdpOutput(builder, "x");
245     efield_[YY].buildMdpOutput(builder, "y");
246     efield_[ZZ].buildMdpOutput(builder, "z");
247 }
248
249 void ElectricField::initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv)
250 {
251     if (isActive())
252     {
253         please_cite(fplog, "Caleman2008a");
254
255         // Optional output file showing the field, see manual.
256         if (opt2bSet("-field", nfile, fnm))
257         {
258             if (bAppendFiles)
259             {
260                 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
261             }
262             else
263             {
264                 fpField_ = xvgropen(opt2fn("-field", nfile, fnm), "Applied electric field",
265                                     "Time (ps)", "E (V/nm)", oenv);
266             }
267         }
268     }
269 }
270
271 void ElectricField::finishOutput()
272 {
273     if (fpField_ != nullptr)
274     {
275         // This is opened sometimes with xvgropen, sometimes with
276         // gmx_fio_fopen, so we use the least common denominator for closing.
277         gmx_fio_fclose(fpField_);
278         fpField_ = nullptr;
279     }
280 }
281
282 real ElectricField::field(int dim, real t) const
283 {
284     return efield_[dim].evaluate(t);
285 }
286
287 bool ElectricField::isActive() const
288 {
289     return (efield_[XX].a() != 0 || efield_[YY].a() != 0 || efield_[ZZ].a() != 0);
290 }
291
292 void ElectricField::printComponents(double t) const
293 {
294     fprintf(fpField_, "%10g  %10g  %10g  %10g\n", t, field(XX, t), field(YY, t), field(ZZ, t));
295 }
296
297 void ElectricField::calculateForces(const ForceProviderInput& forceProviderInput,
298                                     ForceProviderOutput*      forceProviderOutput)
299 {
300     if (isActive())
301     {
302         const t_mdatoms& mdatoms = forceProviderInput.mdatoms_;
303         const double     t       = forceProviderInput.t_;
304         const t_commrec& cr      = forceProviderInput.cr_;
305
306         // NOTE: The non-conservative electric field does not have a virial
307         ArrayRef<RVec> f = forceProviderOutput->forceWithVirial_.force_;
308
309         for (int m = 0; (m < DIM); m++)
310         {
311             const real fieldStrength = FIELDFAC * field(m, t);
312
313             if (fieldStrength != 0)
314             {
315                 // TODO: Check parallellism
316                 for (index i = 0; i != ssize(f); ++i)
317                 {
318                     // NOTE: Not correct with perturbed charges
319                     f[i][m] += mdatoms.chargeA[i] * fieldStrength;
320                 }
321             }
322         }
323         if (MASTER(&cr) && fpField_ != nullptr)
324         {
325             printComponents(t);
326         }
327     }
328 }
329
330 } // namespace
331
332 std::unique_ptr<IMDModule> createElectricFieldModule()
333 {
334     return std::make_unique<ElectricField>();
335 }
336
337 } // namespace gmx