diff --git a/src/Classic/Fields/Interpolation/NDGrid.cpp b/src/Classic/Fields/Interpolation/NDGrid.cpp index 59d551caa7cd24f1e785080607425b7f4da8e154..be17aca95214b6aa10b20c599ac37fb3696ab216 100644 --- a/src/Classic/Fields/Interpolation/NDGrid.cpp +++ b/src/Classic/Fields/Interpolation/NDGrid.cpp @@ -1,4 +1,30 @@ -// MAUS WARNING: THIS IS LEGACY CODE. +/* + * Copyright (c) 2017, Chris Rogers + * All rights reserved. + * Redistribution and use in source and binary forms, with or without + * modification, are permitted provided that the following conditions are met: + * 1. Redistributions of source code must retain the above copyright notice, + * this list of conditions and the following disclaimer. + * 2. Redistributions in binary form must reproduce the above copyright notice, + * this list of conditions and the following disclaimer in the documentation + * and/or other materials provided with the distribution. + * 3. Neither the name of STFC nor the names of its contributors may be used to + * endorse or promote products derived from this software without specific + * prior written permission. + * + * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" + * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE + * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE + * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE + * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR + * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF + * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS + * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN + * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) + * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE + * POSSIBILITY OF SUCH DAMAGE. + */ + #include <math.h> #include <iomanip> #include "Utilities/OpalException.h" @@ -10,8 +36,6 @@ namespace interpolation { NDGrid::NDGrid() : coord_m(), maps_m(), constantSpacing_m(false) { } - - NDGrid::NDGrid(std::vector<int> size, std::vector<const double *> gridCoordinates) : coord_m(), maps_m(), constantSpacing_m(false) { for (unsigned int i=0; i<size.size(); i++) { @@ -141,7 +165,7 @@ Mesh::Iterator& NDGrid::subOne(Mesh::Iterator& lhs) const { return lhs; } -void NDGrid::setConstantSpacing() { +void NDGrid::setConstantSpacing(double tolerance_m) { constantSpacing_m = true; for (unsigned int i = 0; i < coord_m.size(); i++) { for (unsigned int j = 0; j < coord_m[i].size()-1; j++) { diff --git a/src/Classic/Fields/Interpolation/NDGrid.h b/src/Classic/Fields/Interpolation/NDGrid.h index 4900a93d9778ad2c339deeae5ba1797f241e6b87..69434b19a8e68d6e5be67fbed4299cf8c13d938d 100644 --- a/src/Classic/Fields/Interpolation/NDGrid.h +++ b/src/Classic/Fields/Interpolation/NDGrid.h @@ -172,8 +172,8 @@ class NDGrid : public Mesh { * @param xIndex: will be set with the resultant grid index * * Sets xIndex to the index of the largest mesh position that is less than - * x in the given dimension. - * + * x in the given dimension. If x is less than all mesh positions, sets + * xIndex to -1. */ inline void coordLowerBound(const double& x, const int& dimension, int& xIndex) const; @@ -187,7 +187,7 @@ class NDGrid : public Mesh { * Sets xIndex to the index of the largest mesh position that is less than * pos in all dimension. * - * Note that pos and xIndex are not bound checked. + * Note that dimension and xIndex are not bound checked. */ inline void lowerBound(const std::vector<double>& pos, std::vector<int>& xIndex) const; @@ -199,7 +199,7 @@ class NDGrid : public Mesh { /** Reset the coordinates in the dimension to a new value * - * @param dimension: the dimension of the coordates that will be reset + * @param dimension: the dimension of the coordinates that will be reset * @param nCoords: the length of array x * @param x: array of points. Caller owns the memory allocated to x. */ @@ -242,7 +242,7 @@ class NDGrid : public Mesh { * Loops over all grid points and checks the step size is regular, with * tolerance given by tolernance_m */ - void setConstantSpacing(); + void setConstantSpacing(double tolerance_m = 1e-9); /** Return an integer corresponding to the iterator position * @@ -278,7 +278,6 @@ private: std::vector< std::vector<double> > coord_m; std::vector<VectorMap*> maps_m; bool constantSpacing_m; - const double tolerance_m = 1e-9; friend Mesh::Iterator operator++(Mesh::Iterator& lhs, int); friend Mesh::Iterator operator--(Mesh::Iterator& lhs, int); @@ -314,10 +313,16 @@ std::vector<double> NDGrid::coordVector(const int& dimension) const { } void NDGrid::coordLowerBound(const double& x, const int& dimension, int& xIndex) const { - if(constantSpacing_m) { + if (constantSpacing_m) { double x0 = coord_m[dimension][0]; double x1 = coord_m[dimension][1]; - xIndex = static_cast<int>(floor((x - x0)/(x1-x0))); + xIndex = static_cast<int>(floor((x - x0)/(x1-x0))); + int coordSize = static_cast<int>(coord_m[dimension].size()); + if (xIndex >= coordSize) { + xIndex = coordSize-1; + } else if (xIndex < -1) { + xIndex = -1; + } } else { const std::vector<double>& c_t(coord_m[dimension]); xIndex = std::lower_bound(c_t.begin(), c_t.end(), x)-c_t.begin()-1; @@ -325,8 +330,8 @@ void NDGrid::coordLowerBound(const double& x, const int& dimension, int& xIndex) } void NDGrid::lowerBound(const std::vector<double>& pos, std::vector<int>& xIndex) const { - xIndex = std::vector<int>(pos.size()); - for(unsigned int i=0; i<pos.size(); i++) { + xIndex = std::vector<int> (pos.size()); + for(unsigned int i=0; i < pos.size(); i++) { coordLowerBound(pos[i], i, xIndex[i]); } } diff --git a/tests/classic_src/CMakeLists.txt b/tests/classic_src/CMakeLists.txt index fd92d0937f568075f212040489896d4bc56ece8d..435aed0a0795f2f79334f6f97b424cacc241a2f9 100644 --- a/tests/classic_src/CMakeLists.txt +++ b/tests/classic_src/CMakeLists.txt @@ -1,5 +1,7 @@ add_subdirectory (AbsBeamline) add_subdirectory (Algorithms) add_subdirectory (Utilities) +add_subdirectory (Fields) -set (TEST_SRCS_LOCAL ${TEST_SRCS_LOCAL} PARENT_SCOPE) \ No newline at end of file + +set (TEST_SRCS_LOCAL ${TEST_SRCS_LOCAL} PARENT_SCOPE) diff --git a/tests/classic_src/Fields/Interpolation/NDGridTest.cpp b/tests/classic_src/Fields/Interpolation/NDGridTest.cpp index 435bc633d80f8d80fe30a708e98338159be312d1..1f2fa8d37c4e1f9ea276d59005c79223e4264c4b 100644 --- a/tests/classic_src/Fields/Interpolation/NDGridTest.cpp +++ b/tests/classic_src/Fields/Interpolation/NDGridTest.cpp @@ -27,8 +27,10 @@ #include "gtest/gtest.h" #include "Fields/Interpolation/NDGrid.h" +#include "Fields/Interpolation/Mesh.h" + +namespace ndgridtest { -// namespace ndgridtest { class NDGridTest : public ::testing::Test { public: NDGridTest() : grid_m(NULL) { @@ -42,18 +44,24 @@ public: gridCoordinates[1][1] = 5.; gridCoordinates[1][2] = 9.; grid_m = new interpolation::NDGrid(gridCoordinates); + gridCoordinates[1][2] = 10.; // force it to not be regular + grid2_m = new interpolation::NDGrid(gridCoordinates); } void TearDown( ) { delete grid_m; grid_m = NULL; + delete grid2_m; + grid2_m = NULL; } ~NDGridTest() { } -private: interpolation::NDGrid* grid_m; + interpolation::NDGrid* grid2_m; + +private: }; @@ -140,5 +148,160 @@ TEST_F(NDGridTest, CoordVectorTest) { // and newCoordArray delete coords_a; } } -//} // namespace ndgridtest + +TEST_F(NDGridTest, CoordLowerBoundTest) { + // first dimension ... 0., 3.; + int index = -1; + grid_m->coordLowerBound(-1., 0, index); + EXPECT_EQ(index, -1); + index = -2; + grid_m->coordLowerBound(0.0001, 0, index); + EXPECT_EQ(index, 0); + index = -2; + grid_m->coordLowerBound(2.999, 0, index); + EXPECT_EQ(index, 0); + index = -2; + grid_m->coordLowerBound(3.001, 0, index); + EXPECT_EQ(index, 1); + index = -2; + grid_m->coordLowerBound(1000.0001, 0, index); + EXPECT_EQ(index, 1); + + // second dimension ... 1., 5., 9. + index = -2; + grid_m->coordLowerBound(0.999, 1, index); + EXPECT_EQ(index, -1); + index = -2; + grid_m->coordLowerBound(4.999, 1, index); + EXPECT_EQ(index, 0); + index = -2; + grid_m->coordLowerBound(8.999, 1, index); + EXPECT_EQ(index, 1); + index = -2; + grid_m->coordLowerBound(9.0001, 1, index); + EXPECT_EQ(index, 2); + index = -2; + grid_m->coordLowerBound(1000.0001, 1, index); + EXPECT_EQ(index, 2); +} + +TEST_F(NDGridTest, LowerBoundTest) { + // first dimension ... 0., 3.; + // second dimension ... 1., 5., 9. + std::vector<int> index1(2, -2); + std::vector<double> pos1(2, -1.); + grid_m->lowerBound(pos1, index1); + EXPECT_EQ(index1[0], -1); + EXPECT_EQ(index1[1], -1); + + std::vector<int> index2(2, -2); + std::vector<double> pos2(2, 4.); + grid_m->lowerBound(pos2, index2); + EXPECT_EQ(index2[0], 1); + EXPECT_EQ(index2[1], 0); + + std::vector<int> index3(2, -2); + std::vector<double> pos3(2, 100.); + grid_m->lowerBound(pos3, index3); + EXPECT_EQ(index3[0], 1); + EXPECT_EQ(index3[1], 2); +} + +TEST_F(NDGridTest, MinMaxTest) { + EXPECT_EQ(grid_m->min(0), 0.); + EXPECT_EQ(grid_m->min(1), 1.); + EXPECT_EQ(grid_m->max(0), 3.); + EXPECT_EQ(grid_m->max(1), 9.); +} + +TEST_F(NDGridTest, SetCoordTest) { + double xNew[] = {5., 10., 12., 15.}; + grid_m->setCoord(0, 4, xNew); + std::vector<double> xTest = grid_m->coordVector(0); + EXPECT_EQ(xTest.size(), 4); + for (size_t i = 0; i < 4; ++i) { + EXPECT_EQ(xTest[i], xNew[i]); + } +} + +TEST_F(NDGridTest, BeginEndTest) { + ASSERT_EQ(grid_m->begin().getState().size(), 2); + EXPECT_EQ(grid_m->begin().getState()[0], 1); + EXPECT_EQ(grid_m->begin().getState()[1], 1); + ASSERT_EQ(grid_m->end().getState().size(), 2); + EXPECT_EQ(grid_m->end().getState()[0], 3); // one past the last (2, 3) + EXPECT_EQ(grid_m->end().getState()[1], 1); +} + +TEST_F(NDGridTest, GetPositionTest) { + std::vector<double> position(3, -1); + interpolation::Mesh::Iterator it = grid_m->begin(); + grid_m->getPosition(it, &position[0]); + EXPECT_EQ(grid_m->getPositionDimension(), 2); + EXPECT_EQ(position[0], 0.); + EXPECT_EQ(position[1], 1.); + it[0] = 2; + it[1] = 3; + grid_m->getPosition(it, &position[0]); + EXPECT_EQ(position[0], 3.); + EXPECT_EQ(position[1], 9.); +} + +TEST_F(NDGridTest, GetSetConstantSpacingTest) { + EXPECT_TRUE(grid_m->getConstantSpacing()); + grid_m->setConstantSpacing(false); + EXPECT_FALSE(grid_m->getConstantSpacing()); + grid_m->setConstantSpacing(true); + EXPECT_TRUE(grid_m->getConstantSpacing()); + grid_m->setConstantSpacing(false); + grid_m->setConstantSpacing(); + EXPECT_TRUE(grid_m->getConstantSpacing()); + + EXPECT_FALSE(grid2_m->getConstantSpacing()); + grid2_m->setConstantSpacing(true); + EXPECT_TRUE(grid2_m->getConstantSpacing()); + grid2_m->setConstantSpacing(false); + EXPECT_FALSE(grid2_m->getConstantSpacing()); + grid2_m->setConstantSpacing(true); + grid2_m->setConstantSpacing(0.1); + EXPECT_FALSE(grid2_m->getConstantSpacing()); + grid2_m->setConstantSpacing(2.01); + EXPECT_TRUE(grid2_m->getConstantSpacing()); +} + +TEST_F(NDGridTest, ToIntegerTest) { + interpolation::Mesh::Iterator it = grid_m->begin(); + EXPECT_EQ(grid_m->toInteger(it), 0); + it[0] = 2; + it[1] = 3; + EXPECT_EQ(grid_m->toInteger(it), 5); +} + +TEST_F(NDGridTest, GetNearestTest) { + // first dimension ... 0., 3.; + // second dimension ... 1., 5., 9. + + std::vector<double> pos(2, 0.); + interpolation::Mesh::Iterator it; + pos[0] = 1.49; + pos[1] = 2.99; + it = grid_m->getNearest(&pos[0]); + EXPECT_EQ(it[0], 1); + EXPECT_EQ(it[1], 1); + pos[0] = 0.49; + it = grid_m->getNearest(&pos[0]); + EXPECT_EQ(it[0], 1); + EXPECT_EQ(it[1], 1); + pos[1] = 9.49; + it = grid_m->getNearest(&pos[0]); + EXPECT_EQ(it[0], 1); + EXPECT_EQ(it[1], 3); + + pos[0] = 3.49; + it = grid_m->getNearest(&pos[0]); + EXPECT_EQ(it[0], 2); + EXPECT_EQ(it[1], 3); +} + +} // namespace ndgridtest diff --git a/tests/classic_src/Utilities/RingSectionTest.cpp b/tests/classic_src/Utilities/RingSectionTest.cpp index 9e033df30a29c6df40e7352aaa08240b326b3013..028d89e0410634dd50651e9af76838e034f7bf66 100644 --- a/tests/classic_src/Utilities/RingSectionTest.cpp +++ b/tests/classic_src/Utilities/RingSectionTest.cpp @@ -43,7 +43,9 @@ TEST(RingSectionTest, TestConstructDestruct) { EXPECT_EQ(ors.getEndNormal(), vec0); EXPECT_EQ(ors.getComponentPosition(), vec0); EXPECT_EQ(ors.getComponentOrientation(), vec0); - EXPECT_EQ(ors.getVirtualBoundingBox(), std::vector<Vector_t>(4, vec0)); + for (size_t i = 0; i < 4; ++i) + for (size_t j = 0; j < 3; ++j) + EXPECT_EQ(ors.getVirtualBoundingBox().at(i)[j], 0.); RingSection ors_comp; MockComponent comp;