Skip to content
Snippets Groups Projects
Commit d4b14df2 authored by Norihiro Watanabe's avatar Norihiro Watanabe
Browse files

fix minor issues in ElementCoordinatesMappingLocal and its test

parent f472c118
No related branches found
No related tags found
No related merge requests found
...@@ -25,7 +25,7 @@ ElementCoordinatesMappingLocal::ElementCoordinatesMappingLocal( ...@@ -25,7 +25,7 @@ ElementCoordinatesMappingLocal::ElementCoordinatesMappingLocal(
: _coords(global_coords) : _coords(global_coords)
{ {
assert(e.getDimension() <= global_coords.getDimension()); assert(e.getDimension() <= global_coords.getDimension());
for(size_t i = 0; i < e.getNNodes(); i++) for(unsigned i = 0; i < e.getNNodes(); i++)
_point_vec.push_back(MeshLib::Node(*(e.getNode(i)))); _point_vec.push_back(MeshLib::Node(*(e.getNode(i))));
getRotationMatrixToGlobal(e, global_coords, _point_vec, _matR2global); getRotationMatrixToGlobal(e, global_coords, _point_vec, _matR2global);
...@@ -40,18 +40,18 @@ void ElementCoordinatesMappingLocal::rotateToLocal( ...@@ -40,18 +40,18 @@ void ElementCoordinatesMappingLocal::rotateToLocal(
std::vector<MeshLib::Node> &local_pt) const std::vector<MeshLib::Node> &local_pt) const
{ {
// rotate the point coordinates // rotate the point coordinates
const std::size_t global_dim = global_coords.getDimension(); const unsigned global_dim = global_coords.getDimension();
Eigen::VectorXd dx = Eigen::VectorXd::Zero(global_dim); Eigen::VectorXd dx = Eigen::VectorXd::Zero(global_dim);
Eigen::VectorXd x_new = Eigen::VectorXd::Zero(3); Eigen::VectorXd x_new = Eigen::VectorXd::Zero(3);
for(std::size_t i = 0; i < ele.getNNodes(); i++) for(unsigned i = 0; i < ele.getNNodes(); i++)
{ {
for (std::size_t j=0; j<global_dim; j++) for (unsigned j=0; j<global_dim; j++)
dx[j] = vec_pt[i].getCoords()[j]; dx[j] = vec_pt[i].getCoords()[j];
x_new.head(global_dim) = matR2local * dx; x_new.head(global_dim) = matR2local * dx;
local_pt[i] = MeshLib::Node(x_new.data()); local_pt[i] = MeshLib::Node(x_new.data());
} }
}; }
void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal( void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal(
const Element &e, const Element &e,
...@@ -61,7 +61,7 @@ void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal( ...@@ -61,7 +61,7 @@ void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal(
{ {
const std::size_t global_dim = global_coords.getDimension(); const std::size_t global_dim = global_coords.getDimension();
// compute R in x=R*x' where x is original coordinates and x' is local coordinates // compute R in x=R*x' where x are original coordinates and x' are local coordinates
matR = RotationMatrix::Zero(global_dim, global_dim); matR = RotationMatrix::Zero(global_dim, global_dim);
if (global_dim == e.getDimension()) { if (global_dim == e.getDimension()) {
matR = RotationMatrix::Identity(global_dim, global_dim); matR = RotationMatrix::Identity(global_dim, global_dim);
...@@ -82,7 +82,6 @@ void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal( ...@@ -82,7 +82,6 @@ void ElementCoordinatesMappingLocal::getRotationMatrixToGlobal(
MathLib::Vector3 plane_normal; MathLib::Vector3 plane_normal;
double d; double d;
GeoLib::getNewellPlane (pnts, plane_normal, d); GeoLib::getNewellPlane (pnts, plane_normal, d);
//std::cout << "pn=" << plane_normal << std::endl;
// compute a rotation matrix to XY // compute a rotation matrix to XY
MathLib::DenseMatrix<double> matToXY(3,3,0.0); MathLib::DenseMatrix<double> matToXY(3,3,0.0);
GeoLib::computeRotationMatrixToXY(plane_normal, matToXY); GeoLib::computeRotationMatrixToXY(plane_normal, matToXY);
......
...@@ -86,7 +86,7 @@ class TestQuad4 ...@@ -86,7 +86,7 @@ class TestQuad4
static MeshLib::Quad* createXYZ() static MeshLib::Quad* createXYZ()
{ {
MeshLib::Node** nodes = new MeshLib::Node*[e_nnodes]; MeshLib::Node** nodes = new MeshLib::Node*[e_nnodes];
// rotate 40 degree around x axis // rotate 45 degree around x axis
nodes[0] = new MeshLib::Node( 1.0, 0.7071067811865475, 0.7071067811865475); nodes[0] = new MeshLib::Node( 1.0, 0.7071067811865475, 0.7071067811865475);
nodes[1] = new MeshLib::Node(-1.0, 0.7071067811865475, 0.7071067811865475); nodes[1] = new MeshLib::Node(-1.0, 0.7071067811865475, 0.7071067811865475);
nodes[2] = new MeshLib::Node(-1.0, -0.7071067811865475, -0.7071067811865475); nodes[2] = new MeshLib::Node(-1.0, -0.7071067811865475, -0.7071067811865475);
...@@ -118,6 +118,7 @@ class TestQuad4 ...@@ -118,6 +118,7 @@ class TestQuad4
}; };
#ifndef NDEBUG
void debugOutput(MeshLib::Element *ele, MeshLib::ElementCoordinatesMappingLocal &mapping) void debugOutput(MeshLib::Element *ele, MeshLib::ElementCoordinatesMappingLocal &mapping)
{ {
std::cout.precision(12); std::cout.precision(12);
...@@ -136,6 +137,7 @@ void debugOutput(MeshLib::Element *ele, MeshLib::ElementCoordinatesMappingLocal ...@@ -136,6 +137,7 @@ void debugOutput(MeshLib::Element *ele, MeshLib::ElementCoordinatesMappingLocal
std::cout << (matR*v).transpose() << std::endl; std::cout << (matR*v).transpose() << std::endl;
} }
} }
#endif
// check if using the rotation matrix results in the original coordinates // check if using the rotation matrix results in the original coordinates
#define CHECK_COORDS(ele, mapping)\ #define CHECK_COORDS(ele, mapping)\
...@@ -156,7 +158,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineY) ...@@ -156,7 +158,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineY)
//debugOutput(ele, mapping); //debugOutput(ele, mapping);
double exp_R[2*2] = {0, -1, double exp_R[2*2] = {0, -1,
1, 0}; //row major 1, 0};
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
CHECK_COORDS(ele,mapping); CHECK_COORDS(ele,mapping);
...@@ -170,7 +172,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineZ) ...@@ -170,7 +172,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineZ)
auto matR(mapping.getRotationMatrixToGlobal()); auto matR(mapping.getRotationMatrixToGlobal());
//debugOutput(ele, mapping); //debugOutput(ele, mapping);
double exp_R[3*3] = {0, 0, -1, 0, 1, 0, 1, 0, 0}; //row major double exp_R[3*3] = {0, 0, -1, 0, 1, 0, 1, 0, 0};
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
CHECK_COORDS(ele,mapping); CHECK_COORDS(ele,mapping);
...@@ -184,7 +186,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineXY) ...@@ -184,7 +186,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineXY)
//debugOutput(ele, mapping); //debugOutput(ele, mapping);
double exp_R[2*2] = {0.70710678118654757, -0.70710678118654757, double exp_R[2*2] = {0.70710678118654757, -0.70710678118654757,
0.70710678118654757, 0.70710678118654757}; //row major 0.70710678118654757, 0.70710678118654757};
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
CHECK_COORDS(ele,mapping); CHECK_COORDS(ele,mapping);
...@@ -199,7 +201,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineXYZ) ...@@ -199,7 +201,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimLineXYZ)
double exp_R[3*3] = {0.57735026918962584, -0.81649658092772626, 0, double exp_R[3*3] = {0.57735026918962584, -0.81649658092772626, 0,
0.57735026918962584, 0.40824829046386313, -0.70710678118654757, 0.57735026918962584, 0.40824829046386313, -0.70710678118654757,
0.57735026918962584, 0.40824829046386313, 0.70710678118654757}; //row major 0.57735026918962584, 0.40824829046386313, 0.70710678118654757};
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
CHECK_COORDS(ele,mapping); CHECK_COORDS(ele,mapping);
...@@ -215,11 +217,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadXZ) ...@@ -215,11 +217,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadXZ)
// results when using GeoLib::ComputeRotationMatrixToXY() // results when using GeoLib::ComputeRotationMatrixToXY()
double exp_R[3*3] = { 1, 0, 0, double exp_R[3*3] = { 1, 0, 0,
0, 0, -1, 0, 0, -1,
0, 1, 0}; //row major 0, 1, 0};
// // results when using GeoLib::ComputeRotationMatrixToXY2()
// double exp_R[3*3] = { -1, 0, 0,
// 0, 0, -1,
// 0, -1, 0}; //row major
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
...@@ -236,11 +234,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadYZ) ...@@ -236,11 +234,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadYZ)
// results when using GeoLib::ComputeRotationMatrixToXY() // results when using GeoLib::ComputeRotationMatrixToXY()
double exp_R[3*3] = { 0, 0, 1, double exp_R[3*3] = { 0, 0, 1,
0, 1, 0, 0, 1, 0,
-1, 0, 0}; //row major -1, 0, 0};
// // results when using GeoLib::ComputeRotationMatrixToXY2()
// double exp_R[3*3] = { 0, 0, 1,
// -1, 0, 0,
// 0, -1, 0}; //row major
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
...@@ -257,11 +251,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadXYZ) ...@@ -257,11 +251,7 @@ TEST(MeshLib, CoordinatesMappingLocalLowerDimQuadXYZ)
// results when using GeoLib::ComputeRotationMatrixToXY() // results when using GeoLib::ComputeRotationMatrixToXY()
double exp_R[3*3] = { 1, 0, 0, double exp_R[3*3] = { 1, 0, 0,
0, 0.70710678118654757, -0.70710678118654757, 0, 0.70710678118654757, -0.70710678118654757,
0, 0.70710678118654757, 0.70710678118654757}; //row major 0, 0.70710678118654757, 0.70710678118654757};
// // results when using GeoLib::ComputeRotationMatrixToXY2()
// double exp_R[3*3] = { -1, 0, 0,
// 0, -0.70710678118654757, -0.70710678118654757,
// 0, -0.70710678118654757, 0.70710678118654757}; //row major
const double eps(std::numeric_limits<double>::epsilon()); const double eps(std::numeric_limits<double>::epsilon());
ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps); ASSERT_ARRAY_NEAR(exp_R, matR.data(), matR.size(), eps);
......
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