CLHEP Vector

The CLHEP Vector module consists of four classes:

Corresponding header files are:
   CLHEP/Vector/ThreeVector.h
   CLHEP/Vector/Rotation.h
   CLHEP/Vector/LorentzVector.h
   CLHEP/Vector/LorentzRotation.h


Hep3Vector

Hep3Vector is a general three vector class, which can be used for description of different vectors in 3D (see Figure 1).

  
Figure 1: Components of three vector: x,y,z - basic components, - azimuth angle, - polar angle, - magnitude, - transverse component

Using the Hep3Vector class you should remember that it contains only common features of three vectors and does not have member functions which are specific for some particular vector values. For example, it has not translate() function because translation may be applied to vertices (position vectors) but has no sense for free vectors. An implemention of some specific three vector classes you can find in the CLHEP Geometry module.

Declaration

Hep3Vector has been implemented as a vector of three HepDouble variables. You may specify their values when declare a Hep3Vector objects. By default all the components are initialized by zero.

   #include "CLHEP/Vector/ThreeVector.h"
   Hep3Vector v1;             // v1 initialized by (0.,0.,0.)
   Hep3Vector v2(1.,1.,1.);   // v2 initialized by (1.,1.,1.) 
   Hep3Vector v3(v1);         // v3 = v1

There are three predefined global objects of the Hep3Vector class: HepXHat(1.,0.,0.), HepYHat(0.,1.,0.) and HepZHat(0.,0.,1.).

Access to the basic components

You can get values of the basic components either by name using functions x(), y(), z() or by index using operator():

   HepDouble xx = v1.x();     // access to components by name
   HepDouble yy = v1.y();
   HepDouble zz = v1.z();
   xx = v1(0);                // access to componets by index
   yy = v1(1);
   zz = v1(2);

Member functions setX(), setY(), setZ() allow to set the basic components:

   v1.setX(1.);               // set v1 to (1.,2.,3.)
   v1.setY(2.);
   v1.setZ(3.);

Noncartesian coordinates of the vector

The basic components (x,y,z) represent coordinates of the vector in cartesian coodinate system. Other popular coordinate systems are spherical and cylindrical .

Member functions mag(), mag2(), theta(), cosTheta() and phi() allow to get information on components of the vector in spherical coordinate system:

   HepDouble m     = v1.mag();      // get magnitude (rho)
   HepDouble m2    = v1.mag2();     // get magnitude squared

   HepDouble t     = v1.theta();    // get polar angle
   HepDouble cost  = v1.cosTheta(); // get cosine of polar angle
   HepDouble p     = v1.phi();      // get azimuth angle

Member functions setMag(), setTheta(), setPhi() can be used to set components of the vector in spherical coodinate system:

   v1.setTheta(30.*deg); // set theta keeping magnitude and phi unchanged
   v1.setPhi(pi*rad);    // set phi keeping magnitude and theta unchanged
   v1.setMag(10.);       // set magnitude keeping theta and phi unchanged

To get information on the transverse component of the vector (the r-coordinate in the cylindrical coordinate system) and on projection of the transverse component to another vector you can use member functions perp() and perp2():

   HepDouble pp    = v1.perp();     // get transverse component
   HepDouble pp2   = v1.perp2();    // get transverse component squared
   HepDouble ppv2  = v1.perp(v2);   // get transverse conponent w.r.t.
                                    // another vector
   HepDouble pp2v2 = v1.perp2(v2);  // get transverse component w.r.t.
                                    // another vector squared

Arithmetic operators

The Hep3Vector class provides a set of arithmetic operators to add, subtract and scale vectors:

   v3  = -v1;                 // unary minus
   v1  = v2 + v3;             // additions
   v1 += v3;
   v1  = v1 - v2;             // subtractions
   v1 -= v3;
   v1 *= 10.;                 // scalings
   v1  = 5.*v2;

Comparisons

You can test two Hep3Vector objects for equality or inequality:

   if (v1 == v2) {...}        // test for equality
   if (v1 != v3) {...}        // test for inequality

Output to a stream

The following lines:

   Hep3Vector v(1.,2.,3.):
   cout << v << endl;
