Skip to content
Snippets Groups Projects
Commit 86ac5192 authored by Dmitri Naumov's avatar Dmitri Naumov
Browse files

[T] LIE: Add tests for rotation matrix.

The tests are capturing the current behaviour with OGS-surface normals pointing in wrong direction.
parent 6b4eb1f9
No related branches found
No related tags found
No related merge requests found
...@@ -15,6 +15,7 @@ ...@@ -15,6 +15,7 @@
#include <Eigen/Eigen> #include <Eigen/Eigen>
#include "MeshLib/Elements/Line.h" #include "MeshLib/Elements/Line.h"
#include "MeshLib/Elements/Tri.h"
#include "MeshLib/Mesh.h" #include "MeshLib/Mesh.h"
#include "ProcessLib/LIE/Common/Utils.h" #include "ProcessLib/LIE/Common/Utils.h"
...@@ -22,6 +23,19 @@ ...@@ -22,6 +23,19 @@
namespace namespace
{ {
std::unique_ptr<MeshLib::Mesh> createTriangle(
std::array<std::array<double, 3>, 3> const& points)
{
MeshLib::Node** nodes = new MeshLib::Node*[3];
for (int i = 0; i < points.size(); ++i)
nodes[i] = new MeshLib::Node(points[i]);
MeshLib::Element* e = new MeshLib::Tri(nodes);
return std::unique_ptr<MeshLib::Mesh>(new MeshLib::Mesh(
"",
std::vector<MeshLib::Node*>{nodes[0], nodes[1], nodes[2]},
std::vector<MeshLib::Element*>{e}));
}
std::unique_ptr<MeshLib::Mesh> createLine( std::unique_ptr<MeshLib::Mesh> createLine(
std::array<double, 3> const& a, std::array<double, 3> const& b) std::array<double, 3> const& a, std::array<double, 3> const& b)
...@@ -58,6 +72,56 @@ const double eps = std::numeric_limits<double>::epsilon(); ...@@ -58,6 +72,56 @@ const double eps = std::numeric_limits<double>::epsilon();
} }
TEST(LIE, rotationMatrixXYTriangle)
{
auto msh = createTriangle(
{{{{0.0, 0.0, 0.0}}, {{1.0, 0.0, 0.0}}, {{1.0, 1.0, 0.0}}}});
auto e(msh->getElement(0));
Eigen::Vector3d nv;
ProcessLib::LIE::computeNormalVector(*e, 3, nv);
ASSERT_EQ(0., nv[0]);
ASSERT_EQ(0., nv[1]);
ASSERT_EQ(-1., nv[2]);
Eigen::MatrixXd R(3, 3);
ProcessLib::LIE::computeRotationMatrix(*e, nv, 3, R);
ASSERT_NEAR(-1., R(0, 0), eps);
ASSERT_NEAR(0., R(0, 1), eps);
ASSERT_NEAR(0., R(0, 2), eps);
ASSERT_NEAR(0., R(1, 0), eps);
ASSERT_NEAR(1., R(1, 1), eps);
ASSERT_NEAR(0., R(1, 2), eps);
ASSERT_NEAR(0., R(2, 0), eps);
ASSERT_NEAR(0., R(2, 1), eps);
ASSERT_NEAR(-1., R(2, 2), eps);
}
TEST(LIE, rotationMatrixYZTriangle)
{
auto msh = createTriangle(
{{{{0.0, 0.0, 0.0}}, {{0.0, 1.0, 0.0}}, {{0.0, 1.0, 1.0}}}});
auto e(msh->getElement(0));
Eigen::Vector3d nv;
ProcessLib::LIE::computeNormalVector(*e, 3, nv);
ASSERT_EQ(-1., nv[0]);
ASSERT_EQ(0., nv[1]);
ASSERT_EQ(0., nv[2]);
Eigen::MatrixXd R(3, 3);
ProcessLib::LIE::computeRotationMatrix(*e, nv, 3, R);
ASSERT_NEAR(0., R(0, 0), eps);
ASSERT_NEAR(-1., R(0, 1), eps);
ASSERT_NEAR(0., R(0, 2), eps);
ASSERT_NEAR(0., R(1, 0), eps);
ASSERT_NEAR(0., R(1, 1), eps);
ASSERT_NEAR(1., R(1, 2), eps);
ASSERT_NEAR(-1., R(2, 0), eps);
ASSERT_NEAR(0., R(2, 1), eps);
ASSERT_NEAR(0., R(2, 2), eps);
}
TEST(LIE, rotationMatrixX) TEST(LIE, rotationMatrixX)
{ {
auto msh(createX()); auto msh(createX());
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment