Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
Transform3D.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
2// Authors: W. Brown, M. Fischler, L. Moneta 2005
3
4/**********************************************************************
5 * *
6 * Copyright (c) 2005 , LCG ROOT MathLib Team *
7 * *
8 * *
9 **********************************************************************/
10
11// Header file for class Transform3D
12//
13// Created by: Lorenzo Moneta October 21 2005
14//
15//
16#ifndef ROOT_Math_GenVector_Transform3D
17#define ROOT_Math_GenVector_Transform3D 1
18
19
20
22
24
26
28
29
37
38#include <iostream>
39#include <type_traits>
40#include <cmath>
41
42//#include "Math/Vector3Dfwd.h"
43
44
45
46namespace ROOT {
47
48namespace Math {
49
50namespace Impl {
51
52//_________________________________________________________________________________________
53/**
54 Basic 3D Transformation class describing a rotation and then a translation
55 The internal data are a 3D rotation data (represented as a 3x3 matrix) and a 3D vector data.
56 They are represented and held in this class like a 3x4 matrix (a simple array of 12 numbers).
57
58 The class can be constructed from any 3D rotation object
59 (ROOT::Math::Rotation3D, ROOT::Math::AxisAngle, ROOT::Math::Quaternion, etc...) and/or
60 a 3D Vector (ROOT::Math::DislacementVector3D or via ROOT::Math::Translation ) representing a Translation.
61 The Transformation is defined by applying first the rotation and then the translation.
62 A transformation defined by applying first a translation and then a rotation is equivalent to the
63 transformation obtained applying first the rotation and then a translation equivalent to the rotated vector.
64 The operator * can be used to obtain directly such transformations, in addition to combine various
65 transformations.
66 Keep in mind that the operator * (like in the case of rotations ) is not commutative.
67 The operator * is used (in addition to operator() ) to apply a transformations on the vector
68 (DisplacementVector3D and LorentzVector classes) and point (PositionVector3D) classes.
69 In the case of Vector objects the transformation only rotates them and does not translate them.
70 Only Point objects are able to be both rotated and translated.
71
72
73 @ingroup GenVector
74
75 @sa Overview of the @ref GenVector "physics vector library"
76
77*/
78
79template <typename T = double>
81
82public:
83 typedef T Scalar;
84
87
89 kXX = 0, kXY = 1, kXZ = 2, kDX = 3,
90 kYX = 4, kYY = 5, kYZ = 6, kDY = 7,
91 kZX = 8, kZY = 9, kZZ =10, kDZ = 11
92 };
93
94
95
96 /**
97 Default constructor (identy rotation) + zero translation
98 */
100 {
101 SetIdentity();
102 }
103
104 /**
105 Construct given a pair of pointers or iterators defining the
106 beginning and end of an array of 12 Scalars
107 */
108 template<class IT>
109 Transform3D(IT begin, IT end)
110 {
111 SetComponents(begin,end);
112 }
113
114 /**
115 Construct from a rotation and then a translation described by a Vector
116 */
117 Transform3D( const Rotation3D & r, const Vector & v)
118 {
119 AssignFrom( r, v );
120 }
121 /**
122 Construct from a rotation and then a translation described by a Translation3D class
123 */
125
126 /**
127 Construct from a rotation (any rotation object) and then a translation
128 (represented by any DisplacementVector)
129 The requirements on the rotation and vector objects are that they can be transformed in a
130 Rotation3D class and in a Cartesian3D Vector
131 */
132 template <class ARotation, class CoordSystem, class Tag>
134 {
135 AssignFrom( Rotation3D(r), Vector (v.X(),v.Y(),v.Z()) );
136 }
137
138 /**
139 Construct from a rotation (any rotation object) and then a translation
140 represented by a Translation3D class
141 The requirements on the rotation is that it can be transformed in a
142 Rotation3D class
143 */
144 template <class ARotation>
145 Transform3D(const ARotation &r, const Translation3D<T> &t)
146 {
147 AssignFrom( Rotation3D(r), t.Vect() );
148 }
149
150
151#ifdef OLD_VERSION
152 /**
153 Construct from a translation and then a rotation (inverse assignment)
154 */
155 Transform3D( const Vector & v, const Rotation3D & r)
156 {
157 // is equivalent from having first the rotation and then the translation vector rotated
158 AssignFrom( r, r(v) );
159 }
160#endif
161
162 /**
163 Construct from a 3D Rotation only with zero translation
164 */
165 explicit Transform3D( const Rotation3D & r) {
166 AssignFrom(r);
167 }
168
169 // convenience methods for constructing a Transform3D from all the 3D rotations classes
170 // (cannot use templates for conflict with LA)
171
172 explicit Transform3D( const AxisAngle & r) {
174 }
175 explicit Transform3D( const EulerAngles & r) {
177 }
178 explicit Transform3D( const Quaternion & r) {
180 }
181 explicit Transform3D( const RotationZYX & r) {
183 }
184
185 // Constructors from axial rotations
186 // TO DO: implement direct methods for axial rotations without going through Rotation3D
187 explicit Transform3D( const RotationX & r) {
189 }
190 explicit Transform3D( const RotationY & r) {
192 }
193 explicit Transform3D( const RotationZ & r) {
195 }
196
197 /**
198 Construct from a translation only, represented by any DisplacementVector3D
199 and with an identity rotation
200 */
201 template<class CoordSystem, class Tag>
203 AssignFrom(Vector(v.X(),v.Y(),v.Z()));
204 }
205 /**
206 Construct from a translation only, represented by a Cartesian 3D Vector,
207 and with an identity rotation
208 */
209 explicit Transform3D( const Vector & v) {
210 AssignFrom(v);
211 }
212 /**
213 Construct from a translation only, represented by a Translation3D class
214 and with an identity rotation
215 */
216 explicit Transform3D(const Translation3D<T> &t) { AssignFrom(t.Vect()); }
217
218 //#if !defined(__MAKECINT__) && !defined(G__DICTIONARY) // this is ambiguous with double * , double *
219
220
221#ifdef OLD_VERSION
222 /**
223 Construct from a translation (using any type of DisplacementVector )
224 and then a rotation (any rotation object).
225 Requirement on the rotation and vector objects are that they can be transformed in a
226 Rotation3D class and in a Vector
227 */
228 template <class ARotation, class CoordSystem, class Tag>
229 Transform3D(const DisplacementVector3D<CoordSystem,Tag> & v , const ARotation & r)
230 {
231 // is equivalent from having first the rotation and then the translation vector rotated
232 Rotation3D r3d(r);
233 AssignFrom( r3d, r3d( Vector(v.X(),v.Y(),v.Z()) ) );
234 }
235#endif
236
237public:
238 /**
239 Construct transformation from one coordinate system defined by three
240 points (origin + two axis) to
241 a new coordinate system defined by other three points (origin + axis)
242 Scalar version.
243 @param fr0 point defining origin of original reference system
244 @param fr1 point defining first axis of original reference system
245 @param fr2 point defining second axis of original reference system
246 @param to0 point defining origin of transformed reference system
247 @param to1 point defining first axis transformed reference system
248 @param to2 point defining second axis transformed reference system
249 */
250 template <typename SCALAR = T, typename std::enable_if<std::is_arithmetic<SCALAR>::value>::type * = nullptr>
251 Transform3D(const Point &fr0, const Point &fr1, const Point &fr2, const Point &to0, const Point &to1,
252 const Point &to2)
253 {
254 // takes impl. from CLHEP ( E.Chernyaev). To be checked
255
256 Vector x1 = (fr1 - fr0).Unit();
257 Vector y1 = (fr2 - fr0).Unit();
258 Vector x2 = (to1 - to0).Unit();
259 Vector y2 = (to2 - to0).Unit();
260
261 // C H E C K A N G L E S
262
263 const T cos1 = x1.Dot(y1);
264 const T cos2 = x2.Dot(y2);
265
266 if (std::fabs(T(1) - cos1) <= T(0.000001) || std::fabs(T(1) - cos2) <= T(0.000001)) {
267 std::cerr << "Transform3D: Error : zero angle between axes" << std::endl;
268 SetIdentity();
269 } else {
270 if (std::fabs(cos1 - cos2) > T(0.000001)) {
271 std::cerr << "Transform3D: Warning: angles between axes are not equal" << std::endl;
272 }
273
274 // F I N D R O T A T I O N M A T R I X
275
276 Vector z1 = (x1.Cross(y1)).Unit();
277 y1 = z1.Cross(x1);
278
279 Vector z2 = (x2.Cross(y2)).Unit();
280 y2 = z2.Cross(x2);
281
282 T x1x = x1.x();
283 T x1y = x1.y();
284 T x1z = x1.z();
285 T y1x = y1.x();
286 T y1y = y1.y();
287 T y1z = y1.z();
288 T z1x = z1.x();
289 T z1y = z1.y();
290 T z1z = z1.z();
291
292 T x2x = x2.x();
293 T x2y = x2.y();
294 T x2z = x2.z();
295 T y2x = y2.x();
296 T y2y = y2.y();
297 T y2z = y2.z();
298 T z2x = z2.x();
299 T z2y = z2.y();
300 T z2z = z2.z();
301
302 T detxx = (y1y * z1z - z1y * y1z);
303 T detxy = -(y1x * z1z - z1x * y1z);
304 T detxz = (y1x * z1y - z1x * y1y);
305 T detyx = -(x1y * z1z - z1y * x1z);
306 T detyy = (x1x * z1z - z1x * x1z);
307 T detyz = -(x1x * z1y - z1x * x1y);
308 T detzx = (x1y * y1z - y1y * x1z);
309 T detzy = -(x1x * y1z - y1x * x1z);
310 T detzz = (x1x * y1y - y1x * x1y);
311
312 T txx = x2x * detxx + y2x * detyx + z2x * detzx;
313 T txy = x2x * detxy + y2x * detyy + z2x * detzy;
314 T txz = x2x * detxz + y2x * detyz + z2x * detzz;
315 T tyx = x2y * detxx + y2y * detyx + z2y * detzx;
316 T tyy = x2y * detxy + y2y * detyy + z2y * detzy;
317 T tyz = x2y * detxz + y2y * detyz + z2y * detzz;
318 T tzx = x2z * detxx + y2z * detyx + z2z * detzx;
319 T tzy = x2z * detxy + y2z * detyy + z2z * detzy;
320 T tzz = x2z * detxz + y2z * detyz + z2z * detzz;
321
322 // S E T T R A N S F O R M A T I O N
323
324 T dx1 = fr0.x(), dy1 = fr0.y(), dz1 = fr0.z();
325 T dx2 = to0.x(), dy2 = to0.y(), dz2 = to0.z();
326
327 SetComponents(txx, txy, txz, dx2 - txx * dx1 - txy * dy1 - txz * dz1, tyx, tyy, tyz,
328 dy2 - tyx * dx1 - tyy * dy1 - tyz * dz1, tzx, tzy, tzz, dz2 - tzx * dx1 - tzy * dy1 - tzz * dz1);
329 }
330 }
331
332 /**
333 Construct transformation from one coordinate system defined by three
334 points (origin + two axis) to
335 a new coordinate system defined by other three points (origin + axis)
336 Vectorised version.
337 @param fr0 point defining origin of original reference system
338 @param fr1 point defining first axis of original reference system
339 @param fr2 point defining second axis of original reference system
340 @param to0 point defining origin of transformed reference system
341 @param to1 point defining first axis transformed reference system
342 @param to2 point defining second axis transformed reference system
343 */
344 template <typename SCALAR = T, typename std::enable_if<!std::is_arithmetic<SCALAR>::value>::type * = nullptr>
345 Transform3D(const Point &fr0, const Point &fr1, const Point &fr2, const Point &to0, const Point &to1,
346 const Point &to2)
347 {
348 // takes impl. from CLHEP ( E.Chernyaev). To be checked
349
350 Vector x1 = (fr1 - fr0).Unit();
351 Vector y1 = (fr2 - fr0).Unit();
352 Vector x2 = (to1 - to0).Unit();
353 Vector y2 = (to2 - to0).Unit();
354
355 // C H E C K A N G L E S
356
357 const T cos1 = x1.Dot(y1);
358 const T cos2 = x2.Dot(y2);
359
360 const auto m1 = (abs(T(1) - cos1) <= T(0.000001) || abs(T(1) - cos2) <= T(0.000001));
361
362 const auto m2 = (abs(cos1 - cos2) > T(0.000001));
363 if (any_of(m2)) {
364 std::cerr << "Transform3D: Warning: angles between axes are not equal" << std::endl;
365 }
366
367 // F I N D R O T A T I O N M A T R I X
368
369 Vector z1 = (x1.Cross(y1)).Unit();
370 y1 = z1.Cross(x1);
371
372 Vector z2 = (x2.Cross(y2)).Unit();
373 y2 = z2.Cross(x2);
374
375 T x1x = x1.x();
376 T x1y = x1.y();
377 T x1z = x1.z();
378 T y1x = y1.x();
379 T y1y = y1.y();
380 T y1z = y1.z();
381 T z1x = z1.x();
382 T z1y = z1.y();
383 T z1z = z1.z();
384
385 T x2x = x2.x();
386 T x2y = x2.y();
387 T x2z = x2.z();
388 T y2x = y2.x();
389 T y2y = y2.y();
390 T y2z = y2.z();
391 T z2x = z2.x();
392 T z2y = z2.y();
393 T z2z = z2.z();
394
395 T detxx = (y1y * z1z - z1y * y1z);
396 T detxy = -(y1x * z1z - z1x * y1z);
397 T detxz = (y1x * z1y - z1x * y1y);
398 T detyx = -(x1y * z1z - z1y * x1z);
399 T detyy = (x1x * z1z - z1x * x1z);
400 T detyz = -(x1x * z1y - z1x * x1y);
401 T detzx = (x1y * y1z - y1y * x1z);
402 T detzy = -(x1x * y1z - y1x * x1z);
403 T detzz = (x1x * y1y - y1x * x1y);
404
405 T txx = x2x * detxx + y2x * detyx + z2x * detzx;
406 T txy = x2x * detxy + y2x * detyy + z2x * detzy;
407 T txz = x2x * detxz + y2x * detyz + z2x * detzz;
408 T tyx = x2y * detxx + y2y * detyx + z2y * detzx;
409 T tyy = x2y * detxy + y2y * detyy + z2y * detzy;
410 T tyz = x2y * detxz + y2y * detyz + z2y * detzz;
411 T tzx = x2z * detxx + y2z * detyx + z2z * detzx;
412 T tzy = x2z * detxy + y2z * detyy + z2z * detzy;
413 T tzz = x2z * detxz + y2z * detyz + z2z * detzz;
414
415 // S E T T R A N S F O R M A T I O N
416
417 T dx1 = fr0.x(), dy1 = fr0.y(), dz1 = fr0.z();
418 T dx2 = to0.x(), dy2 = to0.y(), dz2 = to0.z();
419
420 SetComponents(txx, txy, txz, dx2 - txx * dx1 - txy * dy1 - txz * dz1, tyx, tyy, tyz,
421 dy2 - tyx * dx1 - tyy * dy1 - tyz * dz1, tzx, tzy, tzz, dz2 - tzx * dx1 - tzy * dy1 - tzz * dz1);
422
423 if (any_of(m1)) {
424 std::cerr << "Transform3D: Error : zero angle between axes" << std::endl;
425 SetIdentity(m1);
426 }
427 }
428
429 // use compiler generated copy ctor, copy assignment and destructor
430
431 /**
432 Construct from a linear algebra matrix of size at least 3x4,
433 which must support operator()(i,j) to obtain elements (0,0) thru (2,3).
434 The 3x3 sub-block is assumed to be the rotation part and the translations vector
435 are described by the 4-th column
436 */
437 template<class ForeignMatrix>
438 explicit Transform3D(const ForeignMatrix & m) {
440 }
441
442 /**
443 Raw constructor from 12 Scalar components
444 */
445 Transform3D(T xx, T xy, T xz, T dx, T yx, T yy, T yz, T dy, T zx, T zy, T zz, T dz)
446 {
447 SetComponents (xx, xy, xz, dx, yx, yy, yz, dy, zx, zy, zz, dz);
448 }
449
450
451 /**
452 Construct from a linear algebra matrix of size at least 3x4,
453 which must support operator()(i,j) to obtain elements (0,0) thru (2,3).
454 The 3x3 sub-block is assumed to be the rotation part and the translations vector
455 are described by the 4-th column
456 */
457 template <class ForeignMatrix>
458 Transform3D<T> &operator=(const ForeignMatrix &m)
459 {
461 return *this;
462 }
463
464
465 // ======== Components ==============
466
467
468 /**
469 Set the 12 matrix components given an iterator to the start of
470 the desired data, and another to the end (12 past start).
471 */
472 template<class IT>
473 void SetComponents(IT begin, IT end) {
474 for (int i = 0; i <12; ++i) {
475 fM[i] = *begin;
476 ++begin;
477 }
478 (void)end;
479 assert (end==begin);
480 }
481
482 /**
483 Get the 12 matrix components into data specified by an iterator begin
484 and another to the end of the desired data (12 past start).
485 */
486 template<class IT>
487 void GetComponents(IT begin, IT end) const {
488 for (int i = 0; i <12; ++i) {
489 *begin = fM[i];
490 ++begin;
491 }
492 (void)end;
493 assert (end==begin);
494 }
495
496 /**
497 Get the 12 matrix components into data specified by an iterator begin
498 */
499 template<class IT>
500 void GetComponents(IT begin) const {
501 std::copy(fM, fM + 12, begin);
502 }
503
504 /**
505 Set components from a linear algebra matrix of size at least 3x4,
506 which must support operator()(i,j) to obtain elements (0,0) thru (2,3).
507 The 3x3 sub-block is assumed to be the rotation part and the translations vector
508 are described by the 4-th column
509 */
510 template<class ForeignMatrix>
511 void
512 SetTransformMatrix (const ForeignMatrix & m) {
513 fM[kXX]=m(0,0); fM[kXY]=m(0,1); fM[kXZ]=m(0,2); fM[kDX]=m(0,3);
514 fM[kYX]=m(1,0); fM[kYY]=m(1,1); fM[kYZ]=m(1,2); fM[kDY]=m(1,3);
515 fM[kZX]=m(2,0); fM[kZY]=m(2,1); fM[kZZ]=m(2,2); fM[kDZ]=m(2,3);
516 }
517
518 /**
519 Get components into a linear algebra matrix of size at least 3x4,
520 which must support operator()(i,j) for write access to elements
521 (0,0) thru (2,3).
522 */
523 template<class ForeignMatrix>
524 void
525 GetTransformMatrix (ForeignMatrix & m) const {
526 m(0,0)=fM[kXX]; m(0,1)=fM[kXY]; m(0,2)=fM[kXZ]; m(0,3)=fM[kDX];
527 m(1,0)=fM[kYX]; m(1,1)=fM[kYY]; m(1,2)=fM[kYZ]; m(1,3)=fM[kDY];
528 m(2,0)=fM[kZX]; m(2,1)=fM[kZY]; m(2,2)=fM[kZZ]; m(2,3)=fM[kDZ];
529 }
530
531
532 /**
533 Set the components from 12 scalars
534 */
535 void SetComponents(T xx, T xy, T xz, T dx, T yx, T yy, T yz, T dy, T zx, T zy, T zz, T dz)
536 {
537 fM[kXX]=xx; fM[kXY]=xy; fM[kXZ]=xz; fM[kDX]=dx;
538 fM[kYX]=yx; fM[kYY]=yy; fM[kYZ]=yz; fM[kDY]=dy;
539 fM[kZX]=zx; fM[kZY]=zy; fM[kZZ]=zz; fM[kDZ]=dz;
540 }
541
542 /**
543 Get the components into 12 scalars
544 */
545 void GetComponents(T &xx, T &xy, T &xz, T &dx, T &yx, T &yy, T &yz, T &dy, T &zx, T &zy, T &zz, T &dz) const
546 {
547 xx=fM[kXX]; xy=fM[kXY]; xz=fM[kXZ]; dx=fM[kDX];
548 yx=fM[kYX]; yy=fM[kYY]; yz=fM[kYZ]; dy=fM[kDY];
549 zx=fM[kZX]; zy=fM[kZY]; zz=fM[kZZ]; dz=fM[kDZ];
550 }
551
552
553 /**
554 Get the rotation and translation vector representing the 3D transformation
555 in any rotation and any vector (the Translation class could also be used)
556 */
557 template<class AnyRotation, class V>
558 void GetDecomposition(AnyRotation &r, V &v) const {
559 GetRotation(r);
561 }
562
563
564 /**
565 Get the rotation and translation vector representing the 3D transformation
566 */
568 GetRotation(r);
570 }
571
572 /**
573 Get the 3D rotation representing the 3D transformation
574 */
576 return Rotation3D( fM[kXX], fM[kXY], fM[kXZ],
577 fM[kYX], fM[kYY], fM[kYZ],
578 fM[kZX], fM[kZY], fM[kZZ] );
579 }
580
581 /**
582 Get the rotation representing the 3D transformation
583 */
584 template <class AnyRotation>
585 AnyRotation Rotation() const {
586 return AnyRotation(Rotation3D(fM[kXX], fM[kXY], fM[kXZ], fM[kYX], fM[kYY], fM[kYZ], fM[kZX], fM[kZY], fM[kZZ]));
587 }
588
589 /**
590 Get the rotation (any type) representing the 3D transformation
591 */
592 template <class AnyRotation>
593 void GetRotation(AnyRotation &r) const {
594 r = Rotation();
595 }
596
597 /**
598 Get the translation representing the 3D transformation in a Cartesian vector
599 */
601
602 /**
603 Get the translation representing the 3D transformation in any vector
604 which implements the SetXYZ method
605 */
606 template <class AnyVector>
607 void GetTranslation(AnyVector &v) const {
608 v.SetXYZ(fM[kDX], fM[kDY], fM[kDZ]);
609 }
610
611
612
613 // operations on points and vectors
614
615 /**
616 Transformation operation for Position Vector in Cartesian coordinate
617 For a Position Vector first a rotation and then a translation is applied
618 */
619 Point operator() (const Point & p) const {
620 return Point ( fM[kXX]*p.X() + fM[kXY]*p.Y() + fM[kXZ]*p.Z() + fM[kDX],
621 fM[kYX]*p.X() + fM[kYY]*p.Y() + fM[kYZ]*p.Z() + fM[kDY],
622 fM[kZX]*p.X() + fM[kZY]*p.Y() + fM[kZZ]*p.Z() + fM[kDZ] );
623 }
624
625
626 /**
627 Transformation operation for Displacement Vectors in Cartesian coordinate
628 For the Displacement Vectors only the rotation applies - no translations
629 */
630 Vector operator() (const Vector & v) const {
631 return Vector( fM[kXX]*v.X() + fM[kXY]*v.Y() + fM[kXZ]*v.Z() ,
632 fM[kYX]*v.X() + fM[kYY]*v.Y() + fM[kYZ]*v.Z() ,
633 fM[kZX]*v.X() + fM[kZY]*v.Y() + fM[kZZ]*v.Z() );
634 }
635
636
637 /**
638 Transformation operation for Position Vector in any coordinate system
639 */
640 template <class CoordSystem>
642 {
643 return PositionVector3D<CoordSystem>(operator()(Point(p)));
644 }
645 /**
646 Transformation operation for Position Vector in any coordinate system
647 */
648 template <class CoordSystem>
650 {
651 return operator()(v);
652 }
653
654 /**
655 Transformation operation for Displacement Vector in any coordinate system
656 */
657 template<class CoordSystem >
658 DisplacementVector3D<CoordSystem> operator() (const DisplacementVector3D <CoordSystem> & v) const {
659 return DisplacementVector3D<CoordSystem>(operator()(Vector(v)));
660 }
661 /**
662 Transformation operation for Displacement Vector in any coordinate system
663 */
664 template <class CoordSystem>
666 {
667 return operator()(v);
668 }
669
670 /**
671 Directly apply the inverse affine transformation on vectors.
672 Avoids having to calculate the inverse as an intermediate result.
673 This is possible since the inverse of a rotation is its transpose.
674 */
676 {
677 return Vector(fM[kXX] * v.X() + fM[kYX] * v.Y() + fM[kZX] * v.Z(),
678 fM[kXY] * v.X() + fM[kYY] * v.Y() + fM[kZY] * v.Z(),
679 fM[kXZ] * v.X() + fM[kYZ] * v.Y() + fM[kZZ] * v.Z());
680 }
681
682 /**
683 Directly apply the inverse affine transformation on points
684 (first inverse translation then inverse rotation).
685 Avoids having to calculate the inverse as an intermediate result.
686 This is possible since the inverse of a rotation is its transpose.
687 */
688 Point ApplyInverse(const Point &p) const
689 {
690 Point tmp(p.X() - fM[kDX], p.Y() - fM[kDY], p.Z() - fM[kDZ]);
691 return Point(fM[kXX] * tmp.X() + fM[kYX] * tmp.Y() + fM[kZX] * tmp.Z(),
692 fM[kXY] * tmp.X() + fM[kYY] * tmp.Y() + fM[kZY] * tmp.Z(),
693 fM[kXZ] * tmp.X() + fM[kYZ] * tmp.Y() + fM[kZZ] * tmp.Z());
694 }
695
696 /**
697 Directly apply the inverse affine transformation on an arbitrary
698 coordinate-system point.
699 Involves casting to Point(p) type.
700 */
701 template <class CoordSystem>
703 {
705 }
706
707 /**
708 Directly apply the inverse affine transformation on an arbitrary
709 coordinate-system vector.
710 Involves casting to Vector(p) type.
711 */
712 template <class CoordSystem>
714 {
716 }
717
718 /**
719 Transformation operation for points between different coordinate system tags
720 */
721 template <class CoordSystem, class Tag1, class Tag2>
723 {
724 const Point xyzNew = operator()(Point(p1.X(), p1.Y(), p1.Z()));
725 p2.SetXYZ( xyzNew.X(), xyzNew.Y(), xyzNew.Z() );
726 }
727
728
729 /**
730 Transformation operation for Displacement Vector of different coordinate systems
731 */
732 template <class CoordSystem, class Tag1, class Tag2>
734 {
735 const Vector xyzNew = operator()(Vector(v1.X(), v1.Y(), v1.Z()));
736 v2.SetXYZ( xyzNew.X(), xyzNew.Y(), xyzNew.Z() );
737 }
738
739 /**
740 Transformation operation for a Lorentz Vector in any coordinate system
741 */
742 template <class CoordSystem >
744 const Vector xyzNew = operator()(Vector(q.Vect()));
745 return LorentzVector<CoordSystem>(xyzNew.X(), xyzNew.Y(), xyzNew.Z(), q.E());
746 }
747 /**
748 Transformation operation for a Lorentz Vector in any coordinate system
749 */
750 template <class CoordSystem>
752 {
753 return operator()(q);
754 }
755
756 /**
757 Transformation on a 3D plane
758 */
759 template <typename TYPE>
761 {
762 // transformations on a 3D plane
763 const auto n = plane.Normal();
764 // take a point on the plane. Use origin projection on the plane
765 // ( -ad, -bd, -cd) if (a**2 + b**2 + c**2 ) = 1
766 const auto d = plane.HesseDistance();
767 Point p(-d * n.X(), -d * n.Y(), -d * n.Z());
768 return Plane3D<TYPE>(operator()(n), operator()(p));
769 }
770
771 /// Multiplication operator for 3D plane
772 template <typename TYPE>
774 {
775 return operator()(plane);
776 }
777
778 // skip transformation for arbitrary vectors - not really defined if point or displacement vectors
779
780 /**
781 multiply (combine) with another transformation in place
782 */
783 inline Transform3D<T> &operator*=(const Transform3D<T> &t);
784
785 /**
786 multiply (combine) two transformations
787 */
788 inline Transform3D<T> operator*(const Transform3D<T> &t) const;
789
790 /**
791 Invert the transformation in place (scalar)
792 */
793 template <typename SCALAR = T, typename std::enable_if<std::is_arithmetic<SCALAR>::value>::type * = nullptr>
794 void Invert()
795 {
796 //
797 // Name: Transform3D::inverse Date: 24.09.96
798 // Author: E.Chernyaev (IHEP/Protvino) Revised:
799 //
800 // Function: Find inverse affine transformation.
801
802 T detxx = fM[kYY] * fM[kZZ] - fM[kYZ] * fM[kZY];
803 T detxy = fM[kYX] * fM[kZZ] - fM[kYZ] * fM[kZX];
804 T detxz = fM[kYX] * fM[kZY] - fM[kYY] * fM[kZX];
805 T det = fM[kXX] * detxx - fM[kXY] * detxy + fM[kXZ] * detxz;
806 if (det == T(0)) {
807 std::cerr << "Transform3D::inverse error: zero determinant" << std::endl;
808 return;
809 }
810 det = T(1) / det;
811 detxx *= det;
812 detxy *= det;
813 detxz *= det;
814 T detyx = (fM[kXY] * fM[kZZ] - fM[kXZ] * fM[kZY]) * det;
815 T detyy = (fM[kXX] * fM[kZZ] - fM[kXZ] * fM[kZX]) * det;
816 T detyz = (fM[kXX] * fM[kZY] - fM[kXY] * fM[kZX]) * det;
817 T detzx = (fM[kXY] * fM[kYZ] - fM[kXZ] * fM[kYY]) * det;
818 T detzy = (fM[kXX] * fM[kYZ] - fM[kXZ] * fM[kYX]) * det;
819 T detzz = (fM[kXX] * fM[kYY] - fM[kXY] * fM[kYX]) * det;
820 SetComponents(detxx, -detyx, detzx, -detxx * fM[kDX] + detyx * fM[kDY] - detzx * fM[kDZ], -detxy, detyy, -detzy,
821 detxy * fM[kDX] - detyy * fM[kDY] + detzy * fM[kDZ], detxz, -detyz, detzz,
822 -detxz * fM[kDX] + detyz * fM[kDY] - detzz * fM[kDZ]);
823 }
824
825 /**
826 Invert the transformation in place (vectorised)
827 */
828 template <typename SCALAR = T, typename std::enable_if<!std::is_arithmetic<SCALAR>::value>::type * = nullptr>
829 void Invert()
830 {
831 //
832 // Name: Transform3D::inverse Date: 24.09.96
833 // Author: E.Chernyaev (IHEP/Protvino) Revised:
834 //
835 // Function: Find inverse affine transformation.
836
837 T detxx = fM[kYY] * fM[kZZ] - fM[kYZ] * fM[kZY];
838 T detxy = fM[kYX] * fM[kZZ] - fM[kYZ] * fM[kZX];
839 T detxz = fM[kYX] * fM[kZY] - fM[kYY] * fM[kZX];
840 T det = fM[kXX] * detxx - fM[kXY] * detxy + fM[kXZ] * detxz;
841 const auto detZmask = (det == T(0));
842 if (any_of(detZmask)) {
843 std::cerr << "Transform3D::inverse error: zero determinant" << std::endl;
844 det(detZmask) = T(1);
845 }
846 det = T(1) / det;
847 detxx *= det;
848 detxy *= det;
849 detxz *= det;
850 T detyx = (fM[kXY] * fM[kZZ] - fM[kXZ] * fM[kZY]) * det;
851 T detyy = (fM[kXX] * fM[kZZ] - fM[kXZ] * fM[kZX]) * det;
852 T detyz = (fM[kXX] * fM[kZY] - fM[kXY] * fM[kZX]) * det;
853 T detzx = (fM[kXY] * fM[kYZ] - fM[kXZ] * fM[kYY]) * det;
854 T detzy = (fM[kXX] * fM[kYZ] - fM[kXZ] * fM[kYX]) * det;
855 T detzz = (fM[kXX] * fM[kYY] - fM[kXY] * fM[kYX]) * det;
856 // Set det=0 cases to 0
857 if (any_of(detZmask)) {
858 detxx(detZmask) = T(0);
859 detxy(detZmask) = T(0);
860 detxz(detZmask) = T(0);
861 detyx(detZmask) = T(0);
862 detyy(detZmask) = T(0);
863 detyz(detZmask) = T(0);
864 detzx(detZmask) = T(0);
865 detzy(detZmask) = T(0);
866 detzz(detZmask) = T(0);
867 }
868 // set final components
869 SetComponents(detxx, -detyx, detzx, -detxx * fM[kDX] + detyx * fM[kDY] - detzx * fM[kDZ], -detxy, detyy, -detzy,
870 detxy * fM[kDX] - detyy * fM[kDY] + detzy * fM[kDZ], detxz, -detyz, detzz,
871 -detxz * fM[kDX] + detyz * fM[kDY] - detzz * fM[kDZ]);
872 }
873
874 /**
875 Return the inverse of the transformation.
876 */
878 {
879 Transform3D<T> t(*this);
880 t.Invert();
881 return t;
882 }
883
884 /**
885 Equality operator. Check equality for each element
886 To do: use T tolerance
887 */
888 bool operator==(const Transform3D<T> &rhs) const
889 {
890 return (fM[0] == rhs.fM[0] && fM[1] == rhs.fM[1] && fM[2] == rhs.fM[2] && fM[3] == rhs.fM[3] &&
891 fM[4] == rhs.fM[4] && fM[5] == rhs.fM[5] && fM[6] == rhs.fM[6] && fM[7] == rhs.fM[7] &&
892 fM[8] == rhs.fM[8] && fM[9] == rhs.fM[9] && fM[10] == rhs.fM[10] && fM[11] == rhs.fM[11]);
893 }
894
895 /**
896 Inequality operator. Check equality for each element
897 To do: use T tolerance
898 */
899 bool operator!=(const Transform3D<T> &rhs) const { return !operator==(rhs); }
900
901protected:
902
903 /**
904 make transformation from first a rotation then a translation
905 */
906 void AssignFrom(const Rotation3D &r, const Vector &v)
907 {
908 // assignment from rotation + translation
909
910 T rotData[9];
911 r.GetComponents(rotData, rotData + 9);
912 // first raw
913 for (int i = 0; i < 3; ++i) fM[i] = rotData[i];
914 // second raw
915 for (int i = 0; i < 3; ++i) fM[kYX + i] = rotData[3 + i];
916 // third raw
917 for (int i = 0; i < 3; ++i) fM[kZX + i] = rotData[6 + i];
918
919 // translation data
920 T vecData[3];
921 v.GetCoordinates(vecData, vecData + 3);
922 fM[kDX] = vecData[0];
923 fM[kDY] = vecData[1];
924 fM[kDZ] = vecData[2];
925 }
926
927 /**
928 make transformation from only rotations (zero translation)
929 */
931 {
932 // assign from only a rotation (null translation)
933 T rotData[9];
934 r.GetComponents(rotData, rotData + 9);
935 for (int i = 0; i < 3; ++i) {
936 for (int j = 0; j < 3; ++j) fM[4 * i + j] = rotData[3 * i + j];
937 // empty vector data
938 fM[4 * i + 3] = T(0);
939 }
940 }
941
942 /**
943 make transformation from only translation (identity rotations)
944 */
945 void AssignFrom(const Vector &v)
946 {
947 // assign from a translation only (identity rotations)
948 fM[kXX] = T(1);
949 fM[kXY] = T(0);
950 fM[kXZ] = T(0);
951 fM[kDX] = v.X();
952 fM[kYX] = T(0);
953 fM[kYY] = T(1);
954 fM[kYZ] = T(0);
955 fM[kDY] = v.Y();
956 fM[kZX] = T(0);
957 fM[kZY] = T(0);
958 fM[kZZ] = T(1);
959 fM[kDZ] = v.Z();
960 }
961
962 /**
963 Set identity transformation (identity rotation , zero translation)
964 */
966 {
967 // set identity ( identity rotation and zero translation)
968 fM[kXX] = T(1);
969 fM[kXY] = T(0);
970 fM[kXZ] = T(0);
971 fM[kDX] = T(0);
972 fM[kYX] = T(0);
973 fM[kYY] = T(1);
974 fM[kYZ] = T(0);
975 fM[kDY] = T(0);
976 fM[kZX] = T(0);
977 fM[kZY] = T(0);
978 fM[kZZ] = T(1);
979 fM[kDZ] = T(0);
980 }
981
982 /**
983 Set identity transformation (identity rotation , zero translation)
984 vectorised version that sets using a mask
985 */
986 template <typename SCALAR = T, typename std::enable_if<!std::is_arithmetic<SCALAR>::value>::type * = nullptr>
987 void SetIdentity(const typename SCALAR::mask_type m)
988 {
989 // set identity ( identity rotation and zero translation)
990 fM[kXX](m) = T(1);
991 fM[kXY](m) = T(0);
992 fM[kXZ](m) = T(0);
993 fM[kDX](m) = T(0);
994 fM[kYX](m) = T(0);
995 fM[kYY](m) = T(1);
996 fM[kYZ](m) = T(0);
997 fM[kDY](m) = T(0);
998 fM[kZX](m) = T(0);
999 fM[kZY](m) = T(0);
1000 fM[kZZ](m) = T(1);
1001 fM[kDZ](m) = T(0);
1002 }
1003
1004private:
1005 T fM[12]; // transformation elements (3x4 matrix)
1006};
1007
1008
1009
1010
1011// inline functions (combination of transformations)
1012
1013template <class T>
1015{
1016 // combination of transformations
1017
1018 SetComponents(fM[kXX]*t.fM[kXX]+fM[kXY]*t.fM[kYX]+fM[kXZ]*t.fM[kZX],
1019 fM[kXX]*t.fM[kXY]+fM[kXY]*t.fM[kYY]+fM[kXZ]*t.fM[kZY],
1020 fM[kXX]*t.fM[kXZ]+fM[kXY]*t.fM[kYZ]+fM[kXZ]*t.fM[kZZ],
1021 fM[kXX]*t.fM[kDX]+fM[kXY]*t.fM[kDY]+fM[kXZ]*t.fM[kDZ]+fM[kDX],
1022
1023 fM[kYX]*t.fM[kXX]+fM[kYY]*t.fM[kYX]+fM[kYZ]*t.fM[kZX],
1024 fM[kYX]*t.fM[kXY]+fM[kYY]*t.fM[kYY]+fM[kYZ]*t.fM[kZY],
1025 fM[kYX]*t.fM[kXZ]+fM[kYY]*t.fM[kYZ]+fM[kYZ]*t.fM[kZZ],
1026 fM[kYX]*t.fM[kDX]+fM[kYY]*t.fM[kDY]+fM[kYZ]*t.fM[kDZ]+fM[kDY],
1027
1028 fM[kZX]*t.fM[kXX]+fM[kZY]*t.fM[kYX]+fM[kZZ]*t.fM[kZX],
1029 fM[kZX]*t.fM[kXY]+fM[kZY]*t.fM[kYY]+fM[kZZ]*t.fM[kZY],
1030 fM[kZX]*t.fM[kXZ]+fM[kZY]*t.fM[kYZ]+fM[kZZ]*t.fM[kZZ],
1031 fM[kZX]*t.fM[kDX]+fM[kZY]*t.fM[kDY]+fM[kZZ]*t.fM[kDZ]+fM[kDZ]);
1032
1033 return *this;
1034}
1035
1036template <class T>
1038{
1039 // combination of transformations
1040
1041 return Transform3D<T>(fM[kXX] * t.fM[kXX] + fM[kXY] * t.fM[kYX] + fM[kXZ] * t.fM[kZX],
1042 fM[kXX] * t.fM[kXY] + fM[kXY] * t.fM[kYY] + fM[kXZ] * t.fM[kZY],
1043 fM[kXX] * t.fM[kXZ] + fM[kXY] * t.fM[kYZ] + fM[kXZ] * t.fM[kZZ],
1044 fM[kXX] * t.fM[kDX] + fM[kXY] * t.fM[kDY] + fM[kXZ] * t.fM[kDZ] + fM[kDX],
1045
1046 fM[kYX] * t.fM[kXX] + fM[kYY] * t.fM[kYX] + fM[kYZ] * t.fM[kZX],
1047 fM[kYX] * t.fM[kXY] + fM[kYY] * t.fM[kYY] + fM[kYZ] * t.fM[kZY],
1048 fM[kYX] * t.fM[kXZ] + fM[kYY] * t.fM[kYZ] + fM[kYZ] * t.fM[kZZ],
1049 fM[kYX] * t.fM[kDX] + fM[kYY] * t.fM[kDY] + fM[kYZ] * t.fM[kDZ] + fM[kDY],
1050
1051 fM[kZX] * t.fM[kXX] + fM[kZY] * t.fM[kYX] + fM[kZZ] * t.fM[kZX],
1052 fM[kZX] * t.fM[kXY] + fM[kZY] * t.fM[kYY] + fM[kZZ] * t.fM[kZY],
1053 fM[kZX] * t.fM[kXZ] + fM[kZY] * t.fM[kYZ] + fM[kZZ] * t.fM[kZZ],
1054 fM[kZX] * t.fM[kDX] + fM[kZY] * t.fM[kDY] + fM[kZZ] * t.fM[kDZ] + fM[kDZ]);
1055}
1056
1057
1058
1059
1060//--- global functions resulting in Transform3D -------
1061
1062
1063// ------ combination of a translation (first) and a rotation ------
1064
1065
1066/**
1067 combine a translation and a rotation to give a transform3d
1068 First the translation then the rotation
1069 */
1070template <class T>
1072{
1073 return Transform3D<T>(r, r(t.Vect()));
1074}
1075template <class T>
1077{
1078 Rotation3D r3(r);
1079 return Transform3D<T>(r3, r3(t.Vect()));
1080}
1081template <class T>
1083{
1084 Rotation3D r3(r);
1085 return Transform3D<T>(r3, r3(t.Vect()));
1086}
1087template <class T>
1089{
1090 Rotation3D r3(r);
1091 return Transform3D<T>(r3, r3(t.Vect()));
1092}
1093template <class T>
1095{
1096 Rotation3D r3(r);
1097 return Transform3D<T>(r3, r3(t.Vect()));
1098}
1099template <class T>
1101{
1102 Rotation3D r3(r);
1103 return Transform3D<T>(r3, r3(t.Vect()));
1104}
1105template <class T>
1107{
1108 Rotation3D r3(r);
1109 return Transform3D<T>(r3, r3(t.Vect()));
1110}
1111template <class T>
1113{
1114 Rotation3D r3(r);
1115 return Transform3D<T>(r3, r3(t.Vect()));
1116}
1117
1118// ------ combination of a rotation (first) and then a translation ------
1119
1120/**
1121 combine a rotation and a translation to give a transform3d
1122 First a rotation then the translation
1123 */
1124template <class T>
1126{
1127 return Transform3D<T>(r, t.Vect());
1128}
1129template <class T>
1131{
1132 return Transform3D<T>(Rotation3D(r), t.Vect());
1133}
1134template <class T>
1136{
1137 return Transform3D<T>(Rotation3D(r), t.Vect());
1138}
1139template <class T>
1141{
1142 return Transform3D<T>(Rotation3D(r), t.Vect());
1143}
1144template <class T>
1146{
1147 return Transform3D<T>(Rotation3D(r), t.Vect());
1148}
1149template <class T>
1151{
1152 return Transform3D<T>(Rotation3D(r), t.Vect());
1153}
1154template <class T>
1156{
1157 return Transform3D<T>(Rotation3D(r), t.Vect());
1158}
1159template <class T>
1161{
1162 return Transform3D<T>(Rotation3D(r), t.Vect());
1163}
1164
1165// ------ combination of a Transform3D and a pure translation------
1166
1167/**
1168 combine a transformation and a translation to give a transform3d
1169 First the translation then the transform3D
1170 */
1171template <class T>
1173{
1174 Rotation3D r = t.Rotation();
1175 return Transform3D<T>(r, r(d.Vect()) + t.Translation().Vect());
1176}
1177
1178/**
1179 combine a translation and a transformation to give a transform3d
1180 First the transformation then the translation
1181 */
1182template <class T>
1184{
1185 return Transform3D<T>(t.Rotation(), t.Translation().Vect() + d.Vect());
1186}
1187
1188// ------ combination of a Transform3D and any rotation------
1189
1190
1191/**
1192 combine a transformation and a rotation to give a transform3d
1193 First the rotation then the transform3D
1194 */
1195template <class T>
1197{
1198 return Transform3D<T>(t.Rotation() * r, t.Translation());
1199}
1200template <class T>
1202{
1203 return Transform3D<T>(t.Rotation() * r, t.Translation());
1204}
1205template <class T>
1207{
1208 return Transform3D<T>(t.Rotation() * r, t.Translation());
1209}
1210template <class T>
1212{
1213 return Transform3D<T>(t.Rotation() * r, t.Translation());
1214}
1215template <class T>
1217{
1218 return Transform3D<T>(t.Rotation() * r, t.Translation());
1219}
1220template <class T>
1222{
1223 return Transform3D<T>(t.Rotation() * r, t.Translation());
1224}
1225template <class T>
1227{
1228 return Transform3D<T>(t.Rotation() * r, t.Translation());
1229}
1230template <class T>
1232{
1233 return Transform3D<T>(t.Rotation() * r, t.Translation());
1234}
1235
1236
1237
1238/**
1239 combine a rotation and a transformation to give a transform3d
1240 First the transformation then the rotation
1241 */
1242template <class T>
1244{
1245 return Transform3D<T>(r * t.Rotation(), r * t.Translation().Vect());
1246}
1247template <class T>
1249{
1250 Rotation3D r3d(r);
1251 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1252}
1253template <class T>
1255{
1256 Rotation3D r3d(r);
1257 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1258}
1259template <class T>
1261{
1262 Rotation3D r3d(r);
1263 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1264}
1265template <class T>
1267{
1268 Rotation3D r3d(r);
1269 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1270}
1271template <class T>
1273{
1274 Rotation3D r3d(r);
1275 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1276}
1277template <class T>
1279{
1280 Rotation3D r3d(r);
1281 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1282}
1283template <class T>
1285{
1286 Rotation3D r3d(r);
1287 return Transform3D<T>(r3d * t.Rotation(), r3d * t.Translation().Vect());
1288}
1289
1290
1291//---I/O functions
1292// TODO - I/O should be put in the manipulator form
1293
1294/**
1295 print the 12 components of the Transform3D
1296 */
1297template <class T>
1298std::ostream &operator<<(std::ostream &os, const Transform3D<T> &t)
1299{
1300 // TODO - this will need changing for machine-readable issues
1301 // and even the human readable form needs formatting improvements
1302
1303 T m[12];
1304 t.GetComponents(m, m + 12);
1305 os << "\n" << m[0] << " " << m[1] << " " << m[2] << " " << m[3];
1306 os << "\n" << m[4] << " " << m[5] << " " << m[6] << " " << m[7];
1307 os << "\n" << m[8] << " " << m[9] << " " << m[10] << " " << m[11] << "\n";
1308 return os;
1309}
1310
1311} // end namespace Impl
1312
1313// typedefs for double and float versions
1316
1317} // end namespace Math
1318
1319} // end namespace ROOT
1320
1321
1322#endif /* ROOT_Math_GenVector_Transform3D */
#define d(i)
Definition RSha256.hxx:102
TBuffer & operator<<(TBuffer &buf, const Tmpl *obj)
Definition TBuffer.h:399
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint xy
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
Option_t Option_t TPoint TPoint const char y1
float * q
AxisAngle class describing rotation represented with direction axis (3D Vector) and an angle of rotat...
Definition AxisAngle.h:42
DefaultCoordinateSystemTag Default tag for identifying any coordinate system.
Class describing a generic displacement vector in 3 dimensions.
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
DisplacementVector3D Cross(const DisplacementVector3D< OtherCoords, Tag > &v) const
Return vector (cross) product of two displacement vectors, as a vector in the coordinate system of th...
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
EulerAngles class describing rotation as three angles (Euler Angles).
Definition EulerAngles.h:45
Class describing a geometrical plane in 3 dimensions.
Definition Plane3D.h:53
Vector Normal() const
Return normal vector to the plane as Cartesian DisplacementVector.
Definition Plane3D.h:158
Scalar HesseDistance() const
Return the Hesse Distance (distance from the origin) of the plane or the d coefficient expressed in n...
Definition Plane3D.h:164
Basic 3D Transformation class describing a rotation and then a translation The internal data are a 3D...
Definition Transform3D.h:80
void GetDecomposition(AnyRotation &r, V &v) const
Get the rotation and translation vector representing the 3D transformation in any rotation and any ve...
void GetComponents(T &xx, T &xy, T &xz, T &dx, T &yx, T &yy, T &yz, T &dy, T &zx, T &zy, T &zz, T &dz) const
Get the components into 12 scalars.
Transform3D< T > Inverse() const
Return the inverse of the transformation.
void GetDecomposition(Rotation3D &r, Vector &v) const
Get the rotation and translation vector representing the 3D transformation.
Transform3D(const Point &fr0, const Point &fr1, const Point &fr2, const Point &to0, const Point &to1, const Point &to2)
Construct transformation from one coordinate system defined by three points (origin + two axis) to a ...
Transform3D(const AxisAngle &r)
Transform3D(const Rotation3D &r)
Construct from a 3D Rotation only with zero translation.
PositionVector3D< CoordSystem > operator*(const PositionVector3D< CoordSystem > &v) const
Transformation operation for Position Vector in any coordinate system.
void GetRotation(AnyRotation &r) const
Get the rotation (any type) representing the 3D transformation.
Translation3D< T > Translation() const
Get the translation representing the 3D transformation in a Cartesian vector.
void SetIdentity()
Set identity transformation (identity rotation , zero translation)
Transform3D(const Translation3D< T > &t)
Construct from a translation only, represented by a Translation3D class and with an identity rotation...
Plane3D< TYPE > operator*(const Plane3D< TYPE > &plane) const
Multiplication operator for 3D plane.
AnyRotation Rotation() const
Get the rotation representing the 3D transformation.
DisplacementVector3D< CoordSystem > ApplyInverse(const DisplacementVector3D< CoordSystem > &p) const
Directly apply the inverse affine transformation on an arbitrary coordinate-system vector.
Transform3D(T xx, T xy, T xz, T dx, T yx, T yy, T yz, T dy, T zx, T zy, T zz, T dz)
Raw constructor from 12 Scalar components.
Transform3D(const RotationZYX &r)
void AssignFrom(const Rotation3D &r)
make transformation from only rotations (zero translation)
PositionVector3D< CoordSystem > operator()(const PositionVector3D< CoordSystem > &p) const
Transformation operation for Position Vector in any coordinate system.
Vector ApplyInverse(const Vector &v) const
Directly apply the inverse affine transformation on vectors.
Transform3D(const RotationX &r)
void Transform(const DisplacementVector3D< CoordSystem, Tag1 > &v1, DisplacementVector3D< CoordSystem, Tag2 > &v2) const
Transformation operation for Displacement Vector of different coordinate systems.
Transform3D(IT begin, IT end)
Construct given a pair of pointers or iterators defining the beginning and end of an array of 12 Scal...
void Invert()
Invert the transformation in place (scalar)
void SetIdentity(const typename SCALAR::mask_type m)
Set identity transformation (identity rotation , zero translation) vectorised version that sets using...
Transform3D< T > & operator=(const ForeignMatrix &m)
Construct from a linear algebra matrix of size at least 3x4, which must support operator()(i,...
PositionVector3D< CoordSystem > ApplyInverse(const PositionVector3D< CoordSystem > &p) const
Directly apply the inverse affine transformation on an arbitrary coordinate-system point.
void GetTranslation(AnyVector &v) const
Get the translation representing the 3D transformation in any vector which implements the SetXYZ meth...
Transform3D(const DisplacementVector3D< CoordSystem, Tag > &v)
Construct from a translation only, represented by any DisplacementVector3D and with an identity rotat...
bool operator!=(const Transform3D< T > &rhs) const
Inequality operator.
void SetComponents(IT begin, IT end)
Set the 12 matrix components given an iterator to the start of the desired data, and another to the e...
void GetComponents(IT begin, IT end) const
Get the 12 matrix components into data specified by an iterator begin and another to the end of the d...
Transform3D()
Default constructor (identy rotation) + zero translation.
Definition Transform3D.h:99
Transform3D(const Rotation3D &r, const Translation3D< T > &t)
Construct from a rotation and then a translation described by a Translation3D class.
Rotation3D Rotation() const
Get the 3D rotation representing the 3D transformation.
Transform3D(const RotationY &r)
void GetComponents(IT begin) const
Get the 12 matrix components into data specified by an iterator begin.
PositionVector3D< Cartesian3D< T >, DefaultCoordinateSystemTag > Point
Definition Transform3D.h:86
Transform3D(const ARotation &r, const Translation3D< T > &t)
Construct from a rotation (any rotation object) and then a translation represented by a Translation3D...
void SetComponents(T xx, T xy, T xz, T dx, T yx, T yy, T yz, T dy, T zx, T zy, T zz, T dz)
Set the components from 12 scalars.
DisplacementVector3D< CoordSystem > operator*(const DisplacementVector3D< CoordSystem > &v) const
Transformation operation for Displacement Vector in any coordinate system.
Transform3D(const Quaternion &r)
bool operator==(const Transform3D< T > &rhs) const
Equality operator.
Transform3D(const RotationZ &r)
Transform3D< T > & operator*=(const Transform3D< T > &t)
multiply (combine) with another transformation in place
void AssignFrom(const Rotation3D &r, const Vector &v)
make transformation from first a rotation then a translation
Transform3D(const Rotation3D &r, const Vector &v)
Construct from a rotation and then a translation described by a Vector.
Transform3D(const EulerAngles &r)
Transform3D(const Vector &v)
Construct from a translation only, represented by a Cartesian 3D Vector, and with an identity rotatio...
DisplacementVector3D< Cartesian3D< T >, DefaultCoordinateSystemTag > Vector
Definition Transform3D.h:85
LorentzVector< CoordSystem > operator*(const LorentzVector< CoordSystem > &q) const
Transformation operation for a Lorentz Vector in any coordinate system.
Plane3D< TYPE > operator()(const Plane3D< TYPE > &plane) const
Transformation on a 3D plane.
Transform3D(const ForeignMatrix &m)
Construct from a linear algebra matrix of size at least 3x4, which must support operator()(i,...
void SetTransformMatrix(const ForeignMatrix &m)
Set components from a linear algebra matrix of size at least 3x4, which must support operator()(i,...
void GetTransformMatrix(ForeignMatrix &m) const
Get components into a linear algebra matrix of size at least 3x4, which must support operator()(i,...
Point operator()(const Point &p) const
Transformation operation for Position Vector in Cartesian coordinate For a Position Vector first a ro...
void AssignFrom(const Vector &v)
make transformation from only translation (identity rotations)
void Transform(const PositionVector3D< CoordSystem, Tag1 > &p1, PositionVector3D< CoordSystem, Tag2 > &p2) const
Transformation operation for points between different coordinate system tags.
Point ApplyInverse(const Point &p) const
Directly apply the inverse affine transformation on points (first inverse translation then inverse ro...
Transform3D(const ARotation &r, const DisplacementVector3D< CoordSystem, Tag > &v)
Construct from a rotation (any rotation object) and then a translation (represented by any Displaceme...
Class describing a 3 dimensional translation.
const Vector & Vect() const
return a const reference to the underline vector representing the translation
Class describing a generic LorentzVector in the 4D space-time, using the specified coordinate system ...
Class describing a generic position vector (point) in 3 dimensions.
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
PositionVector3D< CoordSystem, Tag > & SetXYZ(Scalar a, Scalar b, Scalar c)
set the values of the vector from the cartesian components (x,y,z) (if the vector is held in polar or...
Rotation class with the (3D) rotation represented by a unit quaternion (u, i, j, k).
Definition Quaternion.h:49
Rotation class with the (3D) rotation represented by a 3x3 orthogonal matrix.
Definition Rotation3D.h:67
Rotation class representing a 3D rotation about the X axis by the angle of rotation.
Definition RotationX.h:45
Rotation class representing a 3D rotation about the Y axis by the angle of rotation.
Definition RotationY.h:45
Rotation class with the (3D) rotation represented by angles describing first a rotation of an angle p...
Definition RotationZYX.h:63
Rotation class representing a 3D rotation about the Z axis by the angle of rotation.
Definition RotationZ.h:45
SVector< T, D > Unit(const SVector< T, D > &rhs)
Unit.
Definition Functions.h:382
const Int_t n
Definition legend1.C:16
Namespace for new Math classes and functions.
Transform3D< T > operator*(const Rotation3D &r, const Translation3D< T > &t)
combine a translation and a rotation to give a transform3d First the translation then the rotation
Impl::Transform3D< float > Transform3DF
Impl::Transform3D< double > Transform3D
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
TMarker m
Definition textangle.C:8