produce the following output:
   (1,2,3)

Related vectors

Member functions unit() and orthogonal() return respectively a vector of unit length parallel to the current and a vector orthogonal to the current:

   v2 = v1.unit();            // get unit vector parallel to v1
   v3 = v1.orthogonal();      // get vector orthogonal to v1

Scalar and vector products

To find a scalar product of the current vector with another one you can use either the member function dot() or the operator *:

   HepDouble s  = v1.dot(v2); // scalar product
   s  = v1*v2;

Member function cross() returns a vector product of the current vector with another one:

   v3 = v1.cross(v2)          // vector product

Angle between two vectors

Member function angle() allows to find angle between the current vector and another one:

   HepDouble ang = v1.angle(v2);  // get angle w.r.t. another vector

Rotation around axes

Member functions rotateX(), rotateY() and rotateZ() allow to rotate Hep3Vector around the axes:

   v1.rotateX(30.*deg);       // rotation around axes
   v1.rotateY(halfpi*rad);
   v1.rotateZ(45.*deg);

Rotation around specified vector

You can rotate Hep3Vector around arbitrary axis using the member function rotate():

   v1.rotate(pi/4.*rad,v2);   // rotation around specified vector

Rotation by HepRotation

You can rotate Hep3Vector by an object of the HepRotation class (see next section) using either the member function transform() or the operator *= or the operator * of the HepRotation class:

   HepRotation M;
   v1.transform(M);           // rotation with a rotation matrix
   v1  = M*v1;
   v1 *= M;                   // Attention: v1 = M * v1 !!!

Transformation from rotated frame

Suppose that there is a direction given by two angles and . Let assign a frame , , to the direction in the following way: let coinsides with the direction, lies in the -plane and lies in the xy-plane and be perpendicular to the -plane.

A frame , , originally coinsided with global frame can be transformed to the frame , , by two consequentive rotations, first through angle around y-axis and then trough angle around z-axis. Such a transformation is described by the following matrix:

The same transformation describes transformation from the rotated frame to the original frame (see answers on frequently asked questions for the HepRotation class) and is used in member function rotateUz() for transformation of the current vector from the rotated frame to the original frame:

   v1.rotateUz(direction); // direction must be Hep3Vector of unit length


HepRotation

The HepRotation class describes a rotation of objects of the Hep3Vector class. It has been implemented as a matrix of HepDoubles:

Before to consider the public interface of the HepRotation class, let us to answer on a couple of frequently asked questions.

a) What rotation does the HepRotation class describe?

The HepRotation class describes so called active rotation, i.e. rotation of objects inside static system of coordinates. In case if you want to rotate the frame and want to know the coordinates of objects in the rotated coordinate system then you should apply the inverse rotation to the objects. But if you want to trasform coordinates from the rotated frame to the original frame then you should again apply the direct transformation.

When we speak about rotation around specified axis we mean counterclockwise rotation around positive direction of the axis.

b) How to find a compound rotation of two rotations and ?

The operator * of the HepRotation class has been implemented in a way that follows the mathematical notation of a product of two matrices which describe two consequentive rotations. It means that you should place the second rotation from the left side of the * and the first one from the right side:

   R = R2 * R1;

When you apply several consequentive rotations to a Hep3Vector object you should remember that the operator * is left associative. It means that the notation with parenthesis:

   v = R3 * (R2 * (R1 * v));
is more effective then the notation without parenthesis:
   v = R3 * R2 * R1 * v;
because multiplication of a HepRotation object with a Hep3Vector object is less expensive then multiplication of two HepRotation objects.

Declaration

There is no direct write access to the members of the HepRotation class. It ensures that a HepRotation object always describes a real rotation. By default a HepRotation object is initialized by null rotation, i.e. an identity matrix:

   #include "CLHEP/Vector/Rotation.h"
   HepRotation R;             // R initialized as identity
   HepRotation M(R);          // M = R

Access to the matrix components

Member functions xx(), xy(), xz(), yx(), yy(), yz(), zx(), zy(), zz() allow to get values of the elements of the rotation matrix:

   HepDouble xx, xy, xz, yx, yy, yz, zx, zy, zz;
   xx = R.xx(); xy = R.xy(); xz = R.xz();
   yx = R.yx(); yy = R.yy(); yz = R.yz();
   zx = R.zx(); zy = R.zy(); zz = R.zz();
You can also access the matrix elements using operator():
   xx = R(0,0); xy = R(0,1); xz = R(0,2);
   yx = R(1,0); yy = R(1,1); yz = R(1,2);
   zx = R(2,0); zy = R(2,1); zz = R(2,2);

Comparisons

You can test two HepRotation objects for equality/inequality or test whether a HepRotation object is an identity matrix:

   if (R == M) {...}          // test for equality
   if (R != M) {...}          // test for inequality
   if (R.isIdentity()) {...}  // test for identity

Rotation of vector

Mathematically rotation of a vector is expressed as:

or in more brief form:

The HepRotation class provides an operator * which allows to express a rotation of a Hep3Vector object in a similar form:

   Hep3Vector v(1.,1.,1.);
   v = R * v;

See also the transform() function and operator *= of the Hep3Vector class.

Compound rotation

To find the product of two rotations you can use either operator * or operator *= or member function transform(). Note that the operator *= should be applied to the later rotation.

   HepRotation a,b,c; 

   c  = b * a;                // product of two rotations
   a *= b;                    // Attention: a = a * b  
   c  = a.transform(b);       // a = b * a and then c = a

Rotation around axes

The following matrices describe counterclockwise rotations around the coordinate axes:

Member functions rotateX(), rotateY(), rotateZ() allow to add a rotation around respective axis to the current rotation and return the result.

   a.rotateX(30.*deg);        // rotation around x-axis
   a.rotateY(pi*rad);         // rotation around y-axis
   a.rotateZ(45.*deg);        // rotation around z-axis

Rotation around arbitrary axis

The following formula expresses a counterclockwise rotation through angle around unit vector (x,y,z):

Member function rotate() allows to add a rotation around arbitrary vector (not neccessary a unit one) to the current rotation and returns the result. Note that there are two different signatures of the rotate() function:

   Hep3Vector axis(1.,1.,1.); 
   Hep3Vector *pAxis = &axis; 
   a.rotate(45.*deg, axis);     // rotation around specified vector
   b.rotate(halfpi*rad, pAxis);

Member function getAngleAxis() finds a unit vector (x,y,z) and an angle , which describe the same rotartion as the current:

    HepDouble  angle;
    Hep3Vector axis;
    a.getAngleAxis(angle, axis); // get angle and axis of rotation
To do this the following equations are used:

Rotation of local axes

Suppose that there is a 3D object and its local axes , , originally coincide with global axes. If it is known the positions of local axes , , after rotation then matrix which describes the rotation looks as follows:

Member function rotateAxes() adds a rotation of local axes to the current rotation and returns the result:

   Hep3Vector newX(0.,1.,0.), 
              newY(0.,0.,1.),
              newZ(1.,0.,0.);
   a.rotateAxes(newX, newY, newZ); // rotation of local axes

Member functions thetaX(), thetaY(), thetaZ(), phiX(), phiY(), phiZ() return azimuth and polar angles of the rotated axes:

    HepDouble tx = a.thetaX();
    HepDouble ty = a.thetaY();
    HepDouble tz = a.thetaZ();

    HepDouble px = a.phiX();
    HepDouble py = a.phiY();
    HepDouble pz = a.phiZ();

Inverse rotation

In case of rotation the inverse matrix is equal to the transpose one:

Member function inverse() returns the inverse rotation keeping the current rotation unchanged. Member function invert() inverts the current rotation and returns the result:

    b = a.inverse();          // b is inverse of a, a is unchanged
    b = a.invert();           // invert a and then b = a


HepLorentzVector

HepLorentzVector is a general four-vector class, which can be used either for description of position and time (x,y,z,t) or momentum and energy .

Declaration

