Fix undefined behavior flagged by UBSAN
[gromacs.git] / src / gromacs / applied_forces / tests / electricfield.cpp
blob8e19955ae37654f42471b1793933e171251ec962
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 * Tests for functionality of the electric field module.
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_applied_forces
42 #include "gmxpre.h"
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/mdlib/forcerec.h"
51 #include "gromacs/mdtypes/commrec.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"
71 namespace gmx
73 namespace test
75 namespace
78 /********************************************************************
79 * ElectricFieldTest
82 class ElectricFieldTest : public ::testing::Test
84 public:
85 static void test(int dim, real E0, real omega, real t0, real sigma, real expectedValue)
87 // Make the electric field module
88 auto module = createElectricFieldModule();
90 // Fill the module as if from .mdp inputs
92 const char* dimXYZ[3] = { "x", "y", "z" };
93 GMX_RELEASE_ASSERT(dim >= 0 && dim < DIM, "Dimension should be 0, 1 or 2");
95 KeyValueTreeBuilder mdpValues;
96 mdpValues.rootObject().addValue(formatString("electric-field-%s", dimXYZ[dim]),
97 formatString("%g %g %g %g", E0, omega, t0, sigma));
99 KeyValueTreeTransformer transform;
100 transform.rules()->addRule().keyMatchType("/", StringCompareType::CaseAndDashInsensitive);
101 module->mdpOptionProvider()->initMdpTransform(transform.rules());
102 auto result = transform.transform(mdpValues.build(), nullptr);
103 Options moduleOptions;
104 module->mdpOptionProvider()->initMdpOptions(&moduleOptions);
105 assignOptionsFromKeyValueTree(&moduleOptions, result.object(), nullptr);
108 // Prepare a ForceProviderInput
109 t_mdatoms md;
110 std::vector<real> chargeA{ 1 };
111 md.homenr = ssize(chargeA);
112 md.chargeA = chargeA.data();
113 t_commrec cr;
114 matrix boxDummy = { { 0, 0, 0 }, { 0, 0, 0 }, { 0, 0, 0 } };
115 ForceProviderInput forceProviderInput({}, md, 0.0, boxDummy, cr);
117 // Prepare a ForceProviderOutput
118 PaddedVector<RVec> f = { { 0, 0, 0 } };
119 ForceWithVirial forceWithVirial(f, true);
120 gmx_enerdata_t enerdDummy(1, 0);
121 ForceProviderOutput forceProviderOutput(&forceWithVirial, &enerdDummy);
123 // Use the ForceProviders to calculate forces
124 ForceProviders forceProviders;
125 module->initForceProviders(&forceProviders);
126 forceProviders.calculateForces(forceProviderInput, &forceProviderOutput);
128 FloatingPointTolerance tolerance(relativeToleranceAsFloatingPoint(1.0, 0.005));
129 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);
148 } // namespace
149 } // namespace test
150 } // namespace gmx