Transform3D.h
// -*- C++ -*-
// CLASSDOC OFF
// $Id: Transform3D.h,v 1.5 1999/04/18 13:36:44 evc Exp $
// ---------------------------------------------------------------------------
// CLASSDOC ON
//
// Hep geometrical 3D Transformation class
//
// This file is part of Geant4 (simulation toolkit for HEP).
//
// Author: Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch>
//
// ******************************************
// * *
// * Transform *
// * / / \ \ *
// * -------- / \ -------- *
// * / / \ \ *
// * Rotate Translate Reflect Scale *
// * / | \ / | \ / | \ / | \ *
// * X Y Z X Y Z X Y Z X Y Z *
// * *
// ******************************************
//
// Identity transformation:
// HepTransform3D::Identity - global identity transformation;
// any constructor without parameters, e.g. HepTransform3D();
// m.setIdentity() - set "m" to identity;
//
// General transformations:
// HepTransform3D(m,v) - transformation given by HepRotation "m"
// and Hep3Vector "v";
// HepTransform3D(a0,a1,a2, b0,b1,b2) - transformation given by initial
// and transformed positions of three points;
// Rotations:
// HepRotate3D(m) - rotation given by HepRotation "m";
// HepRotate3D(ang,v) - rotation through the angle "ang" around
// vector "v";
// HepRotate3D(ang,p1,p2) - rotation through the angle "ang"
// counterclockwise around the axis given by
// two points p1->p2;
// HepRotate3D(a1,a2, b1,b2) - rotation around the origin defined by initial
// and transformed positions of two points;
// HepRotateX3D(ang) - rotation around X-axis;
// HepRotateY3D(ang) - rotation around Y-axis;
// HepRotateZ3D(ang) - rotation around Z-axis;
//
// Translations:
// HepTranslate3D(v) - translation given by Hep3Vector "v";
// HepTranslate3D(dx,dy,dz) - translation on vector (dx,dy,dz);
// HepTraslateX3D(dx) - translation along X-axis;
// HepTraslateY3D(dy) - translation along Y-axis;
// HepTraslateZ3D(dz) - translation along Z-axis;
//
// Reflections:
// HepReflect3D(a,b,c,d) - reflection in the plane a*x+b*y+c*z+d=0;
// HepReflect3D(normal,p) - reflection in the plane going through "p"
// and whose normal is equal to "normal";
// HepReflectX3D(a) - reflect X in the plane x=a (default a=0);
// HepReflectY3D(a) - reflect Y in the plane y=a (default a=0);
// HepReflectZ3D(a) - reflect Z in the plane z=a (default a=0);
//
// Scalings:
// HepScale3D(sx,sy,sz) - general scaling with factors "sx","sy","sz"
// along X, Y and Z;
// HepScale3D(s) - scaling with constant factor "s" along all
// directions;
// HepScaleX3D(sx) - scale X;
// HepScaleY3D(sy) - scale Y;
// HepScaleZ3D(sz) - scale Z;
//
// Inverse transformation:
// m.inverse() or - returns inverse transformation;
// m^-1
//
// Compound transformation:
// m3 = m2 * m1 - it is relatively slow in comparison with
// transformation of a vector. Use parenthesis
// to avoid this operation (see remark below);
// Transformation of point:
// p2 = m * p1
//
// Transformation of vector:
// v2 = m * v1
//
// Transformation of normal:
// n2 = m * n1
//
// The following table explains how different transformations affect
// point, vector and normal. "+" means affect, "-" means do not affect,
// "*" meas affect but in different way than "+"
//
// Point Vector Normal
// -------------+-------+-------+-------
// Rotation ! + ! + ! +
// Translation ! + ! - ! -
// Reflection ! + ! + ! *
// Scaling ! + ! + ! *
// -------------+-------+-------+-------
//
// Example of the usage:
//
// HepTransform3D m1, m2, m3;
// HepVector3D v2, v1(0,0,0);
//
// m1 = HepRotate3D(angle, HepVector3D(1,1,1));
// m2 = HepTranslate3D(dx,dy,dz);
// m3 = m1^-1;
//
// v2 = m3*(m2*(m1*v1));
//
// Remark: For the reason that the operator * is left associative the notation
// v2 = m3*(m2*(m1*v1)) is much more effective then the notation
// v2 = m3*m2*m1*v1. In the first case three operations "Transform*Vector"
// are executed, in the second case two operations "Transform*Transform" and
// one "Transform*Vector" are perfomed. "Transform*Transform" is roughly
// 3 times slower than "Transform*Vector".
//
// History:
// 24.09.96 E.Chernyaev - initial version
//
// 26.02.97 E.Chernyaev
// - added global Identity by request of John Allison
// (to avoid problems with compilation on HP)
// - added getRotation and getTranslation
#ifndef HepTransform3D_hh
#define HepTransform3D_hh
#include "CLHEP/config/CLHEP.h"
class Hep3Vector;
class HepRotation;
class HepPoint3D;
class HepVector3D;
class HepNormal3D;
class HepTransform3D {
protected:
HepDouble xx, xy, xz, dx, // 4x3 Transformation Matrix
yx, yy, yz, dy,
zx, zy, zz, dz;
// Protected constructor
HepTransform3D(HepDouble XX, HepDouble XY, HepDouble XZ, HepDouble DX,
HepDouble YX, HepDouble YY, HepDouble YZ, HepDouble DY,
HepDouble ZX, HepDouble ZY, HepDouble ZZ, HepDouble DZ)
: xx(XX), xy(XY), xz(XZ), dx(DX),
yx(YX), yy(YY), yz(YZ), dy(DY),
zx(ZX), zy(ZY), zz(ZZ), dz(DZ) {}
// Set transformation matrix
void setTransform(HepDouble XX, HepDouble XY, HepDouble XZ, HepDouble DX,
HepDouble YX, HepDouble YY, HepDouble YZ, HepDouble DY,
HepDouble ZX, HepDouble ZY, HepDouble ZZ, HepDouble DZ) {
xx = XX; xy = XY; xz = XZ; dx = DX;
yx = YX; yy = YY; yz = YZ; dy = DY;
zx = ZX; zy = ZY; zz = ZZ; dz = DZ;
}
public:
// Global identity transformation
static const HepTransform3D Identity;
// For RogueWave...
HepBoolean operator == (const HepTransform3D& transform) const {
return this == &transform?true:false;
}
// Constructor: identity
HepTransform3D()
: xx(1), xy(0), xz(0), dx(0),
yx(0), yy(1), yz(0), dy(0),
zx(0), zy(0), zz(1), dz(0) {}
// Constructor: rotation and then translation
inline HepTransform3D(const HepRotation &m, const Hep3Vector &v);
// Constructor: transformation of basis
HepTransform3D
(const HepPoint3D &fr0, const HepPoint3D &fr1, const HepPoint3D &fr2,
const HepPoint3D &to0, const HepPoint3D &to1, const HepPoint3D &to2);
// Copy constructor
HepTransform3D(const HepTransform3D & m)
: xx(m.xx), xy(m.xy), xz(m.xz), dx(m.dx),
yx(m.yx), yy(m.yy), yz(m.yz), dy(m.dy),
zx(m.zx), zy(m.zy), zz(m.zz), dz(m.dz) {}
// Assignment
HepTransform3D& operator=(const HepTransform3D &m) {
setTransform(m.xx, m.xy, m.xz, m.dx,
m.yx, m.yy, m.yz, m.dy,
m.zx, m.zy, m.zz, m.dz);
return *this;
}
// Set identity
void setIdentity() {
xy = xz = dx = yx = yz = dy = zx = zy = dz = 0; xx = yy = zz = 1;
}
// Inverse transformation
HepTransform3D inverse() const;
// Inverse transformation: m^-1
HepTransform3D operator^(HepInt) const { return inverse(); }
// Compound transformations
inline HepTransform3D operator*(const HepTransform3D &b) const;
// Transformation of point
inline HepPoint3D operator*(const HepPoint3D &p) const;
// Transformation of vector
inline HepVector3D operator*(const HepVector3D &v) const;
// Transformation of normal
inline HepNormal3D operator*(const HepNormal3D &n) const;
// Extract the rotation matrix
inline HepRotation getRotation() const;
// Extract the translation vector
inline Hep3Vector getTranslation() const;
};
// R O T A T I O N S
class HepRotate3D : public HepTransform3D {
public:
HepRotate3D() : HepTransform3D() {}
inline HepRotate3D(const HepRotation &m);
HepRotate3D(HepDouble a, const HepPoint3D &p1, const HepPoint3D &p2);
inline HepRotate3D(HepDouble a, const HepVector3D &v);
inline HepRotate3D(const HepPoint3D &fr1, const HepPoint3D &fr2,
const HepPoint3D &to1, const HepPoint3D &to2);
};
class HepRotateX3D : public HepRotate3D {
public:
HepRotateX3D() : HepRotate3D() {}
HepRotateX3D(HepDouble a) {
HepDouble cosa = cos(a), sina = sin(a);
setTransform(1,0,0,0, 0,cosa,-sina,0, 0,sina,cosa,0);
}
};
class HepRotateY3D : public HepRotate3D {
public:
HepRotateY3D() : HepRotate3D() {}
HepRotateY3D(HepDouble a) {
HepDouble cosa = cos(a), sina = sin(a);
setTransform(cosa,0,sina,0, 0,1,0,0, -sina,0,cosa,0);
}
};
class HepRotateZ3D : public HepRotate3D {
public:
HepRotateZ3D() : HepRotate3D() {}
HepRotateZ3D(HepDouble a) {
HepDouble cosa = cos(a), sina = sin(a);
setTransform(cosa,-sina,0,0, sina,cosa,0,0, 0,0,1,0);
}
};
// T R A N S L A T I O N S
class HepTranslate3D : public HepTransform3D {
public:
HepTranslate3D() : HepTransform3D() {}
inline HepTranslate3D(const Hep3Vector &v);
HepTranslate3D(HepDouble x, HepDouble y, HepDouble z)
: HepTransform3D(1,0,0,x, 0,1,0,y, 0,0,1,z) {}
};
class HepTranslateX3D : public HepTranslate3D {
public:
HepTranslateX3D() : HepTranslate3D() {}
HepTranslateX3D(HepDouble x) : HepTranslate3D(x, 0, 0) {}
};
class HepTranslateY3D : public HepTranslate3D {
public:
HepTranslateY3D() : HepTranslate3D() {}
HepTranslateY3D(HepDouble y) : HepTranslate3D(0, y, 0) {}
};
class HepTranslateZ3D : public HepTranslate3D {
public:
HepTranslateZ3D() : HepTranslate3D() {}
HepTranslateZ3D(HepDouble z) : HepTranslate3D(0, 0, z) {}
};
// R E F L E C T I O N S
class HepReflect3D : public HepTransform3D {
protected:
HepReflect3D(HepDouble XX, HepDouble XY, HepDouble XZ, HepDouble DX,
HepDouble YX, HepDouble YY, HepDouble YZ, HepDouble DY,
HepDouble ZX, HepDouble ZY, HepDouble ZZ, HepDouble DZ)
: HepTransform3D(XX,XY,XZ,DX, YX,YY,YZ,DY, ZX,ZY,ZZ,DZ) {}
public:
HepReflect3D() : HepTransform3D() {}
HepReflect3D(HepDouble a, HepDouble b, HepDouble c, HepDouble d);
inline HepReflect3D(const HepNormal3D &normal, const HepPoint3D &point);
};
class HepReflectX3D : public HepReflect3D {
public:
HepReflectX3D(HepDouble x=0) : HepReflect3D(-1,0,0,x+x, 0,1,0,0, 0,0,1,0) {}
};
class HepReflectY3D : public HepReflect3D {
public:
HepReflectY3D(HepDouble y=0) : HepReflect3D(1,0,0,0, 0,-1,0,y+y, 0,0,1,0) {}
};
class HepReflectZ3D : public HepReflect3D {
public:
HepReflectZ3D(HepDouble z=0) : HepReflect3D(1,0,0,0, 0,1,0,0, 0,0,-1,z+z) {}
};
// S C A L I N G S
class HepScale3D : public HepTransform3D {
public:
HepScale3D() : HepTransform3D() {}
HepScale3D(HepDouble x, HepDouble y, HepDouble z)
: HepTransform3D(x,0,0,0, 0,y,0,0, 0,0,z,0) {}
HepScale3D(HepDouble s)
: HepTransform3D(s,0,0,0, 0,s,0,0, 0,0,s,0) {}
};
class HepScaleX3D : public HepScale3D {
public:
HepScaleX3D() : HepScale3D() {}
HepScaleX3D(HepDouble x) : HepScale3D(x, 1, 1) {}
};
class HepScaleY3D : public HepScale3D {
public:
HepScaleY3D() : HepScale3D() {}
HepScaleY3D(HepDouble y) : HepScale3D(1, y, 1) {}
};
class HepScaleZ3D : public HepScale3D {
public:
HepScaleZ3D() : HepScale3D() {}
HepScaleZ3D(HepDouble z) : HepScale3D(1, 1, z) {}
};
// I N L I N E S F O R T R A N S F O R M A T I O N
#include "CLHEP/Vector/ThreeVector.h"
#include "CLHEP/Vector/Rotation.h"
#include "CLHEP/Geometry/Point3D.h"
#include "CLHEP/Geometry/Vector3D.h"
#include "CLHEP/Geometry/Normal3D.h"
inline HepTransform3D::HepTransform3D(const HepRotation &m, const Hep3Vector &v)
{
Hep3Vector w;
w = m * Hep3Vector(1,0,0); xx = w.x(); yx = w.y(); zx = w.z();
w = m * Hep3Vector(0,1,0); xy = w.x(); yy = w.y(); zy = w.z();
w = m * Hep3Vector(0,0,1); xz = w.x(); yz = w.y(); zz = w.z();
dx = v.x(); dy = v.y(); dz = v.z();
}
inline HepTransform3D HepTransform3D::operator*(const HepTransform3D &b) const {
return HepTransform3D
(xx*b.xx+xy*b.yx+xz*b.zx, xx*b.xy+xy*b.yy+xz*b.zy,
xx*b.xz+xy*b.yz+xz*b.zz, xx*b.dx+xy*b.dy+xz*b.dz+dx,
yx*b.xx+yy*b.yx+yz*b.zx, yx*b.xy+yy*b.yy+yz*b.zy,
yx*b.xz+yy*b.yz+yz*b.zz, yx*b.dx+yy*b.dy+yz*b.dz+dy,
zx*b.xx+zy*b.yx+zz*b.zx, zx*b.xy+zy*b.yy+zz*b.zy,
zx*b.xz+zy*b.yz+zz*b.zz, zx*b.dx+zy*b.dy+zz*b.dz+dz);
}
inline HepPoint3D HepTransform3D::operator*(const HepPoint3D &p) const {
HepDouble x = p.x(), y = p.y(), z = p.z();
return HepPoint3D(x*xx+y*xy+z*xz+dx, x*yx+y*yy+z*yz+dy, x*zx+y*zy+z*zz+dz);
}
inline HepVector3D HepTransform3D::operator*(const HepVector3D &v) const {
HepDouble x = v.x(), y = v.y(), z = v.z();
return HepPoint3D(x*xx+y*xy+z*xz, x*yx+y*yy+z*yz, x*zx+y*zy+z*zz);
}
inline HepNormal3D HepTransform3D::operator*(const HepNormal3D &n) const {
HepDouble x = n.x(), y = n.y(), z = n.z();
return HepNormal3D(x*(yy*zz-yz*zy) + y*(yz*zx-yx*zz) + z*(yx*zy-yy*zx),
x*(zy*xz-zz*xy) + y*(zz*xx-zx*xz) + z*(zx*xy-zy*xx),
x*(xy*yz-xz*yy) + y*(xz*yx-xx*yz) + z*(xx*yy-xy*yx));
}
inline HepRotation HepTransform3D::getRotation() const {
HepRotation m;
return m.rotateAxes(Hep3Vector(xx,yx,zx),
Hep3Vector(xy,yy,zy),
Hep3Vector(xz,yz,zz));
}
inline Hep3Vector HepTransform3D::getTranslation() const {
return Hep3Vector(dx,dy,dz);
}
// I N L I N E S F O R R O T A T I O N
inline HepRotate3D::HepRotate3D(const HepRotation &m) {
Hep3Vector w;
w = m * Hep3Vector(1,0,0); xx = w.x(); yx = w.y(); zx = w.z();
w = m * Hep3Vector(0,1,0); xy = w.x(); yy = w.y(); zy = w.z();
w = m * Hep3Vector(0,0,1); xz = w.x(); yz = w.y(); zz = w.z();
dx = 0; dy = 0; dz = 0;
}
inline HepRotate3D::HepRotate3D(HepDouble a, const HepVector3D &v) {
*this = HepRotate3D(a, HepPoint3D(0,0,0), HepPoint3D(v.x(),v.y(),v.z()));
}
inline HepRotate3D::HepRotate3D(const HepPoint3D &fr1, const HepPoint3D &fr2,
const HepPoint3D &to1, const HepPoint3D &to2)
: HepTransform3D(HepPoint3D(0,0,0),fr1,fr2, HepPoint3D(0,0,0),to1,to2) {}
// I N L I N E S F O R T R A N S L A T I O N
inline HepTranslate3D::HepTranslate3D(const Hep3Vector &v)
: HepTransform3D(1,0,0,v.x(), 0,1,0,v.y(), 0,0,1,v.z()) {}
// I N L I N E S F O R R E F L E C T I O N
inline HepReflect3D::HepReflect3D(const HepNormal3D &n, const HepPoint3D &p) {
*this = HepReflect3D(n.x(), n.y(), n.z(), -n*p);
}
#endif /* HepTransform3D_hh */
Generated by GNU enscript 1.6.1.