Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / math / tests / invertmatrix.cpp
blob4021cb65945f4cda3d11a024b2f21825e5cca0ea
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016, 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 matrix inversion routines
39 * \todo Test error conditions when they throw exceptions
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_math
44 #include "gmxpre.h"
46 #include "gromacs/math/invertmatrix.h"
48 #include <vector>
50 #include <gtest/gtest.h>
52 #include "gromacs/math/vec.h"
54 #include "testutils/testasserts.h"
56 namespace
59 using gmx::invertMatrix;
60 using gmx::invertBoxMatrix;
61 using gmx::test::defaultRealTolerance;
63 TEST(InvertMatrixTest, IdentityIsImpotent)
65 matrix in = {{1, 0, 0}, {0, 1, 0}, {0, 0, 1}};
66 matrix out;
68 invertMatrix(in, out);
70 EXPECT_REAL_EQ_TOL(out[XX][XX], in[XX][XX], defaultRealTolerance());
71 EXPECT_REAL_EQ_TOL(out[XX][YY], in[XX][YY], defaultRealTolerance());
72 EXPECT_REAL_EQ_TOL(out[XX][ZZ], in[XX][ZZ], defaultRealTolerance());
73 EXPECT_REAL_EQ_TOL(out[YY][XX], in[YY][XX], defaultRealTolerance());
74 EXPECT_REAL_EQ_TOL(out[YY][YY], in[YY][YY], defaultRealTolerance());
75 EXPECT_REAL_EQ_TOL(out[YY][ZZ], in[YY][ZZ], defaultRealTolerance());
76 EXPECT_REAL_EQ_TOL(out[ZZ][XX], in[ZZ][XX], defaultRealTolerance());
77 EXPECT_REAL_EQ_TOL(out[ZZ][YY], in[ZZ][YY], defaultRealTolerance());
78 EXPECT_REAL_EQ_TOL(out[ZZ][ZZ], in[ZZ][ZZ], defaultRealTolerance());
81 TEST(InvertMatrixTest, ComputesInverse)
83 matrix in = {{1, 2, 3}, {-1, real(2.5), 1}, {10, -2, real(1.2)}};
84 matrix out;
85 matrix expected = {{real(-0.12019230769230768),
86 real(0.20192307692307693),
87 real(0.13221153846153844)},
88 {real(-0.26923076923076916),
89 real(0.69230769230769229),
90 real(0.096153846153846145)},
91 {real(0.55288461538461531),
92 real(-0.52884615384615374),
93 real(-0.10817307692307691)}};
95 invertMatrix(in, out);
97 EXPECT_REAL_EQ_TOL(out[XX][XX], expected[XX][XX], defaultRealTolerance());
98 EXPECT_REAL_EQ_TOL(out[XX][YY], expected[XX][YY], defaultRealTolerance());
99 EXPECT_REAL_EQ_TOL(out[XX][ZZ], expected[XX][ZZ], defaultRealTolerance());
100 EXPECT_REAL_EQ_TOL(out[YY][XX], expected[YY][XX], defaultRealTolerance());
101 EXPECT_REAL_EQ_TOL(out[YY][YY], expected[YY][YY], defaultRealTolerance());
102 EXPECT_REAL_EQ_TOL(out[YY][ZZ], expected[YY][ZZ], defaultRealTolerance());
103 EXPECT_REAL_EQ_TOL(out[ZZ][XX], expected[ZZ][XX], defaultRealTolerance());
104 EXPECT_REAL_EQ_TOL(out[ZZ][YY], expected[ZZ][YY], defaultRealTolerance());
105 EXPECT_REAL_EQ_TOL(out[ZZ][ZZ], expected[ZZ][ZZ], defaultRealTolerance());
108 TEST(InvertBoxMatrixTest, IdentityIsImpotent)
110 matrix in = {{1, 0, 0}, {0, 1, 0}, {0, 0, 1}};
112 invertBoxMatrix(in, in);
114 EXPECT_REAL_EQ_TOL(in[XX][XX], in[XX][XX], defaultRealTolerance());
115 EXPECT_REAL_EQ_TOL(in[XX][YY], in[XX][YY], defaultRealTolerance());
116 EXPECT_REAL_EQ_TOL(in[XX][ZZ], in[XX][ZZ], defaultRealTolerance());
117 EXPECT_REAL_EQ_TOL(in[YY][XX], in[YY][XX], defaultRealTolerance());
118 EXPECT_REAL_EQ_TOL(in[YY][YY], in[YY][YY], defaultRealTolerance());
119 EXPECT_REAL_EQ_TOL(in[YY][ZZ], in[YY][ZZ], defaultRealTolerance());
120 EXPECT_REAL_EQ_TOL(in[ZZ][XX], in[ZZ][XX], defaultRealTolerance());
121 EXPECT_REAL_EQ_TOL(in[ZZ][YY], in[ZZ][YY], defaultRealTolerance());
122 EXPECT_REAL_EQ_TOL(in[ZZ][ZZ], in[ZZ][ZZ], defaultRealTolerance());
125 TEST(InvertBoxMatrixTest, ComputesInverseInPlace)
127 matrix in = {{1, 0, 0}, {-1, real(2.5), 0}, {10, -2, real(1.2)}};
128 matrix expected = {{1, 0, 0},
129 {real(0.4), real(0.4), 0},
130 {real(-23.0/3.0), real(2.0/3.0), real(5.0/6.0)}};
132 invertBoxMatrix(in, in);
134 EXPECT_REAL_EQ_TOL(expected[XX][XX], in[XX][XX], defaultRealTolerance());
135 EXPECT_REAL_EQ_TOL(expected[XX][YY], in[XX][YY], defaultRealTolerance());
136 EXPECT_REAL_EQ_TOL(expected[XX][ZZ], in[XX][ZZ], defaultRealTolerance());
137 EXPECT_REAL_EQ_TOL(expected[YY][XX], in[YY][XX], defaultRealTolerance());
138 EXPECT_REAL_EQ_TOL(expected[YY][YY], in[YY][YY], defaultRealTolerance());
139 EXPECT_REAL_EQ_TOL(expected[YY][ZZ], in[YY][ZZ], defaultRealTolerance());
140 EXPECT_REAL_EQ_TOL(expected[ZZ][XX], in[ZZ][XX], defaultRealTolerance());
141 EXPECT_REAL_EQ_TOL(expected[ZZ][YY], in[ZZ][YY], defaultRealTolerance());
142 EXPECT_REAL_EQ_TOL(expected[ZZ][ZZ], in[ZZ][ZZ], defaultRealTolerance());
145 } // namespace