2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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
44 #include "electricfield.h"
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"
79 * \brief Declaration of storage unit for fields
81 class ElectricFieldData
84 ElectricFieldData() : a_(0), omega_(0), t0_(0), sigma_(0)
89 * Adds an option section to specify parameters for this field component.
91 void initMdpOptions(IOptionsContainerWithSections *options, const char *sectionName)
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_));
100 * Creates mdp parameters for this field component.
102 void buildMdpOutput(KeyValueTreeObjectBuilder *builder, const std::string &name) const
104 builder->addUniformArray<real>("electric-field-" + name, {a_, omega_, t0_, sigma_});
107 /*! \brief Evaluates this field component at given time.
109 * \param[in] t The time to evualate at
110 * \return The electric field
112 real evaluate(real t) const
116 return a_ * (std::cos(omega_*(t-t0_))
117 * std::exp(-square(t-t0_)/(2.0*square(sigma_))));
121 return a_ * std::cos(omega_*t);
125 //! Return the amplitude
126 real a() const { return a_; }
129 //! Coeffient (V / nm)
133 //! Central time point (ps) for pulse
135 //! Width of pulse (ps, if zero there is no pulse)
140 * \brief Describe time dependent electric field
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.
146 class ElectricField final : public IMDModule,
147 public IMdpOptionProvider, public IMDOutputProvider,
148 public IForceProvider
151 ElectricField() : fpField_(nullptr) {}
154 IMdpOptionProvider *mdpOptionProvider() override { return this; }
155 IMDOutputProvider *outputProvider() override { return this; }
156 void initForceProviders(ForceProviders *forceProviders) override
160 forceProviders->addForceProvider(this);
164 // From IMdpOptionProvider
165 void initMdpTransform(IKeyValueTreeTransformRules *transform) override;
166 void initMdpOptions(IOptionsContainerWithSections *options) override;
167 void buildMdpOutput(KeyValueTreeObjectBuilder *builder) const override;
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;
174 // From IForceProvider
175 //! \copydoc IForceProvider::calculateForces()
176 void calculateForces(const ForceProviderInput &forceProviderInput,
177 ForceProviderOutput *forceProviderOutput) override;
180 //! Return whether or not to apply a field
181 bool isActive() const;
183 /*! \brief Return the field strength
185 * \param[in] dim The spatial direction
186 * \param[in] t The time (ps)
187 * \return The field strength in V/nm units
189 real field(int dim, real t) const;
191 /*! \brief Print the field components to a file
193 * \param[in] t The time
194 * Will throw and exit with fatal error if file is not open.
196 void printComponents(double t) const;
198 //! The field strength in each dimension
199 ElectricFieldData efield_[DIM];
200 //! File pointer for electric field
204 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
205 void convertParameters(gmx::KeyValueTreeObjectBuilder *builder,
206 const std::string &value)
208 const std::vector<std::string> sxt = splitString(value);
215 GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
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]));
223 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules *rules)
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);
233 void ElectricField::initMdpOptions(IOptionsContainerWithSections *options)
235 auto section = options->addSection(OptionSection("electric-field"));
236 efield_[XX].initMdpOptions(§ion, "x");
237 efield_[YY].initMdpOptions(§ion, "y");
238 efield_[ZZ].initMdpOptions(§ion, "z");
241 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder *builder) const
243 const char *const comment[] = {
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."
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");
256 void ElectricField::initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
257 bool bAppendFiles, const gmx_output_env_t *oenv)
261 please_cite(fplog, "Caleman2008a");
263 // Optional outpuf file showing the field, see manual.
264 if (opt2bSet("-field", nfile, fnm))
268 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
272 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
273 "Applied electric field", "Time (ps)",
280 void ElectricField::finishOutput()
282 if (fpField_ != nullptr)
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_);
291 real ElectricField::field(int dim, real t) const
293 return efield_[dim].evaluate(t);
296 bool ElectricField::isActive() const
298 return (efield_[XX].a() != 0 ||
299 efield_[YY].a() != 0 ||
300 efield_[ZZ].a() != 0);
303 void ElectricField::printComponents(double t) const
305 fprintf(fpField_, "%10g %10g %10g %10g\n", t,
306 field(XX, t), field(YY, t), field(ZZ, t));
309 void ElectricField::calculateForces(const ForceProviderInput &forceProviderInput,
310 ForceProviderOutput *forceProviderOutput)
314 const t_mdatoms &mdatoms = forceProviderInput.mdatoms_;
315 const double t = forceProviderInput.t_;
316 const t_commrec &cr = forceProviderInput.cr_;
318 // NOTE: The non-conservative electric field does not have a virial
319 rvec *f = as_rvec_array(forceProviderOutput->forceWithVirial_.force_.data());
321 for (int m = 0; (m < DIM); m++)
323 real Ext = FIELDFAC*field(m, t);
327 // TODO: Check parallellism
328 for (int i = 0; i < mdatoms.homenr; ++i)
330 // NOTE: Not correct with perturbed charges
331 f[i][m] += mdatoms.chargeA[i]*Ext;
335 if (MASTER(&cr) && fpField_ != nullptr)
344 std::unique_ptr<IMDModule> createElectricFieldModule()
346 return std::unique_ptr<IMDModule>(new ElectricField());