Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / applied-forces / electricfield.cpp
blobb75d479293e36d44de84566bcd9e3a24780c1ade
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Declares data structure and utilities for electric fields
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_applied_forces
42 #include "gmxpre.h"
44 #include "electricfield.h"
46 #include <cmath>
48 #include "gromacs/commandline/filenm.h"
49 #include "gromacs/fileio/gmxfio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/mdtypes/forcerec.h"
55 #include "gromacs/mdtypes/iforceprovider.h"
56 #include "gromacs/mdtypes/imdmodule.h"
57 #include "gromacs/mdtypes/imdoutputprovider.h"
58 #include "gromacs/mdtypes/imdpoptionprovider.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/options/basicoptions.h"
61 #include "gromacs/options/ioptionscontainerwithsections.h"
62 #include "gromacs/options/optionsection.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/keyvaluetreebuilder.h"
65 #include "gromacs/utility/keyvaluetreetransform.h"
66 #include "gromacs/utility/pleasecite.h"
67 #include "gromacs/utility/strconvert.h"
68 #include "gromacs/utility/stringutil.h"
70 namespace gmx
73 namespace
76 /*! \internal
77 * \brief Declaration of storage unit for fields
79 class ElectricFieldData
81 public:
82 ElectricFieldData() : a_(0), omega_(0), t0_(0), sigma_(0)
86 /*! \brief
87 * Adds an option section to specify parameters for this field component.
89 void initMdpOptions(IOptionsContainerWithSections *options, const char *sectionName)
91 auto section = options->addSection(OptionSection(sectionName));
92 section.addOption(RealOption("E0").store(&a_));
93 section.addOption(RealOption("omega").store(&omega_));
94 section.addOption(RealOption("t0").store(&t0_));
95 section.addOption(RealOption("sigma").store(&sigma_));
98 /*! \brief Evaluates this field component at given time.
100 * \param[in] t The time to evualate at
101 * \return The electric field
103 real evaluate(real t) const
105 if (sigma_ > 0)
107 return a_ * (std::cos(omega_*(t-t0_))
108 * std::exp(-square(t-t0_)/(2.0*square(sigma_))));
110 else
112 return a_ * std::cos(omega_*t);
116 /*! \brief Initiate the field values
118 * \param[in] a Amplitude
119 * \param[in] omega Frequency
120 * \param[in] t0 Peak of the pulse
121 * \param[in] sigma Width of the pulse
123 void setField(real a, real omega, real t0, real sigma)
125 a_ = a;
126 omega_ = omega;
127 t0_ = t0;
128 sigma_ = sigma;
131 //! Return the amplitude
132 real a() const { return a_; }
133 //! Return the frequency
134 real omega() const { return omega_; }
135 //! Return the time for the peak of the pulse
136 real t0() const { return t0_; }
137 //! Return the width of the pulse (0 means inifinite)
138 real sigma() const { return sigma_; }
140 private:
141 //! Coeffient (V / nm)
142 real a_;
143 //! Frequency (1/ps)
144 real omega_;
145 //! Central time point (ps) for pulse
146 real t0_;
147 //! Width of pulse (ps, if zero there is no pulse)
148 real sigma_;
151 /*! \internal
152 * \brief Describe time dependent electric field
154 * Class that implements a force to be evaluated in mdrun.
155 * The electric field can be pulsed and oscillating, simply
156 * oscillating, or static, in each of X,Y,Z directions.
158 class ElectricField : public IMDModule,
159 public IMdpOptionProvider, public IMDOutputProvider,
160 public IForceProvider
162 public:
163 ElectricField() : fpField_(nullptr) {}
165 // From IMDModule
166 virtual IMdpOptionProvider *mdpOptionProvider() { return this; }
167 virtual IMDOutputProvider *outputProvider() { return this; }
168 virtual IForceProvider *forceProvider() { return this; }
170 // From IMdpOptionProvider
171 virtual void initMdpTransform(IKeyValueTreeTransformRules *transform);
172 virtual void initMdpOptions(IOptionsContainerWithSections *options);
174 // From IMDOutputProvider
175 virtual void initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
176 bool bAppendFiles, const gmx_output_env_t *oenv);
177 virtual void finishOutput();
179 // From IForceProvider
180 virtual void initForcerec(t_forcerec *fr);
181 //! \copydoc IForceProvider::calculateForces()
182 virtual void calculateForces(const t_commrec *cr,
183 const t_mdatoms *mdatoms,
184 PaddedRVecVector *force,
185 double t);
187 private:
188 //! Return whether or not to apply a field
189 bool isActive() const;
191 /*! \brief Return the field strength
193 * \param[in] dim The spatial direction
194 * \param[in] t The time (ps)
195 * \return The field strength in V/nm units
197 real field(int dim, real t) const;
199 /*! \brief Return amplitude of field
201 * \param[in] dim Direction of the field (XX, YY, ZZ)
202 * \return Amplitude of the field
204 real a(int dim) const { return efield_[dim].a(); }
205 /*! \brief Return frequency of field (1/ps)
207 * \param[in] dim Direction of the field (XX, YY, ZZ)
208 * \return Frequency of the field
210 real omega(int dim) const { return efield_[dim].omega(); }
211 /*! \brief Return time of pulse peak
213 * \param[in] dim Direction of the field (XX, YY, ZZ)
214 * \return Time of pulse peak
216 real t0(int dim) const { return efield_[dim].t0(); }
217 /*! \brief Return width of the pulse
219 * \param[in] dim Direction of the field (XX, YY, ZZ)
220 * \return Width of the pulse
222 real sigma(int dim) const { return efield_[dim].sigma(); }
224 /*! \brief Print the field components to a file
226 * \param[in] t The time
227 * Will throw and exit with fatal error if file is not open.
229 void printComponents(double t) const;
231 //! The field strength in each dimension
232 ElectricFieldData efield_[DIM];
233 //! File pointer for electric field
234 FILE *fpField_;
237 //! Converts static parameters from mdp format to E0.
238 real convertStaticParameters(const std::string &value)
240 // TODO: Better context for the exceptions here (possibly
241 // also convert them to warning_errors or such).
242 const std::vector<std::string> sx = splitString(value);
243 if (sx.empty())
245 return 0.0;
247 const int n = fromString<int>(sx[0]);
248 if (n <= 0)
250 return 0.0;
252 if (n != 1)
254 GMX_THROW(InvalidInputError("Only one electric field term supported for each dimension"));
256 if (sx.size() != 3)
258 GMX_THROW(InvalidInputError("Expected exactly one electric field amplitude value"));
260 return fromString<real>(sx[1]);
263 //! Converts dynamic parameters from mdp format to (omega, t0, sigma).
264 void convertDynamicParameters(gmx::KeyValueTreeObjectBuilder *builder,
265 const std::string &value)
267 const std::vector<std::string> sxt = splitString(value);
268 if (sxt.empty())
270 return;
272 const int n = fromString<int>(sxt[0]);
273 switch (n)
275 case 1:
276 if (sxt.size() != 3)
278 GMX_THROW(InvalidInputError("Please specify 1 omega 0 for non-pulsed fields"));
280 builder->addValue<real>("omega", fromString<real>(sxt[1]));
281 break;
282 case 3:
283 if (sxt.size() != 7)
285 GMX_THROW(InvalidInputError("Please specify 1 omega 0 t0 0 sigma 0 for pulsed fields"));
287 builder->addValue<real>("omega", fromString<real>(sxt[1]));
288 builder->addValue<real>("t0", fromString<real>(sxt[3]));
289 builder->addValue<real>("sigma", fromString<real>(sxt[5]));
290 break;
291 default:
292 GMX_THROW(InvalidInputError("Incomprehensible input for electric field"));
296 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules *rules)
298 // TODO This responsibility should be handled by the caller,
299 // e.g. embedded in the rules, somehow.
300 std::string prefix = "/applied-forces/";
302 rules->addRule().from<std::string>("/E-x").to<real>(prefix + "electric-field/x/E0")
303 .transformWith(&convertStaticParameters);
304 rules->addRule().from<std::string>("/E-xt").toObject(prefix + "electric-field/x")
305 .transformWith(&convertDynamicParameters);
306 rules->addRule().from<std::string>("/E-y").to<real>(prefix + "electric-field/y/E0")
307 .transformWith(&convertStaticParameters);
308 rules->addRule().from<std::string>("/E-yt").toObject(prefix + "electric-field/y")
309 .transformWith(&convertDynamicParameters);
310 rules->addRule().from<std::string>("/E-z").to<real>(prefix + "electric-field/z/E0")
311 .transformWith(&convertStaticParameters);
312 rules->addRule().from<std::string>("/E-zt").toObject(prefix + "electric-field/z")
313 .transformWith(&convertDynamicParameters);
316 void ElectricField::initMdpOptions(IOptionsContainerWithSections *options)
318 //CTYPE ("Format is E0 (V/nm), omega (1/ps), t0 (ps), sigma (ps) ");
319 auto section = options->addSection(OptionSection("electric-field"));
320 efield_[XX].initMdpOptions(&section, "x");
321 efield_[YY].initMdpOptions(&section, "y");
322 efield_[ZZ].initMdpOptions(&section, "z");
325 void ElectricField::initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
326 bool bAppendFiles, const gmx_output_env_t *oenv)
328 if (isActive())
330 please_cite(fplog, "Caleman2008a");
332 // Optional outpuf file showing the field, see manual.
333 if (opt2bSet("-field", nfile, fnm))
335 if (bAppendFiles)
337 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
339 else
341 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
342 "Applied electric field", "Time (ps)",
343 "E (V/nm)", oenv);
349 void ElectricField::finishOutput()
351 if (fpField_ != nullptr)
353 // This is opened sometimes with xvgropen, sometimes with
354 // gmx_fio_fopen, so we use the least common denominator for closing.
355 gmx_fio_fclose(fpField_);
356 fpField_ = nullptr;
360 void ElectricField::initForcerec(t_forcerec *fr)
362 if (isActive())
364 fr->bF_NoVirSum = TRUE;
365 fr->efield = this;
369 real ElectricField::field(int dim, real t) const
371 return efield_[dim].evaluate(t);
374 bool ElectricField::isActive() const
376 return (efield_[XX].a() != 0 ||
377 efield_[YY].a() != 0 ||
378 efield_[ZZ].a() != 0);
381 void ElectricField::printComponents(double t) const
383 fprintf(fpField_, "%10g %10g %10g %10g\n", t,
384 field(XX, t), field(YY, t), field(ZZ, t));
387 void ElectricField::calculateForces(const t_commrec *cr,
388 const t_mdatoms *mdatoms,
389 PaddedRVecVector *force,
390 double t)
392 if (isActive())
394 rvec *f = as_rvec_array(force->data());
396 for (int m = 0; (m < DIM); m++)
398 real Ext = FIELDFAC*field(m, t);
400 if (Ext != 0)
402 // TODO: Check parallellism
403 for (int i = 0; i < mdatoms->homenr; ++i)
405 // NOTE: Not correct with perturbed charges
406 f[i][m] += mdatoms->chargeA[i]*Ext;
410 if (MASTER(cr) && fpField_ != nullptr)
412 printComponents(t);
417 } // namespace
419 std::unique_ptr<IMDModule> createElectricFieldModule()
421 return std::unique_ptr<IMDModule>(new ElectricField());
424 } // namespace gmx