HepLorentzVector has been implemented as a pair of Hep3Vector and HepDouble, thus in total there are four HepDouble variables. By default all the components are initialized by zero.

   #include "CLHEP/Vector/LorentzVector.h"
   HepLorentzVector v1;       // v1 initialized by (0.,0.,0.,0.)
   HepLorentzVector v2(1.,1.,1.,1.);
   HepLorentzVector v3(v1);
   HepLorentzVector v4(Hep3Vector(1.,2.,3.),4.);

Conversion to Hep3Vector

HepLorentzVector has a conversion operator to Hep3Vector. It returns a reference to the vector component of the HepLorentzVector object. So, HepLorentzVector can be used everywhere Hep3Vector expected.

   Hep3Vector v3;
   HepLorentzVector v4;
    ...
   HepDouble a = v3.dot(v4);       
   Hep3Vector(v4).rotateZ(30.*deg); // rotate vector component of v4

Access to the components

There are two sets of access functions to the componets of a HepLorentzVector object: x(), y(), z(), t() and px(), py(), pz(), e(). Both sets return the same values but the first set is more relevant for use where HepLorentzVector describes a combination of position and time and the second set is more relevant where HepLorentzVector describes momentum and energy:

   HepDouble xx = v.x();      // get position and time
   HepDouble yy = v.y();
   HepDouble zz = v.z();
   HepDouble tt = v.t();

   HepDouble px = v.px();     // get momentum and energy
   HepDouble py = v.py();
   HepDouble pz = v.pz();
   HepDouble ee = v.e();

The components of HepLorentzVector can be also accessed by index:

   xx = v(0);                 // access to components by index
   yy = v(1);
   zz = v(2);
   tt = v(3);

You can use the vect() member function to get the vector component of HepLorentzVector:

   Hep3Vector p = v.vect();  // get vector component

For setting components it also can be used two sets of member functions: setX(), setX(), setZ(), setT() and setPx(), setPy(), setPz(), setE().

   v.setX(1.);                // set position and time
   v.setY(2.);
   v.setZ(3.);
   v.setT(4.);
 
   v.setPx(1.);               // set momentum and energy
   v.setPy(2.);      
   v.setPz(3.);
   v.setE(4.);

The vector component of HepLorentzVector can be set by the setVect() function:

   v.setVect(Hep3Vector());  // set vector component

Vector component in noncartesian coordinate systems

Member functions rho(), theta(), cosTheta() and phi() allow to get information on parameters of the vector component in spherical coordinate system:

   HepDouble m     = v4.rho();      // get length of the vector component
   HepDouble t     = v4.theta();    // get polar angle
   HepDouble cost  = v4.cosTheta(); // get cosine of polar angle
   HepDouble p     = v4.phi();      // get azimuth angle

Member functions setRho(), setTheta(), setPhi() can be used to change the vector component in spherical coodinate system, the spacial component remains unchanged:

   v4.setRho(10.);       // set length keeping theta and phi unchanged
   v4.setTheta(30.*deg); // set theta keeping length and phi unchanged
   v4.setPhi(pi*rad);    // set phi keeping length and theta unchanged

To get information about the r-coordinate of the vector component in the cylindrical coordinate system and on projection of it to specified three vector you can use member functions perp() and perp2():

   HepDouble pp    = v4.perp();     // get transverse component
   HepDouble pp2   = v4.perp2();    // get transverse component squared
   HepDouble ppv2  = v4.perp(v3);   // get transverse conponent w.r.t.
                                    // another vector
   HepDouble pp2v2 = v4.perp2(v4);  // get transverse component w.r.t.
                                    // another vector squared

Arithmetic operators

The HepLorentzVector class provides arithmetic operators to add and subtract four-vectors:

   v3  = -v1;                 // unary minus
   v1  = v2 + v3;             // additions
   v1 += v3;
   v1  = v1 - v2;             // subtractions
   v1 -= v3;

Comparisons

You can test two HepLorentzVector objects for equality or inequality:

   if (v1 == v2) {...}        // test for equality
   if (v1 != v3) {...}        // test for inequality

Output to a stream

The following lines:

   HepLorentzVector v(1.,2.,3.,4.):
   cout << v << endl;
produce the following output:
   (1,2,3,4)

Magnitude/Invariant mass

The magnitude squared of the four-vector is calculated as:

In case if is negative the magnitude is calculated as:

The member functions mag2() and mag() return respectively the magnitude squared and the magnitude:

   HepDouble s2 = v.mag2();
   HepDouble s  = v.mag();

Since in case of momentum and energy the magnitude of four-vector has meaning of invariat mass the HepLorentzVector class provides the member functions m2() and m() which are similar to mag2() and mag() but have more relevant names:

   HepDouble s2 = v.m2();
   HepDouble s  = v.m();

Scalar product

For calculation of scalar product of two four-vectors the following formula is used:

To find a scalar product of the current four-vector with another one you can use either the member function dot() or the operator *:

   HepDouble s  = v1.dot(v2); // scalar product
   s  = v1*v2;

Angle between two vectors

Member function angle() returns the angle between the vector component and another three vector:

   HepDouble ang = v4.angle(v3);  // get angle w.r.t. another vector

Light-cone components

Member functions plus() and minus() return respectively positive and negative light-cone components:

  HepDouble pcone = v.plus();    // pcone = t + z;
  HepDouble mcone = v.minus();   // mcone = t - z;

Lorentz boost

A boost in a general direction can be parameterized with three parameters which can be taken as the components of a three vector . Using this, an arbitrary active Lorentz boost transformation (from the rod frame to the original frame) can be written as:

where and

The member function boost() makes Lorentz boost transformation of HepLorentzVector from the rod frame to the original frame:

   v.boost(bx,by,bz);       // Lorentz boost
   v.boost(Hep3Vector(bx,by,bz));

The member fuction boostVector() returns a three vector of the spacial components divided by the time component:

   Hep3Vector beta = v.boostVector(); // beta = (x/t,y/t,z/t)

Rotation around axes

Member functions rotateX(), rotateY() and rotateZ() allow to rotate the vector component of HepLorentzVector around the axes:

   v1.rotateX(30.*deg);       // rotation around axes
   v1.rotateY(halfpi*rad);
   v1.rotateZ(45.*deg);

Rotation around specified vector

You can rotate the vector component of HepLorentzVector around arbitrary axis using the member function rotate():

   v1.rotate(pi/4.*rad,v2);   // rotation around specified vector

Transformation from rotated frame

The member function rotateUz transforms the vector component from rotated frame assigned to the unit vector to the original frame (for more details see similar function of the Hep3Vector class):

   v4.rotateUz(direction); // direction must be Hep3Vector of unit length

Rotation by HepRotation

You can rotate the vector component of HepLorentzVector by an object of the HepRotation class using either the member function transform() or the operator *=:

   HepRotation M;
   v1.transform(M);           // rotation with a rotation matrix
   v1 *= M;                   // Attention: v1 = M * v1 !!!

Transformation by HepLorentzRotation

You can transform HepLorentzVector by an object of the HepLorentzRotation class (see next section) using either the member function transform() or the operator *= or the operator * of the HepLorentzRotation class:

   HepLorentzRotation L;
   v1.transform(L);           // transformation with a Lorentz rotation matrix
   v1  = L*v1;
   v1 *= L;                   // Attention: v1 = L * v1 !!!


HepLorentzRotation

The HepLorentzRotation class describes Lorentz transformations of Lorentz vector. Lorentz transformations include Lorentz boosts and rotations (see HepRotation). You may think of HepLorentzRotation as of matrix of HepDoubles:

Declaration

By default a HepLorentzRotation object is initialized by the identity matrix. It can be also initialized by other HepLorentzRotation, by pure rotation or by Lorentz boost:

   #include "CLHEP/Vector/LorentzRotation.h"
   HepLorentzRotation L;             // L initialized as identity
   HepLorentzRotation M(L);          // M = L
   HepLorentzRotation LR(R);         // R is HepRotation
   HepLorentzRotation LB1(bx,by,bz); // Lorentz boost
   HepLorentzRotation LB2(Hep3Vector(bx,by,bz)); // LB2 = LB1
Matrix for Lorentz boost looks as follows:

where is a boost vector, and .

Access to the matrix components

