2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019, The GROMACS development team.
5 * Copyright (c) 2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * Defines data structure and utilities for electric fields
40 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 * \ingroup module_applied_forces
45 #include "electricfield.h"
52 #include "gromacs/commandline/filenm.h"
53 #include "gromacs/fileio/gmxfio.h"
54 #include "gromacs/fileio/xvgr.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/forceoutput.h"
58 #include "gromacs/mdtypes/iforceprovider.h"
59 #include "gromacs/mdtypes/imdmodule.h"
60 #include "gromacs/mdtypes/imdoutputprovider.h"
61 #include "gromacs/mdtypes/imdpoptionprovider.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"
78 * \brief Describes an applied electric field in a coordinate
81 * Can compute the applied field strength at a time, and supports
82 * operations to read and form corresponding .mdp contents.
84 class ElectricFieldDimension
88 * Adds an option section to specify parameters for this field component.
90 void initMdpOptions(IOptionsContainerWithSections* options, const char* sectionName)
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_));
99 * Creates mdp parameters for this field component.
101 void buildMdpOutput(KeyValueTreeObjectBuilder* builder, const std::string& name) const
103 builder->addUniformArray<real>("electric-field-" + name, { a_, omega_, t0_, sigma_ });
106 /*! \brief Evaluates this field component at given time.
108 * \param[in] t The time to evualate at
109 * \return The electric field
111 real evaluate(real t) const
115 return a_ * (std::cos(omega_ * (t - t0_)) * std::exp(-square(t - t0_) / (2.0 * square(sigma_))));
119 return a_ * std::cos(omega_ * t);
123 //! Return the amplitude
124 real a() const { return a_; }
127 //! Coeffient (V / nm)
131 //! Central time point (ps) for pulse
133 //! Width of pulse (ps, if zero there is no pulse)
138 * \brief Describe time dependent electric field
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.
144 class ElectricField final : public IMDModule, public IMdpOptionProvider, public IMDOutputProvider, public IForceProvider
147 ElectricField() : fpField_(nullptr) {}
150 IMdpOptionProvider* mdpOptionProvider() override { return this; }
151 IMDOutputProvider* outputProvider() override { return this; }
152 void initForceProviders(ForceProviders* forceProviders) override
156 forceProviders->addForceProvider(this);
160 // From IMdpOptionProvider
161 void initMdpTransform(IKeyValueTreeTransformRules* transform) override;
162 void initMdpOptions(IOptionsContainerWithSections* options) override;
163 void buildMdpOutput(KeyValueTreeObjectBuilder* builder) const override;
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;
169 // From IForceProvider
170 //! \copydoc IForceProvider::calculateForces()
171 void calculateForces(const ForceProviderInput& forceProviderInput,
172 ForceProviderOutput* forceProviderOutput) override;
174 void subscribeToSimulationSetupNotifications(MDModulesNotifiers* /* notifiers */) override {}
175 void subscribeToPreProcessingNotifications(MDModulesNotifiers* /* notifiers */) override {}
178 //! Return whether or not to apply a field
179 bool isActive() const;
181 /*! \brief Return the field strength
183 * \param[in] dim The spatial direction
184 * \param[in] t The time (ps)
185 * \return The field strength in V/nm units
187 real field(int dim, real t) const;
189 /*! \brief Print the field components to a file
191 * \param[in] t The time
192 * Will throw and exit with fatal error if file is not open.
194 void printComponents(double t) const;
196 //! The components of the applied electric field in each coordinate dimension
197 ElectricFieldDimension efield_[DIM];
198 //! File pointer for electric field
202 //! Converts dynamic parameters from new mdp format to (E0, omega, t0, sigma).
203 void convertParameters(gmx::KeyValueTreeObjectBuilder* builder, const std::string& value)
205 const std::vector<std::string> sxt = splitString(value);
212 GMX_THROW(InvalidInputError("Please specify E0 omega t0 sigma for electric fields"));
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]));
220 void ElectricField::initMdpTransform(IKeyValueTreeTransformRules* rules)
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);
227 void ElectricField::initMdpOptions(IOptionsContainerWithSections* options)
229 auto section = options->addSection(OptionSection("electric-field"));
230 efield_[XX].initMdpOptions(§ion, "x");
231 efield_[YY].initMdpOptions(§ion, "y");
232 efield_[ZZ].initMdpOptions(§ion, "z");
235 void ElectricField::buildMdpOutput(KeyValueTreeObjectBuilder* builder) const
237 std::string comment =
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");
249 void ElectricField::initOutput(FILE* fplog, int nfile, const t_filenm fnm[], bool bAppendFiles, const gmx_output_env_t* oenv)
253 please_cite(fplog, "Caleman2008a");
255 // Optional output file showing the field, see manual.
256 if (opt2bSet("-field", nfile, fnm))
260 fpField_ = gmx_fio_fopen(opt2fn("-field", nfile, fnm), "a+");
264 fpField_ = xvgropen(opt2fn("-field", nfile, fnm),
265 "Applied electric field",
274 void ElectricField::finishOutput()
276 if (fpField_ != nullptr)
278 // This is opened sometimes with xvgropen, sometimes with
279 // gmx_fio_fopen, so we use the least common denominator for closing.
280 gmx_fio_fclose(fpField_);
285 real ElectricField::field(int dim, real t) const
287 return efield_[dim].evaluate(t);
290 bool ElectricField::isActive() const
292 return (efield_[XX].a() != 0 || efield_[YY].a() != 0 || efield_[ZZ].a() != 0);
295 void ElectricField::printComponents(double t) const
297 fprintf(fpField_, "%10g %10g %10g %10g\n", t, field(XX, t), field(YY, t), field(ZZ, t));
300 void ElectricField::calculateForces(const ForceProviderInput& forceProviderInput,
301 ForceProviderOutput* forceProviderOutput)
305 const double t = forceProviderInput.t_;
306 const t_commrec& cr = forceProviderInput.cr_;
308 // NOTE: The non-conservative electric field does not have a virial
309 ArrayRef<RVec> f = forceProviderOutput->forceWithVirial_.force_;
311 auto chargeA = forceProviderInput.chargeA_;
312 for (int m = 0; (m < DIM); m++)
314 const real fieldStrength = gmx::c_fieldfac * field(m, t);
316 if (fieldStrength != 0)
318 // TODO: Check parallellism
319 for (int i = 0; i < forceProviderInput.homenr_; ++i)
321 // NOTE: Not correct with perturbed charges
322 f[i][m] += chargeA[i] * fieldStrength;
326 if (MASTER(&cr) && fpField_ != nullptr)
335 std::unique_ptr<IMDModule> createElectricFieldModule()
337 return std::make_unique<ElectricField>();