Make a shared header of 3D transformations

Affine transformations of homogeneous coordinates using 4x4 matrices are
quite common in visualization. Create a new math header file in the base
vtkm namespace that has common functions for such coordinates.

Much of this implementation was taken from the rendering matrix helpers.
This commit is contained in:
Kenneth Moreland 2016-06-08 12:39:47 -06:00
parent d55402e998
commit cdeeda67bb
7 changed files with 448 additions and 136 deletions

@ -36,6 +36,7 @@ set(headers
Range.h
StaticAssert.h
TopologyElementTag.h
Transform3D.h
TypeListTag.h
Types.h
TypeTraits.h

216
vtkm/Transform3D.h Normal file

@ -0,0 +1,216 @@
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2016 Sandia Corporation.
// Copyright 2016 UT-Battelle, LLC.
// Copyright 2016 Los Alamos National Security.
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
// Under the terms of Contract DE-AC52-06NA25396 with Los Alamos National
// Laboratory (LANL), the U.S. Government retains certain rights in
// this software.
//
//=============================================================================
#ifndef vtk_m_Transform3D_h
#define vtk_m_Transform3D_h
// This header file contains a collection of math functions useful in the
// linear transformation of homogeneous points for rendering in 3D.
#include <vtkm/Matrix.h>
#include <vtkm/VectorAnalysis.h>
namespace vtkm {
/// \brief Transform a 3D point by a transformation matrix.
///
/// Given a 4x4 transformation matrix and a 3D point, returns the point
/// transformed by the given matrix in homogeneous coordinates.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Vec<T,3> Transform3DPoint(const vtkm::Matrix<T,4,4> &matrix,
const vtkm::Vec<T,3> &point)
{
vtkm::Vec<T,4> homogeneousPoint(point[0], point[1], point[2], T(1));
homogeneousPoint = vtkm::MatrixMultiply(matrix, homogeneousPoint);
T invW = T(1)/homogeneousPoint[3];
return vtkm::Vec<T,3>(
homogeneousPoint[0]*invW,
homogeneousPoint[1]*invW,
homogeneousPoint[2]*invW);
}
/// \brief Transform a 3D point by a transformation matrix.
///
/// Given a 4x4 transformation matrix and a 3D point, returns the point
/// transformed by the given matrix in homogeneous coordinates.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Vec<T,3> Transform3DVector(const vtkm::Matrix<T,4,4> &matrix,
const vtkm::Vec<T,3> &vector)
{
vtkm::Vec<T,4> homogeneousVector(vector[0], vector[1], vector[2], T(0));
homogeneousVector = vtkm::MatrixMultiply(matrix, homogeneousVector);
return vtkm::Vec<T,3>(
homogeneousVector[0],
homogeneousVector[1],
homogeneousVector[2]);
}
/// \brief Returns a scale matrix.
///
/// Given a scale factor for the x, y, and z directions, returns a
/// transformation matrix for those scales.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4>
Transform3DScale(const T &scaleX, const T &scaleY, const T &scaleZ)
{
vtkm::Matrix<T,4,4> scaleMatrix(T(0));
scaleMatrix(0,0) = scaleX;
scaleMatrix(1,1) = scaleY;
scaleMatrix(2,2) = scaleZ;
scaleMatrix(3,3) = T(1);
return scaleMatrix;
}
/// \brief Returns a scale matrix.
///
/// Given a scale factor for the x, y, and z directions (defined in a Vec),
/// returns a transformation matrix for those scales.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DScale(const vtkm::Vec<T,3> &scaleVec)
{
return vtkm::Transform3DScale(scaleVec[0], scaleVec[1], scaleVec[2]);
}
/// \brief Returns a scale matrix.
///
/// Given a uniform scale factor, returns a transformation matrix for those
/// scales.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DScale(const T &scale)
{
return vtkm::Transform3DScale(scale, scale, scale);
}
/// \brief Returns a translation matrix.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DTranslate(const T &x, const T &y, const T &z)
{
vtkm::Matrix<T,4,4> translateMatrix;
vtkm::MatrixIdentity(translateMatrix);
translateMatrix(0,3) = x;
translateMatrix(1,3) = y;
translateMatrix(2,3) = z;
return translateMatrix;
}
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DTranslate(const vtkm::Vec<T,3> &v)
{
return vtkm::Transform3DTranslate(v[0], v[1], v[2]);
}
/// \brief Returns a rotation matrix.
///
/// Given an angle (in radians) and an axis of rotation, returns a
/// transformation matrix that rotates around the given axis. The rotation
/// follows the right-hand rule, so if the vector points toward the user, the
/// rotation will be counterclockwise.
///
/// Note that, unlike with OpenGL, the angle is given in radians, not degrees.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DRotate(T angleRadians,
const vtkm::Vec<T,3> &axisOfRotation)
{
const vtkm::Vec<T,3> normAxis = vtkm::Normal(axisOfRotation);
T sinAngle = vtkm::Sin(angleRadians);
T cosAngle = vtkm::Cos(angleRadians);
vtkm::Matrix<T,4,4> matrix;
matrix(0,0) = normAxis[0]*normAxis[0]*(1-cosAngle) + cosAngle;
matrix(0,1) = normAxis[0]*normAxis[1]*(1-cosAngle) - normAxis[2]*sinAngle;
matrix(0,2) = normAxis[0]*normAxis[2]*(1-cosAngle) + normAxis[1]*sinAngle;
matrix(0,3) = T(0);
matrix(1,0) = normAxis[1]*normAxis[0]*(1-cosAngle) + normAxis[2]*sinAngle;
matrix(1,1) = normAxis[1]*normAxis[1]*(1-cosAngle) + cosAngle;
matrix(1,2) = normAxis[1]*normAxis[2]*(1-cosAngle) - normAxis[0]*sinAngle;
matrix(1,3) = T(0);
matrix(2,0) = normAxis[2]*normAxis[0]*(1-cosAngle) - normAxis[1]*sinAngle;
matrix(2,1) = normAxis[2]*normAxis[1]*(1-cosAngle) + normAxis[0]*sinAngle;
matrix(2,2) = normAxis[2]*normAxis[2]*(1-cosAngle) + cosAngle;
matrix(2,3) = T(0);
matrix(3,0) = T(0);
matrix(3,1) = T(0);
matrix(3,2) = T(0);
matrix(3,3) = T(1);
return matrix;
}
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DRotate(T angleRadians, T x, T y, T z)
{
return vtkm::Transform3DRotate(angleRadians, vtkm::Vec<T,3>(x,y,z));
}
/// \brief Returns a rotation matrix.
///
/// Returns a transformation matrix that rotates around the x axis.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DRotateX(T angleRadians)
{
return vtkm::Transform3DRotate(angleRadians, T(1), T(0), T(0));
}
/// \brief Returns a rotation matrix.
///
/// Returns a transformation matrix that rotates around the y axis.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DRotateY(T angleRadians)
{
return vtkm::Transform3DRotate(angleRadians, T(0), T(1), T(0));
}
/// \brief Returns a rotation matrix.
///
/// Returns a transformation matrix that rotates around the z axis.
///
template<typename T>
VTKM_EXEC_CONT_EXPORT
vtkm::Matrix<T,4,4> Transform3DRotateZ(T angleRadians)
{
return vtkm::Transform3DRotate(angleRadians, T(0), T(0), T(1));
}
} // namespace vtkm
#endif //vtk_m_Transform3D_h

@ -22,6 +22,7 @@
#include <vtkm/Bounds.h>
#include <vtkm/Math.h>
#include <vtkm/Matrix.h>
#include <vtkm/Transform3D.h>
#include <vtkm/Range.h>
#include <vtkm/VectorAnalysis.h>
#include <vtkm/rendering/MatrixHelpers.h>
@ -79,8 +80,8 @@ class Camera
matrix(3,3) = 0.f;
vtkm::Matrix<vtkm::Float32,4,4> T, Z;
T = MatrixHelpers::TranslateMatrix(this->XPan, this->YPan, 0);
Z = MatrixHelpers::ScaleMatrix(this->Zoom, this->Zoom, 1);
T = vtkm::Transform3DTranslate(this->XPan, this->YPan, 0.f);
Z = vtkm::Transform3DScale(this->Zoom, this->Zoom, 1.f);
matrix = vtkm::MatrixMultiply(Z, vtkm::MatrixMultiply(T, matrix));
return matrix;
}
@ -249,18 +250,6 @@ public:
}
}
VTKM_CONT_EXPORT
vtkm::Vec<vtkm::Float32, 3>
MultVector(const vtkm::Matrix<vtkm::Float32,4,4> &matrix, vtkm::Vec<vtkm::Float32, 3> &v) const
{
vtkm::Vec<vtkm::Float32,4> v4(v[0],v[1],v[2], 1);
v4 = vtkm::MatrixMultiply(matrix, v4);
v[0] = v4[0];
v[1] = v4[1];
v[2] = v4[2];
return v;
}
/// \brief The mode of the camera (2D or 3D).
///
/// \c vtkm::Camera can be set to a 2D or 3D mode. 2D mode is used for
@ -496,11 +485,11 @@ public:
//Translate matrix
vtkm::Matrix<vtkm::Float32,4,4> translate =
MatrixHelpers::TranslateMatrix(-this->Camera3D.LookAt);
vtkm::Transform3DTranslate(-this->Camera3D.LookAt);
//Translate matrix
vtkm::Matrix<vtkm::Float32,4,4> inverseTranslate =
MatrixHelpers::TranslateMatrix(this->Camera3D.LookAt);
vtkm::Transform3DTranslate(this->Camera3D.LookAt);
vtkm::Matrix<vtkm::Float32,4,4> view = this->CreateViewMatrix();
view(0,3) = 0;
@ -516,9 +505,12 @@ public:
inverseView, vtkm::MatrixMultiply(
rotate, vtkm::MatrixMultiply(
view,translate))));
this->Camera3D.Position = MultVector(fullTransform, this->Camera3D.Position);
this->Camera3D.LookAt = MultVector(fullTransform, this->Camera3D.LookAt);
this->Camera3D.ViewUp = MultVector(fullTransform, this->Camera3D.ViewUp);
this->Camera3D.Position =
vtkm::Transform3DPoint(fullTransform, this->Camera3D.Position);
this->Camera3D.LookAt =
vtkm::Transform3DPoint(fullTransform, this->Camera3D.LookAt);
this->Camera3D.ViewUp =
vtkm::Transform3DVector(fullTransform, this->Camera3D.ViewUp);
}
/// \brief Set up the camera to look at geometry

@ -109,119 +109,6 @@ struct MatrixHelpers
return matrix;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> ScaleMatrix(const vtkm::Vec<vtkm::Float32,3> &v)
{
return ScaleMatrix(v[0], v[1], v[2]);
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> ScaleMatrix(const vtkm::Float32 &s)
{
return ScaleMatrix(s,s,s);
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> ScaleMatrix(const vtkm::Float32 &x,
const vtkm::Float32 &y,
const vtkm::Float32 &z)
{
vtkm::Matrix<vtkm::Float32,4,4> scaleMatrix(0.0f);
scaleMatrix(0,0) = x;
scaleMatrix(1,1) = y;
scaleMatrix(2,2) = z;
scaleMatrix(3,3) = 1.0f;
return scaleMatrix;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> TranslateMatrix(const vtkm::Vec<vtkm::Float32,3> &v)
{
return TranslateMatrix(v[0], v[1], v[2]);
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> TranslateMatrix(const vtkm::Float32 &x,
const vtkm::Float32 &y,
const vtkm::Float32 &z)
{
vtkm::Matrix<vtkm::Float32,4,4> translateMatrix;
vtkm::MatrixIdentity(translateMatrix);
translateMatrix(0,3) = x;
translateMatrix(1,3) = y;
translateMatrix(2,3) = z;
return translateMatrix;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> RotateXMatrix(vtkm::Float32 angleRadians)
{
vtkm::Matrix<vtkm::Float32,4,4> M;
M(0,0) = 1.f;
M(0,1) = 0.f;
M(0,2) = 0.f;
M(0,3) = 0.f;
M(1,0) = 0.f;
M(1,1) = cosf(angleRadians);
M(1,2) = - sinf(angleRadians);
M(1,3) = 0.f;
M(2,0) = 0.f;
M(2,1) = sinf(angleRadians);
M(2,2) = cosf(angleRadians);
M(2,3) = 0.f;
M(3,0) = 0.f;
M(3,1) = 0.f;
M(3,2) = 0.f;
M(3,3) = 1.f;
return M;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> RotateYMatrix(vtkm::Float32 angleRadians)
{
vtkm::Matrix<vtkm::Float32,4,4> M;
M(0,0) = cosf(angleRadians);
M(0,1) = 0.f;
M(0,2) = sinf(angleRadians);
M(0,3) = 0.f;
M(1,0) = 0.f;
M(1,1) = 1.f;
M(1,2) = 0.f;
M(1,3) = 0.f;
M(2,0) = - sinf(angleRadians);
M(2,1) = 0.f;
M(2,2) = cosf(angleRadians);
M(2,3) = 0.f;
M(3,0) = 0.f;
M(3,1) = 0.f;
M(3,2) = 0.f;
M(3,3) = 1.f;
return M;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> RotateZMatrix(vtkm::Float32 angleRadians)
{
vtkm::Matrix<vtkm::Float32,4,4> M;
M(0,0) = cosf(angleRadians);
M(0,1) = - sinf(angleRadians);
M(0,2) = 0.f;
M(0,3) = 0.f;
M(1,0) = sinf(angleRadians);
M(1,1) = cosf(angleRadians);
M(1,2) = 0.f;
M(1,3) = 0.f;
M(2,0) = 0.f;
M(2,1) = 0.f;
M(2,2) = 1.f;
M(2,3) = 0.f;
M(3,0) = 0.f;
M(3,1) = 0.f;
M(3,2) = 0.f;
M(3,3) = 1.f;
return M;
}
static VTKM_CONT_EXPORT
vtkm::Matrix<vtkm::Float32,4,4> TrackballMatrix(vtkm::Float32 p1x,
vtkm::Float32 p1y,

@ -182,13 +182,13 @@ public:
vtkm::Matrix<vtkm::Float32, 4, 4> T;
T = MatrixHelpers::TranslateMatrix(psx,psy,-psz);
T = vtkm::Transform3DTranslate(psx,psy,-psz);
vtkm::Float32 WindowAspect =
vtkm::Float32(canvas.GetWidth()) / vtkm::Float32(canvas.GetHeight());
vtkm::Matrix<vtkm::Float32, 4, 4> SW;
SW = MatrixHelpers::ScaleMatrix(1.f/WindowAspect, 1, 1);
SW = vtkm::Transform3DScale(1.f/WindowAspect, 1.f, 1.f);
vtkm::Matrix<vtkm::Float32, 4, 4> SV;
vtkm::MatrixIdentity(SV);
@ -198,11 +198,11 @@ public:
camera.GetRealViewport(canvas.GetWidth(),canvas.GetHeight(),vl,vr,vb,vt);
vtkm::Float32 xs = (vr-vl);
vtkm::Float32 ys = (vt-vb);
SV = MatrixHelpers::ScaleMatrix(2.f/xs, 2.f/ys, 1);
SV = vtkm::Transform3DScale(2.f/xs, 2.f/ys, 1.f);
}
vtkm::Matrix<vtkm::Float32, 4, 4> R;
R = MatrixHelpers::RotateZMatrix(Angle * 3.14159265f / 180.f);
R = vtkm::Transform3DRotateZ(Angle * 3.14159265f / 180.f);
vtkm::Vec<vtkm::Float32,4> origin4(0,0,0,1);
vtkm::Vec<vtkm::Float32,4> right4(1,0,0,0);

@ -41,6 +41,7 @@ set(unit_tests
UnitTestPair.cxx
UnitTestRange.cxx
UnitTestTesting.cxx
UnitTestTransform3D.cxx
UnitTestTypeListTag.cxx
UnitTestTypes.cxx
UnitTestTypeTraits.cxx

@ -0,0 +1,215 @@
//=============================================================================
//
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
//
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//
// Copyright 2016 Sandia Corporation.
// Copyright 2016 UT-Battelle, LLC.
// Copyright 2016 Los Alamos National Security.
//
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
// the U.S. Government retains certain rights in this software.
// Under the terms of Contract DE-AC52-06NA25396 with Los Alamos National
// Laboratory (LANL), the U.S. Government retains certain rights in
// this software.
//
//=============================================================================
#include <vtkm/Transform3D.h>
#include <vtkm/testing/Testing.h>
VTKM_THIRDPARTY_PRE_INCLUDE
#include <boost/random/mersenne_twister.hpp>
#include <boost/random/uniform_real_distribution.hpp>
VTKM_THIRDPARTY_POST_INCLUDE
namespace {
boost::mt19937 g_RandomGenerator;
template<typename T>
struct TransformTests
{
boost::random::uniform_real_distribution<T> RandomDistribution;
TransformTests()
: RandomDistribution(0.0f, 1.0f) { }
T RandomNum() const { return this->RandomDistribution(g_RandomGenerator); }
typedef vtkm::Vec<T,3> Vec;
typedef vtkm::Matrix<T,4,4> Transform;
Vec RandomVector() const
{
Vec vec(this->RandomNum(), this->RandomNum(), this->RandomNum());
return T(2)*vec - Vec(1);
}
void CheckTranslate() const
{
std::cout << "--- Checking translate" << std::endl;
Vec startPoint = this->RandomVector();
std::cout << " Starting point: " << startPoint << std::endl;
Vec translateAmount = this->RandomVector();
std::cout << " Translation amount: " << translateAmount << std::endl;
Transform translate = vtkm::Transform3DTranslate(translateAmount);
Vec translated1 = vtkm::Transform3DPoint(translate, startPoint);
std::cout << " First translation: " << translated1 << std::endl;
VTKM_TEST_ASSERT(test_equal(translated1, startPoint+translateAmount),
"Bad translation.");
Vec translated2 = vtkm::Transform3DPoint(translate, translated1);
std::cout << " Second translation: " << translated2 << std::endl;
VTKM_TEST_ASSERT(test_equal(translated2, startPoint+T(2)*translateAmount),
"Bad translation.");
// Vectors should be invarient to translation.
translated1 = vtkm::Transform3DVector(translate, startPoint);
std::cout << " Translated vector: " << translated1 << std::endl;
VTKM_TEST_ASSERT(test_equal(translated1, startPoint), "Bad translation.");
}
void CheckScale() const
{
std::cout << "--- Checking scale" << std::endl;
Vec startPoint = this->RandomVector();
std::cout << " Starting point: " << startPoint << std::endl;
Vec scaleAmount = this->RandomVector();
std::cout << " Scale amount: " << scaleAmount << std::endl;
Transform scale = vtkm::Transform3DScale(scaleAmount);
Vec scaled1 = vtkm::Transform3DPoint(scale, startPoint);
std::cout << " First scale: " << scaled1 << std::endl;
VTKM_TEST_ASSERT(test_equal(scaled1, startPoint*scaleAmount),
"Bad scale.");
Vec scaled2 = vtkm::Transform3DPoint(scale, scaled1);
std::cout << " Second scale: " << scaled2 << std::endl;
VTKM_TEST_ASSERT(test_equal(scaled2, startPoint*scaleAmount*scaleAmount),
"Bad scale.");
// Vectors should scale the same as points.
scaled1 = vtkm::Transform3DVector(scale, startPoint);
std::cout << " Scaled vector: " << scaled1 << std::endl;
VTKM_TEST_ASSERT(test_equal(scaled1, startPoint*scaleAmount),
"Bad scale.");
}
void CheckRotate() const
{
std::cout << "--- Checking rotate" << std::endl;
Vec startPoint = this->RandomVector();
std::cout << " Starting point: " << startPoint << std::endl;
const T ninetyDegrees = T(vtkm::Pi_2());
std::cout << "--Rotate 90 degrees around X" << std::endl;
Transform rotateX = vtkm::Transform3DRotateX(ninetyDegrees);
Vec rotated1 = vtkm::Transform3DPoint(rotateX, startPoint);
std::cout << " First rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(startPoint[0],-startPoint[2],startPoint[1])),
"Bad rotate.");
Vec rotated2 = vtkm::Transform3DPoint(rotateX, rotated1);
std::cout << " Second rotate: " << rotated2 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated2, Vec(startPoint[0],-startPoint[1],-startPoint[2])),
"Bad rotate.");
// Vectors should rotate the same as points.
rotated1 = vtkm::Transform3DVector(rotateX, startPoint);
std::cout << " Vector rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(startPoint[0],-startPoint[2],startPoint[1])),
"Bad rotate.");
std::cout << "--Rotate 90 degrees around Y" << std::endl;
Transform rotateY = vtkm::Transform3DRotateY(ninetyDegrees);
rotated1 = vtkm::Transform3DPoint(rotateY, startPoint);
std::cout << " First rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(startPoint[2],startPoint[1],-startPoint[0])),
"Bad rotate.");
rotated2 = vtkm::Transform3DPoint(rotateY, rotated1);
std::cout << " Second rotate: " << rotated2 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated2, Vec(-startPoint[0],startPoint[1],-startPoint[2])),
"Bad rotate.");
// Vectors should rotate the same as points.
rotated1 = vtkm::Transform3DVector(rotateY, startPoint);
std::cout << " Vector rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(startPoint[2],startPoint[1],-startPoint[0])),
"Bad rotate.");
std::cout << "--Rotate 90 degrees around Z" << std::endl;
Transform rotateZ = vtkm::Transform3DRotateZ(ninetyDegrees);
rotated1 = vtkm::Transform3DPoint(rotateZ, startPoint);
std::cout << " First rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(-startPoint[1],startPoint[0],startPoint[2])),
"Bad rotate.");
rotated2 = vtkm::Transform3DPoint(rotateZ, rotated1);
std::cout << " Second rotate: " << rotated2 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated2, Vec(-startPoint[0],-startPoint[1],startPoint[2])),
"Bad rotate.");
// Vectors should rotate the same as points.
rotated1 = vtkm::Transform3DVector(rotateZ, startPoint);
std::cout << " Vector rotate: " << rotated1 << std::endl;
VTKM_TEST_ASSERT(
test_equal(rotated1, Vec(-startPoint[1],startPoint[0],startPoint[2])),
"Bad rotate.");
}
};
struct TryTransformsFunctor
{
template<typename T>
void operator()(T) const
{
TransformTests<T> tests;
tests.CheckTranslate();
tests.CheckScale();
tests.CheckRotate();
}
};
void TestTransforms()
{
vtkm::UInt32 seed = static_cast<vtkm::UInt32>(time(NULL));
std::cout << "Seed: " << seed << std::endl;
g_RandomGenerator.seed(seed);
vtkm::testing::Testing::TryTypes(TryTransformsFunctor(),
vtkm::TypeListTagFieldScalar());
}
} // anonymous namespace
int UnitTestTransform3D(int, char *[])
{
return vtkm::testing::Testing::Run(TestTransforms);
}