Member functions xx(), xy(), xz(), xt(), yx(), yy(), yz(), yt(), zx(), zy(), zz(), zt(), tx(), ty(), tz(), tt() allow to get values of the elements of the HepLorentzRotation matrix:

   HepDouble xx, xy, xz, xt, yx, yy, yz, yt, zx, zy, zz, zt, tz, ty, tz, tt;
   xx = L.xx(); xy = L.xy(); xz = L.xz(); xt = L.xt();
   yx = L.yx(); yy = L.yy(); yz = L.yz(); yt = L.yt();
   zx = L.zx(); zy = L.zy(); zz = L.zz(); yt = L.zt();
   tx = L.zx(); ty = L.zy(); tz = L.zz(); tt = L.tt();
You can also access the matrix elements using operator():
   xx = L(0,0); xy = L(0,1); xz = L(0,2); xt = L(0,3); 
   yx = L(1,0); yy = L(1,1); yz = L(1,2); yt = L(1,3);
   zx = L(2,0); zy = L(2,1); zz = L(2,2); zt = L(2,3);
   tx = L(3,0); ty = L(3,1); tz = L(3,2); tt = L(3,3);

Comparisons

You can test two HepLorentzRotation objects for equality/inequality or test whether a HepLorentzRotation object is an identity matrix:

   if (L == M) {...}          // test for equality
   if (L != M) {...}          // test for inequality
   if (L.isIdentity()) {...}  // test for identity

Transformation of Lorentz vector

To apply HepLorentzRotation to HepLorentzVector you can use either the vectorMultiplication() function or operator *:

   HepLorentzVector v;
    ... 
   v = L.vectorMultiplication(v);
   v = L * v;
See also the transform() function and the operator *= of the HepLorentzVector class.

Compound transformation

There are four possibilities to find the product of two HepLorentzRotation transformations: operator *, operator *= and functions transform() and matrixMultiplication(). Note that the operator *= should be applied to the later transformation.

   HepLorentzRotation a,b,c; 
   c  = b * a;                     // product of two transformations
   c  = a.matrixMultiplication(b); // a is unchanged
   a *= b;                         // Attention: a = a * b  
   c  = a.transform(b);            // a = b * a and then c = a

Lorentz boost

Member function boost() allow to add a Lorentz boost to the current transformation and returns the result. Note that there are two different signatures of the boost() function:

   HepDouble bx,by,bz;
    ...
   a.boost(bx,by,bz);         // Lorentz boost
   a.boost(Hep3Vector(bx,by,bz));

Rotation around axes

Member functions rotateX(), rotateY(), rotateZ() allow to add a rotation around respective axis to the current transformation and return the result.

   a.rotateX(30.*deg);        // rotation around x-axis
   a.rotateY(pi*rad);         // rotation around y-axis
   a.rotateZ(45.*deg);        // rotation around z-axis

Rotation around arbitrary axis

Member function rotate() allows to add a rotation around arbitrary Hep3Vector to the current transformation and returns the result. Note that there are two different signatures of the rotate() function:

   Hep3Vector axis(1.,1.,1.); 
   Hep3Vector *pAxis = &axis; 
   a.rotate(45.*deg, axis);     // rotation around specified vector
   b.rotate(halfpi*rad, pAxis);

Inverse transformion

Matrix for the inverse transformation of HepLorentzRotation looks as follows:

Member function inverse() returns the inverse transformation keeping current HepLorentzRotation unchanged. Member function invert() inverts current HepLorentzRotation and returns the result:

    b = a.inverse();          // b is inverse of a, a is unchanged
    b = a.invert();           // invert a and then b = a


About this document ...

This document was generated using the LaTeX2HTML translator Version 96.1-g (July 19, 1996) Copyright © 1993, 1994, 1995, 1996, Nikos Drakos, Computer Based Learning Unit, University of Leeds.

The command line arguments were:
latex2html -split 0 -no_navigation -t "CLHEP Vector" vector.tex.

The translation was initiated by Evgueni Tcherniaev on Thu Feb 26 09:06:16 MET 1998


Evgueni Tcherniaev
Thu Feb 26 09:06:16 MET 1998