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 * Tests for functionality of the "angle" trajectory analysis module.
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_applied_forces
44 #include "gromacs/applied_forces/electricfield.h"
46 #include <gtest/gtest.h>
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/math/paddedvector.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/forcerec.h"
52 #include "gromacs/mdtypes/enerdata.h"
53 #include "gromacs/mdtypes/forceoutput.h"
54 #include "gromacs/mdtypes/iforceprovider.h"
55 #include "gromacs/mdtypes/imdmodule.h"
56 #include "gromacs/mdtypes/imdpoptionprovider.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/options/options.h"
60 #include "gromacs/options/treesupport.h"
61 #include "gromacs/utility/arrayref.h"
62 #include "gromacs/utility/keyvaluetreebuilder.h"
63 #include "gromacs/utility/keyvaluetreetransform.h"
64 #include "gromacs/utility/real.h"
65 #include "gromacs/utility/smalloc.h"
66 #include "gromacs/utility/stringcompare.h"
67 #include "gromacs/utility/stringutil.h"
69 #include "testutils/testasserts.h"
74 /********************************************************************
78 class ElectricFieldTest : public ::testing::Test
88 gmx::test::FloatingPointTolerance tolerance(
89 gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
90 auto module(gmx::createElectricFieldModule());
93 const char *dimXYZ[3] = { "x", "y", "z" };
94 GMX_RELEASE_ASSERT(dim >= 0 && dim < DIM, "Dimension should be 0, 1 or 2");
96 gmx::KeyValueTreeBuilder mdpValues;
97 mdpValues.rootObject().addValue(gmx::formatString("electric-field-%s", dimXYZ[dim]),
98 gmx::formatString("%g %g %g %g", E0, omega, t0, sigma));
100 gmx::KeyValueTreeTransformer transform;
101 transform.rules()->addRule()
102 .keyMatchType("/", gmx::StringCompareType::CaseAndDashInsensitive);
103 module->mdpOptionProvider()->initMdpTransform(transform.rules());
104 auto result = transform.transform(mdpValues.build(), nullptr);
105 gmx::Options moduleOptions;
106 module->mdpOptionProvider()->initMdpOptions(&moduleOptions);
107 gmx::assignOptionsFromKeyValueTree(&moduleOptions, result.object(), nullptr);
109 ForceProviders forceProviders;
110 module->initForceProviders(&forceProviders);
113 PaddedVector<gmx::RVec> f = { {0, 0, 0} };
114 gmx::ForceWithVirial forceWithVirial(f, true);
116 snew(md.chargeA, md.homenr);
119 t_commrec *cr = init_commrec();
120 matrix boxDummy = { {0, 0, 0}, {0, 0, 0}, {0, 0, 0} };
121 gmx_enerdata_t enerdDummy;
123 gmx::ForceProviderInput forceProviderInput(gmx::EmptyArrayRef {}, md, 0.0, boxDummy, *cr);
124 gmx::ForceProviderOutput forceProviderOutput(&forceWithVirial, &enerdDummy);
125 forceProviders.calculateForces(forceProviderInput, &forceProviderOutput);
128 EXPECT_REAL_EQ_TOL(f[0][dim], expectedValue, tolerance);
133 TEST_F(ElectricFieldTest, Static)
135 test(0, 1, 0, 0, 0, 96.4853363);
138 TEST_F(ElectricFieldTest, Oscillating)
140 test(0, 1, 5, 0.2, 0, 96.4853363);
143 TEST_F(ElectricFieldTest, Pulsed)
145 test(0, 1, 5, 0.5, 1, -68.215782);