Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoMatrix.cxx
Go to the documentation of this file.
1// @(#)root/geom:$Id$
2// Author: Andrei Gheata 25/10/01
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TGeoMatrix
13\ingroup Geometry_classes
14
15Geometrical transformation package.
16
17 All geometrical transformations handled by the modeller are provided as a
18built-in package. This was designed to minimize memory requirements and
19optimize performance of point/vector master-to-local and local-to-master
20computation. We need to have in mind that a transformation in TGeo has 2
21major use-cases. The first one is for defining the placement of a volume
22with respect to its container reference frame. This frame will be called
23'master' and the frame of the positioned volume - 'local'. If T is a
24transformation used for positioning volume daughters, then:
25
26~~~ {.cpp}
27 MASTER = T * LOCAL
28~~~
29
30 Therefore a local-to-master conversion will be performed by using T, while
31a master-to-local by using its inverse. The second use case is the computation
32of the global transformation of a given object in the geometry. Since the
33geometry is built as 'volumes-inside-volumes', this global transformation
34represent the pile-up of all local transformations in the corresponding
35branch. The conversion from the global reference frame and the given object
36is also called master-to-local, but it is handled by the manager class.
37 A general homogenous transformation is defined as a 4x4 matrix embedding
38a rotation, a translation and a scale. The advantage of this description
39is that each basic transformation can be represented as a homogenous matrix,
40composition being performed as simple matrix multiplication.
41
42 Rotation: Inverse rotation:
43
44~~~ {.cpp}
45 r11 r12 r13 0 r11 r21 r31 0
46 r21 r22 r23 0 r12 r22 r32 0
47 r31 r32 r33 0 r13 r23 r33 0
48 0 0 0 1 0 0 0 1
49~~~
50
51 Translation: Inverse translation:
52
53~~~ {.cpp}
54 1 0 0 tx 1 0 0 -tx
55 0 1 0 ty 0 1 0 -ty
56 0 0 1 tz 0 0 1 -tz
57 0 0 0 1 0 0 0 1
58~~~
59
60 Scale: Inverse scale:
61
62~~~ {.cpp}
63 sx 0 0 0 1/sx 0 0 0
64 0 sy 0 0 0 1/sy 0 0
65 0 0 sz 0 0 0 1/sz 0
66 0 0 0 1 0 0 0 1
67~~~
68
69 where:
70 - `rij` are the 3x3 rotation matrix components,
71 - `tx`, `ty`, `tz` are the translation components
72 - `sx`, `sy`, `sz` are arbitrary scale constants on each axis,
73
74 The disadvantage in using this approach is that computation for 4x4 matrices
75is expensive. Even combining two translation would become a multiplication
76of their corresponding matrices, which is quite an undesired effect. On the
77other hand, it is not a good idea to store a translation as a block of 16
78numbers. We have therefore chosen to implement each basic transformation type
79as a class deriving from the same basic abstract class and handling its specific
80data and point/vector transformation algorithms.
81
82\image html geom_transf.jpg
83
84### The base class TGeoMatrix defines abstract metods for:
85
86#### translation, rotation and scale getters. Every derived class stores only
87 its specific data, e.g. a translation stores an array of 3 doubles and a
88 rotation an array of 9. However, asking which is the rotation array of a
89 TGeoTranslation through the base TGeoMatrix interface is a legal operation.
90 The answer in this case is a pointer to a global constant array representing
91 an identity rotation.
92
93~~~ {.cpp}
94 Double_t *TGeoMatrix::GetTranslation()
95 Double_t *TGeoMatrix::GetRotation()
96 Double_t *TGeoMatrix::GetScale()
97~~~
98
99#### MasterToLocal() and LocalToMaster() point and vector transformations :
100
101~~~ {.cpp}
102 void TGeoMatrix::MasterToLocal(const Double_t *master, Double_t *local)
103 void TGeoMatrix::LocalToMaster(const Double_t *local, Double_t *master)
104 void TGeoMatrix::MasterToLocalVect(const Double_t *master, Double_t *local)
105 void TGeoMatrix::LocalToMasterVect(const Double_t *local, Double_t *master)
106~~~
107
108 These allow correct conversion also for reflections.
109
110#### Transformation type getters :
111
112~~~ {.cpp}
113 Bool_t TGeoMatrix::IsIdentity()
114 Bool_t TGeoMatrix::IsTranslation()
115 Bool_t TGeoMatrix::IsRotation()
116 Bool_t TGeoMatrix::IsScale()
117 Bool_t TGeoMatrix::IsCombi() (translation + rotation)
118 Bool_t TGeoMatrix::IsGeneral() (translation + rotation + scale)
119~~~
120
121 Combinations of basic transformations are represented by specific classes
122deriving from TGeoMatrix. In order to define a matrix as a combination of several
123others, a special class TGeoHMatrix is provided. Here is an example of matrix
124creation :
125
126### Matrix creation example:
127
128~~~ {.cpp}
129 root[0] TGeoRotation r1,r2;
130 r1.SetAngles(90,0,30); // rotation defined by Euler angles
131 r2.SetAngles(90,90,90,180,0,0); // rotation defined by GEANT3 angles
132 TGeoTranslation t1(-10,10,0);
133 TGeoTranslation t2(10,-10,5);
134 TGeoCombiTrans c1(t1,r1);
135 TGeoCombiTrans c2(t2,r2);
136 TGeoHMatrix h = c1 * c2; // composition is done via TGeoHMatrix class
137 root[7] TGeoHMatrix *ph = new TGeoHMatrix(hm); // this is the one we want to
138 // use for positioning a volume
139 root[8] ph->Print();
140 ...
141 pVolume->AddNode(pVolDaughter,id,ph) // now ph is owned by the manager
142~~~
143
144### Rule for matrix creation:
145 Unless explicitly used for positioning nodes (TGeoVolume::AddNode()) all
146matrices deletion have to be managed by users. Matrices passed to geometry
147have to be created by using new() operator and their deletion is done by
148TGeoManager class.
149
150### Available geometrical transformations
151
152#### TGeoTranslation
153Represent a (dx,dy,dz) translation. Data members:
154 Double_t fTranslation[3]. Translations can be added/subtracted.
155
156~~~ {.cpp}
157 TGeoTranslation t1;
158 t1->SetTranslation(-5,10,4);
159 TGeoTranslation *t2 = new TGeoTranslation(4,3,10);
160 t2->Subtract(&t1);
161~~~
162
163#### Rotations
164 Represent a pure rotation. Data members: Double_t fRotationMatrix[3*3].
165 Rotations can be defined either by Euler angles, either, by GEANT3 angles :
166
167~~~ {.cpp}
168 TGeoRotation *r1 = new TGeoRotation();
169 r1->SetAngles(phi, theta, psi); // all angles in degrees
170~~~
171
172 This represent the composition of : first a rotation about Z axis with
173 angle phi, then a rotation with theta about the rotated X axis, and
174 finally a rotation with psi about the new Z axis.
175
176~~~ {.cpp}
177 r1->SetAngles(th1,phi1, th2,phi2, th3,phi3)
178~~~
179
180 This is a rotation defined in GEANT3 style. Theta and phi are the spherical
181 angles of each axis of the rotated coordinate system with respect to the
182 initial one. This construction allows definition of malformed rotations,
183 e.g. not orthogonal. A check is performed and an error message is issued
184 in this case.
185
186 Specific utilities : determinant, inverse.
187
188#### Scale transformations
189 Represent a scale shrinking/enlargement. Data
190 members :Double_t fScale[3]. Not fully implemented yet.
191
192#### Combined transformations
193Represent a rotation followed by a translation.
194Data members: Double_t fTranslation[3], TGeoRotation *fRotation.
195
196~~~ {.cpp}
197 TGeoRotation *rot = new TGeoRotation("rot",10,20,30);
198 TGeoTranslation trans;
199 ...
200 TGeoCombiTrans *c1 = new TGeoCombiTrans(trans, rot);
201 TGeoCombiTrans *c2 = new TGeoCombiTrans("somename",10,20,30,rot)
202~~~
203
204
205#### TGeoGenTrans
206Combined transformations including a scale. Not implemented.
207
208#### TGeoIdentity
209A generic singleton matrix representing a identity transformation
210 NOTE: identified by the global variable gGeoIdentity.
211*/
212
213#include <iostream>
214#include "TObjArray.h"
215
216#include "TGeoManager.h"
217#include "TGeoMatrix.h"
218#include "TMath.h"
219
221const Int_t kN3 = 3*sizeof(Double_t);
222const Int_t kN9 = 9*sizeof(Double_t);
223
224// statics and globals
225
227
228////////////////////////////////////////////////////////////////////////////////
229/// dummy constructor
230
232{
234}
235
236////////////////////////////////////////////////////////////////////////////////
237/// copy constructor
238
240 :TNamed(other)
241{
243}
244
245////////////////////////////////////////////////////////////////////////////////
246/// Constructor
247
249 :TNamed(name, "")
250{
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// Destructor
255
257{
258 if (IsRegistered() && gGeoManager) {
259 if (!gGeoManager->IsCleaning()) {
261 Warning("dtor", "Registered matrix %s was removed", GetName());
262 }
263 }
264}
265
266////////////////////////////////////////////////////////////////////////////////
267/// Returns true if no rotation or the rotation is about Z axis
268
270{
271 if (IsIdentity()) return kTRUE;
272 const Double_t *rot = GetRotationMatrix();
273 if (TMath::Abs(rot[6])>1E-9) return kFALSE;
274 if (TMath::Abs(rot[7])>1E-9) return kFALSE;
275 if ((1.-TMath::Abs(rot[8]))>1E-9) return kFALSE;
276 return kTRUE;
277}
278
279////////////////////////////////////////////////////////////////////////////////
280/// Get total size in bytes of this
281
283{
284 Int_t count = 4+28+strlen(GetName())+strlen(GetTitle()); // fId + TNamed
285 if (IsTranslation()) count += 12;
286 if (IsScale()) count += 12;
287 if (IsCombi() || IsGeneral()) count += 4 + 36;
288 return count;
289}
290
291////////////////////////////////////////////////////////////////////////////////
292/// Provide a pointer name containing uid.
293
294const char *TGeoMatrix::GetPointerName() const
295{
296 static TString name;
297 name.Form("pMatrix%d", GetUniqueID());
298 return name.Data();
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// The homogenous matrix associated with the transformation is used for
303/// piling up's and visualization. A homogenous matrix is a 4*4 array
304/// containing the translation, the rotation and the scale components
305/// ~~~ {.cpp}
306/// | R00*sx R01 R02 dx |
307/// | R10 R11*sy R12 dy |
308/// | R20 R21 R22*sz dz |
309/// | 0 0 0 1 |
310/// ~~~
311/// where Rij is the rotation matrix, (sx, sy, sz) is the scale
312/// transformation and (dx, dy, dz) is the translation.
313
315{
316 Double_t *hmatrix = hmat;
317 const Double_t *mat = GetRotationMatrix();
318 for (Int_t i=0; i<3; i++) {
319 memcpy(hmatrix, mat, kN3);
320 mat += 3;
321 hmatrix += 3;
322 *hmatrix = 0.0;
323 hmatrix++;
324 }
325 memcpy(hmatrix, GetTranslation(), kN3);
326 hmatrix = hmat;
327 if (IsScale()) {
328 for (Int_t i=0; i<3; i++) {
329 *hmatrix *= GetScale()[i];
330 hmatrix += 5;
331 }
332 }
333 hmatrix[15] = 1.;
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
338
339void TGeoMatrix::LocalToMaster(const Double_t *local, Double_t *master) const
340{
341 if (IsIdentity()) {
342 memcpy(master, local, kN3);
343 return;
344 }
345 Int_t i;
346 const Double_t *tr = GetTranslation();
347 if (!IsRotation()) {
348 for (i=0; i<3; i++) master[i] = tr[i] + local[i];
349 return;
350 }
351 const Double_t *rot = GetRotationMatrix();
352 for (i=0; i<3; i++) {
353 master[i] = tr[i]
354 + local[0]*rot[3*i]
355 + local[1]*rot[3*i+1]
356 + local[2]*rot[3*i+2];
357 }
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// convert a vector by multiplying its column vector (x, y, z, 1) to matrix inverse
362
363void TGeoMatrix::LocalToMasterVect(const Double_t *local, Double_t *master) const
364{
365 if (!IsRotation()) {
366 memcpy(master, local, kN3);
367 return;
368 }
369 const Double_t *rot = GetRotationMatrix();
370 for (Int_t i=0; i<3; i++) {
371 master[i] = local[0]*rot[3*i]
372 + local[1]*rot[3*i+1]
373 + local[2]*rot[3*i+2];
374 }
375}
376
377////////////////////////////////////////////////////////////////////////////////
378/// convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
379
380void TGeoMatrix::LocalToMasterBomb(const Double_t *local, Double_t *master) const
381{
382 if (IsIdentity()) {
383 memcpy(master, local, kN3);
384 return;
385 }
386 Int_t i;
387 const Double_t *tr = GetTranslation();
388 Double_t bombtr[3] = {0.,0.,0.};
389 gGeoManager->BombTranslation(tr, &bombtr[0]);
390 if (!IsRotation()) {
391 for (i=0; i<3; i++) master[i] = bombtr[i] + local[i];
392 return;
393 }
394 const Double_t *rot = GetRotationMatrix();
395 for (i=0; i<3; i++) {
396 master[i] = bombtr[i]
397 + local[0]*rot[3*i]
398 + local[1]*rot[3*i+1]
399 + local[2]*rot[3*i+2];
400 }
401}
402
403////////////////////////////////////////////////////////////////////////////////
404/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
405
406void TGeoMatrix::MasterToLocal(const Double_t *master, Double_t *local) const
407{
408 if (IsIdentity()) {
409 memcpy(local, master, kN3);
410 return;
411 }
412 const Double_t *tr = GetTranslation();
413 Double_t mt0 = master[0]-tr[0];
414 Double_t mt1 = master[1]-tr[1];
415 Double_t mt2 = master[2]-tr[2];
416 if (!IsRotation()) {
417 local[0] = mt0;
418 local[1] = mt1;
419 local[2] = mt2;
420 return;
421 }
422 const Double_t *rot = GetRotationMatrix();
423 local[0] = mt0*rot[0] + mt1*rot[3] + mt2*rot[6];
424 local[1] = mt0*rot[1] + mt1*rot[4] + mt2*rot[7];
425 local[2] = mt0*rot[2] + mt1*rot[5] + mt2*rot[8];
426}
427
428////////////////////////////////////////////////////////////////////////////////
429/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
430
431void TGeoMatrix::MasterToLocalVect(const Double_t *master, Double_t *local) const
432{
433 if (!IsRotation()) {
434 memcpy(local, master, kN3);
435 return;
436 }
437 const Double_t *rot = GetRotationMatrix();
438 for (Int_t i=0; i<3; i++) {
439 local[i] = master[0]*rot[i]
440 + master[1]*rot[i+3]
441 + master[2]*rot[i+6];
442 }
443}
444
445////////////////////////////////////////////////////////////////////////////////
446/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
447
448void TGeoMatrix::MasterToLocalBomb(const Double_t *master, Double_t *local) const
449{
450 if (IsIdentity()) {
451 memcpy(local, master, kN3);
452 return;
453 }
454 const Double_t *tr = GetTranslation();
455 Double_t bombtr[3] = {0.,0.,0.};
456 Int_t i;
457 gGeoManager->UnbombTranslation(tr, &bombtr[0]);
458 if (!IsRotation()) {
459 for (i=0; i<3; i++) local[i] = master[i]-bombtr[i];
460 return;
461 }
462 const Double_t *rot = GetRotationMatrix();
463 for (i=0; i<3; i++) {
464 local[i] = (master[0]-bombtr[0])*rot[i]
465 + (master[1]-bombtr[1])*rot[i+3]
466 + (master[2]-bombtr[2])*rot[i+6];
467 }
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// Normalize a vector.
472
474{
475 Double_t normfactor = vect[0]*vect[0] + vect[1]*vect[1] + vect[2]*vect[2];
476 if (normfactor <= 1E-10) return;
477 normfactor = 1./TMath::Sqrt(normfactor);
478 vect[0] *= normfactor;
479 vect[1] *= normfactor;
480 vect[2] *= normfactor;
481}
482
483////////////////////////////////////////////////////////////////////////////////
484/// print the matrix in 4x4 format
485
487{
488 const Double_t *rot = GetRotationMatrix();
489 const Double_t *tr = GetTranslation();
490 printf("matrix %s - tr=%d rot=%d refl=%d scl=%d shr=%d reg=%d own=%d\n", GetName(),(Int_t)IsTranslation(),
492 (Int_t)IsOwned());
493 printf("%10.6f%12.6f%12.6f Tx = %10.6f\n", rot[0], rot[1], rot[2], tr[0]);
494 printf("%10.6f%12.6f%12.6f Ty = %10.6f\n", rot[3], rot[4], rot[5], tr[1]);
495 printf("%10.6f%12.6f%12.6f Tz = %10.6f\n", rot[6], rot[7], rot[8], tr[2]);
496 if (IsScale()) {
497 const Double_t *scl = GetScale();
498 printf("Sx=%10.6fSy=%12.6fSz=%12.6f\n", scl[0], scl[1], scl[2]);
499 }
500}
501
502////////////////////////////////////////////////////////////////////////////////
503/// Multiply by a reflection respect to YZ.
504
506{
507}
508
509////////////////////////////////////////////////////////////////////////////////
510/// Multiply by a reflection respect to ZX.
511
513{
514}
515
516////////////////////////////////////////////////////////////////////////////////
517/// Multiply by a reflection respect to XY.
518
520{
521}
522
523////////////////////////////////////////////////////////////////////////////////
524/// Register the matrix in the current manager, which will become the owner.
525
527{
528 if (!gGeoManager) {
529 Warning("RegisterYourself", "cannot register without geometry");
530 return;
531 }
532 if (!IsRegistered()) {
535 }
536}
537
538////////////////////////////////////////////////////////////////////////////////
539/// If no name was supplied in the ctor, the type of transformation is checked.
540/// A letter will be prepended to the name :
541/// - t - translation
542/// - r - rotation
543/// - s - scale
544/// - c - combi (translation + rotation)
545/// - g - general (tr+rot+scale)
546/// The index of the transformation in gGeoManager list of transformations will
547/// be appended.
548
550{
551 if (!gGeoManager) return;
552 if (strlen(GetName())) return;
553 char type = 'n';
554 if (IsTranslation()) type = 't';
555 if (IsRotation()) type = 'r';
556 if (IsScale()) type = 's';
557 if (IsCombi()) type = 'c';
558 if (IsGeneral()) type = 'g';
560 Int_t index = 0;
561 if (matrices) index =matrices->GetEntriesFast() - 1;
563 SetName(name);
564}
565
566/** \class TGeoTranslation
567\ingroup Geometry_classes
568
569Class describing translations. A translation is
570basically an array of 3 doubles matching the positions 12, 13
571and 14 in the homogenous matrix description.
572*/
573
575
576////////////////////////////////////////////////////////////////////////////////
577/// Default constructor
578
580{
581 for (Int_t i=0; i<3; i++) fTranslation[i] = 0;
582}
583
584////////////////////////////////////////////////////////////////////////////////
585/// Copy ctor.
586
588 :TGeoMatrix(other)
589{
590 SetTranslation(other);
591}
592
593////////////////////////////////////////////////////////////////////////////////
594/// Ctor. based on a general matrix
595
597 :TGeoMatrix(other)
598{
601 SetTranslation(other);
602}
603
604////////////////////////////////////////////////////////////////////////////////
605/// Default constructor defining the translation
606
608 :TGeoMatrix("")
609{
610 if (dx || dy || dz) SetBit(kGeoTranslation);
611 SetTranslation(dx, dy, dz);
612}
613
614////////////////////////////////////////////////////////////////////////////////
615/// Default constructor defining the translation
616
619{
620 if (dx || dy || dz) SetBit(kGeoTranslation);
621 SetTranslation(dx, dy, dz);
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// Assignment from a general matrix
626
628{
629 if (&matrix == this) return *this;
630 Bool_t registered = TestBit(kGeoRegistered);
631 TNamed::operator=(matrix);
632 SetTranslation(matrix);
633 SetBit(kGeoRegistered,registered);
636 return *this;
637}
638
639////////////////////////////////////////////////////////////////////////////////
640/// Translation composition
641
643{
644 const Double_t *tr = right.GetTranslation();
645 fTranslation[0] += tr[0];
646 fTranslation[1] += tr[1];
647 fTranslation[2] += tr[2];
649 return *this;
650}
651
653{
654 TGeoTranslation t = *this;
655 t *= right;
656 return t;
657}
658
660{
661 TGeoHMatrix t = *this;
662 t *= right;
663 return t;
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Is-equal operator
668
670{
671 if (&other == this) return kTRUE;
672 const Double_t *tr = GetTranslation();
673 const Double_t *otr = other.GetTranslation();
674 for (auto i=0; i<3; i++)
675 if (TMath::Abs(tr[i]-otr[i])>1.E-10) return kFALSE;
676 return kTRUE;
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Return a temporary inverse of this.
681
683{
685 h = *this;
687 Double_t tr[3];
688 tr[0] = -fTranslation[0];
689 tr[1] = -fTranslation[1];
690 tr[2] = -fTranslation[2];
691 h.SetTranslation(tr);
692 return h;
693}
694
695////////////////////////////////////////////////////////////////////////////////
696/// Adding a translation to this one
697
699{
700 const Double_t *trans = other->GetTranslation();
701 for (Int_t i=0; i<3; i++)
702 fTranslation[i] += trans[i];
703}
704
705////////////////////////////////////////////////////////////////////////////////
706/// Make a clone of this matrix.
707
709{
710 TGeoMatrix *matrix = new TGeoTranslation(*this);
711 return matrix;
712}
713
714////////////////////////////////////////////////////////////////////////////////
715/// Rotate about X axis of the master frame with angle expressed in degrees.
716
718{
719 Warning("RotateX", "Not implemented. Use TGeoCombiTrans instead");
720}
721
722////////////////////////////////////////////////////////////////////////////////
723/// Rotate about Y axis of the master frame with angle expressed in degrees.
724
726{
727 Warning("RotateY", "Not implemented. Use TGeoCombiTrans instead");
728}
729
730////////////////////////////////////////////////////////////////////////////////
731/// Rotate about Z axis of the master frame with angle expressed in degrees.
732
734{
735 Warning("RotateZ", "Not implemented. Use TGeoCombiTrans instead");
736}
737
738////////////////////////////////////////////////////////////////////////////////
739/// Subtracting a translation from this one
740
742{
743 const Double_t *trans = other->GetTranslation();
744 for (Int_t i=0; i<3; i++)
745 fTranslation[i] -= trans[i];
746}
747
748////////////////////////////////////////////////////////////////////////////////
749/// Set translation components
750
752{
753 fTranslation[0] = dx;
754 fTranslation[1] = dy;
755 fTranslation[2] = dz;
756 if (dx || dy || dz) SetBit(kGeoTranslation);
758}
759
760////////////////////////////////////////////////////////////////////////////////
761/// Set translation components
762
764{
766 const Double_t *transl = other.GetTranslation();
767 memcpy(fTranslation, transl, kN3);
768}
769
770////////////////////////////////////////////////////////////////////////////////
771/// convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
772
773void TGeoTranslation::LocalToMaster(const Double_t *local, Double_t *master) const
774{
775 const Double_t *tr = GetTranslation();
776 for (Int_t i=0; i<3; i++)
777 master[i] = tr[i] + local[i];
778}
779
780////////////////////////////////////////////////////////////////////////////////
781/// convert a vector to MARS
782
783void TGeoTranslation::LocalToMasterVect(const Double_t *local, Double_t *master) const
784{
785 memcpy(master, local, kN3);
786}
787
788////////////////////////////////////////////////////////////////////////////////
789/// convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
790
791void TGeoTranslation::LocalToMasterBomb(const Double_t *local, Double_t *master) const
792{
793 const Double_t *tr = GetTranslation();
794 Double_t bombtr[3] = {0.,0.,0.};
795 gGeoManager->BombTranslation(tr, &bombtr[0]);
796 for (Int_t i=0; i<3; i++)
797 master[i] = bombtr[i] + local[i];
798}
799
800////////////////////////////////////////////////////////////////////////////////
801/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
802
803void TGeoTranslation::MasterToLocal(const Double_t *master, Double_t *local) const
804{
805 const Double_t *tr = GetTranslation();
806 for (Int_t i=0; i<3; i++)
807 local[i] = master[i]-tr[i];
808}
809
810////////////////////////////////////////////////////////////////////////////////
811/// convert a vector from MARS to local
812
813void TGeoTranslation::MasterToLocalVect(const Double_t *master, Double_t *local) const
814{
815 memcpy(local, master, kN3);
816}
817
818////////////////////////////////////////////////////////////////////////////////
819/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
820
821void TGeoTranslation::MasterToLocalBomb(const Double_t *master, Double_t *local) const
822{
823 const Double_t *tr = GetTranslation();
824 Double_t bombtr[3] = {0.,0.,0.};
825 gGeoManager->UnbombTranslation(tr, &bombtr[0]);
826 for (Int_t i=0; i<3; i++)
827 local[i] = master[i]-bombtr[i];
828}
829
830////////////////////////////////////////////////////////////////////////////////
831/// Save a primitive as a C++ statement(s) on output stream "out".
832
833void TGeoTranslation::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
834{
835 if (TestBit(kGeoSavePrimitive)) return;
836 out << " // Translation: " << GetName() << std::endl;
837 out << " dx = " << fTranslation[0] << ";" << std::endl;
838 out << " dy = " << fTranslation[1] << ";" << std::endl;
839 out << " dz = " << fTranslation[2] << ";" << std::endl;
840 out << " TGeoTranslation *" << GetPointerName() << " = new TGeoTranslation(\"" << GetName() << "\",dx,dy,dz);" << std::endl;
842}
843
844/** \class TGeoRotation
845\ingroup Geometry_classes
846Class describing rotations. A rotation is a 3*3 array
847Column vectors has to be orthogonal unit vectors.
848*/
849
851
852////////////////////////////////////////////////////////////////////////////////
853/// Default constructor.
854
856{
857 for (Int_t i=0; i<9; i++) {
858 if (i%4) fRotationMatrix[i] = 0;
859 else fRotationMatrix[i] = 1.0;
860 }
861}
862
863////////////////////////////////////////////////////////////////////////////////
864/// Copy ctor.
865
867 :TGeoMatrix(other)
868{
869 SetRotation(other);
870}
871
872////////////////////////////////////////////////////////////////////////////////
873/// Copy ctor.
874
876 :TGeoMatrix(other)
877{
880 SetRotation(other);
881}
882
883////////////////////////////////////////////////////////////////////////////////
884/// Named rotation constructor
885
888{
889 for (Int_t i=0; i<9; i++) {
890 if (i%4) fRotationMatrix[i] = 0;
891 else fRotationMatrix[i] = 1.0;
892 }
893}
894
895////////////////////////////////////////////////////////////////////////////////
896/// Default rotation constructor with Euler angles. Phi is the rotation angle about
897/// Z axis and is done first, theta is the rotation about new X and is done
898/// second, psi is the rotation angle about new Z and is done third. All angles are in
899/// degrees.
900
903{
904 SetAngles(phi, theta, psi);
905}
906
907////////////////////////////////////////////////////////////////////////////////
908/// Rotation constructor a la GEANT3. Angles theta(i), phi(i) are the polar and azimuthal
909/// angles of the (i) axis of the rotated system with respect to the initial non-rotated
910/// system.
911/// Example : the identity matrix (no rotation) is composed by
912/// theta1=90, phi1=0, theta2=90, phi2=90, theta3=0, phi3=0
913/// SetBit(kGeoRotation);
914
915TGeoRotation::TGeoRotation(const char *name, Double_t theta1, Double_t phi1, Double_t theta2, Double_t phi2,
916 Double_t theta3, Double_t phi3)
918{
919 SetAngles(theta1, phi1, theta2, phi2, theta3, phi3);
920}
921
922////////////////////////////////////////////////////////////////////////////////
923/// Assignment from a general matrix
924
926{
927 if (&other == this) return *this;
928 Bool_t registered = TestBit(kGeoRegistered);
929 TNamed::operator=(other);
930 SetRotation(other);
931 SetBit(kGeoRegistered,registered);
934 return *this;
935}
936
937////////////////////////////////////////////////////////////////////////////////
938/// Composition
939
941{
942 if (!right.IsRotation()) return *this;
943 MultiplyBy(&right, true);
944 return *this;
945}
946
948{
949 TGeoRotation r = *this;
950 r *= right;
951 return r;
952}
953
955{
956 TGeoHMatrix t = *this;
957 t *= right;
958 return t;
959}
960
961////////////////////////////////////////////////////////////////////////////////
962/// Is-equal operator
963
965{
966 if (&other == this) return kTRUE;
967 const Double_t *rot = GetRotationMatrix();
968 const Double_t *orot = other.GetRotationMatrix();
969 for (auto i=0; i<9; i++)
970 if (TMath::Abs(rot[i]-orot[i])>1.E-10) return kFALSE;
971 return kTRUE;
972}
973
974////////////////////////////////////////////////////////////////////////////////
975/// Return a temporary inverse of this.
976
978{
980 h = *this;
982 Double_t newrot[9];
983 newrot[0] = fRotationMatrix[0];
984 newrot[1] = fRotationMatrix[3];
985 newrot[2] = fRotationMatrix[6];
986 newrot[3] = fRotationMatrix[1];
987 newrot[4] = fRotationMatrix[4];
988 newrot[5] = fRotationMatrix[7];
989 newrot[6] = fRotationMatrix[2];
990 newrot[7] = fRotationMatrix[5];
991 newrot[8] = fRotationMatrix[8];
992 h.SetRotation(newrot);
993 return h;
994}
995
996////////////////////////////////////////////////////////////////////////////////
997/// Perform orthogonality test for rotation.
998
1000{
1001 const Double_t *r = fRotationMatrix;
1002 Double_t cij;
1003 for (Int_t i=0; i<2; i++) {
1004 for (Int_t j=i+1; j<3; j++) {
1005 // check columns
1006 cij = TMath::Abs(r[i]*r[j]+r[i+3]*r[j+3]+r[i+6]*r[j+6]);
1007 if (cij>1E-4) return kFALSE;
1008 // check rows
1009 cij = TMath::Abs(r[3*i]*r[3*j]+r[3*i+1]*r[3*j+1]+r[3*i+2]*r[3*j+2]);
1010 if (cij>1E-4) return kFALSE;
1011 }
1012 }
1013 return kTRUE;
1014}
1015
1016////////////////////////////////////////////////////////////////////////////////
1017/// reset data members
1018
1020{
1023}
1024
1025////////////////////////////////////////////////////////////////////////////////
1026/// Perform a rotation about Z having the sine/cosine of the rotation angle.
1027
1029{
1030 fRotationMatrix[0] = sincos[1];
1031 fRotationMatrix[1] = -sincos[0];
1032 fRotationMatrix[3] = sincos[0];
1033 fRotationMatrix[4] = sincos[1];
1035}
1036
1037////////////////////////////////////////////////////////////////////////////////
1038/// Returns rotation angle about Z axis in degrees. If the rotation is a pure
1039/// rotation about Z, fixX parameter does not matter, otherwise its meaning is:
1040/// - fixX = true : result is the phi angle of the projection of the rotated X axis in the un-rotated XY
1041/// - fixX = false : result is the phi angle of the projection of the rotated Y axis - 90 degrees
1042
1044{
1045 Double_t phi;
1046 if (fixX) phi = 180.*TMath::ATan2(-fRotationMatrix[1],fRotationMatrix[4])/TMath::Pi();
1047 else phi = 180.*TMath::ATan2(fRotationMatrix[3], fRotationMatrix[0])/TMath::Pi();
1048 return phi;
1049}
1050
1051////////////////////////////////////////////////////////////////////////////////
1052/// convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
1053
1054void TGeoRotation::LocalToMaster(const Double_t *local, Double_t *master) const
1055{
1056 const Double_t *rot = GetRotationMatrix();
1057 for (Int_t i=0; i<3; i++) {
1058 master[i] = local[0]*rot[3*i]
1059 + local[1]*rot[3*i+1]
1060 + local[2]*rot[3*i+2];
1061 }
1062}
1063
1064////////////////////////////////////////////////////////////////////////////////
1065/// convert a point by multiplying its column vector (x, y, z, 1) to matrix
1066
1067void TGeoRotation::MasterToLocal(const Double_t *master, Double_t *local) const
1068{
1069 const Double_t *rot = GetRotationMatrix();
1070 for (Int_t i=0; i<3; i++) {
1071 local[i] = master[0]*rot[i]
1072 + master[1]*rot[i+3]
1073 + master[2]*rot[i+6];
1074 }
1075}
1076
1077////////////////////////////////////////////////////////////////////////////////
1078/// Make a clone of this matrix.
1079
1081{
1082 TGeoMatrix *matrix = new TGeoRotation(*this);
1083 return matrix;
1084}
1085
1086////////////////////////////////////////////////////////////////////////////////
1087/// Rotate about X axis of the master frame with angle expressed in degrees.
1088
1090{
1093 Double_t c = TMath::Cos(phi);
1094 Double_t s = TMath::Sin(phi);
1095 Double_t v[9];
1096 v[0] = fRotationMatrix[0];
1097 v[1] = fRotationMatrix[1];
1098 v[2] = fRotationMatrix[2];
1099 v[3] = c*fRotationMatrix[3]-s*fRotationMatrix[6];
1100 v[4] = c*fRotationMatrix[4]-s*fRotationMatrix[7];
1101 v[5] = c*fRotationMatrix[5]-s*fRotationMatrix[8];
1102 v[6] = s*fRotationMatrix[3]+c*fRotationMatrix[6];
1103 v[7] = s*fRotationMatrix[4]+c*fRotationMatrix[7];
1104 v[8] = s*fRotationMatrix[5]+c*fRotationMatrix[8];
1105
1106 memcpy(fRotationMatrix, v, kN9);
1107}
1108
1109////////////////////////////////////////////////////////////////////////////////
1110/// Rotate about Y axis of the master frame with angle expressed in degrees.
1111
1113{
1116 Double_t c = TMath::Cos(phi);
1117 Double_t s = TMath::Sin(phi);
1118 Double_t v[9];
1119 v[0] = c*fRotationMatrix[0]+s*fRotationMatrix[6];
1120 v[1] = c*fRotationMatrix[1]+s*fRotationMatrix[7];
1121 v[2] = c*fRotationMatrix[2]+s*fRotationMatrix[8];
1122 v[3] = fRotationMatrix[3];
1123 v[4] = fRotationMatrix[4];
1124 v[5] = fRotationMatrix[5];
1125 v[6] = -s*fRotationMatrix[0]+c*fRotationMatrix[6];
1126 v[7] = -s*fRotationMatrix[1]+c*fRotationMatrix[7];
1127 v[8] = -s*fRotationMatrix[2]+c*fRotationMatrix[8];
1128
1129 memcpy(fRotationMatrix, v, kN9);
1130}
1131
1132////////////////////////////////////////////////////////////////////////////////
1133/// Rotate about Z axis of the master frame with angle expressed in degrees.
1134
1136{
1139 Double_t c = TMath::Cos(phi);
1140 Double_t s = TMath::Sin(phi);
1141 Double_t v[9];
1142 v[0] = c*fRotationMatrix[0]-s*fRotationMatrix[3];
1143 v[1] = c*fRotationMatrix[1]-s*fRotationMatrix[4];
1144 v[2] = c*fRotationMatrix[2]-s*fRotationMatrix[5];
1145 v[3] = s*fRotationMatrix[0]+c*fRotationMatrix[3];
1146 v[4] = s*fRotationMatrix[1]+c*fRotationMatrix[4];
1147 v[5] = s*fRotationMatrix[2]+c*fRotationMatrix[5];
1148 v[6] = fRotationMatrix[6];
1149 v[7] = fRotationMatrix[7];
1150 v[8] = fRotationMatrix[8];
1151
1152 memcpy(&fRotationMatrix[0],v,kN9);
1153}
1154
1155////////////////////////////////////////////////////////////////////////////////
1156/// Multiply by a reflection respect to YZ.
1157
1159{
1160 if (leftside) {
1164 } else {
1168 }
1171}
1172
1173////////////////////////////////////////////////////////////////////////////////
1174/// Multiply by a reflection respect to ZX.
1175
1177{
1178 if (leftside) {
1182 } else {
1186 }
1189}
1190
1191////////////////////////////////////////////////////////////////////////////////
1192/// Multiply by a reflection respect to XY.
1193
1195{
1196 if (leftside) {
1200 } else {
1204 }
1207}
1208
1209////////////////////////////////////////////////////////////////////////////////
1210/// Save a primitive as a C++ statement(s) on output stream "out".
1211
1212void TGeoRotation::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
1213{
1214 if (TestBit(kGeoSavePrimitive)) return;
1215 out << " // Rotation: " << GetName() << std::endl;
1216 Double_t th1,ph1,th2,ph2,th3,ph3;
1217 GetAngles(th1,ph1,th2,ph2,th3,ph3);
1218 out << " thx = " << th1 << "; phx = " << ph1 << ";" << std::endl;
1219 out << " thy = " << th2 << "; phy = " << ph2 << ";" << std::endl;
1220 out << " thz = " << th3 << "; phz = " << ph3 << ";" << std::endl;
1221 out << " TGeoRotation *" << GetPointerName() << " = new TGeoRotation(\"" << GetName() << "\",thx,phx,thy,phy,thz,phz);" << std::endl;
1223}
1224
1225////////////////////////////////////////////////////////////////////////////////
1226/// Copy rotation elements from other rotation matrix.
1227
1229{
1230 SetBit(kGeoRotation, other.IsRotation());
1232}
1233
1234////////////////////////////////////////////////////////////////////////////////
1235/// Set matrix elements according to Euler angles. Phi is the rotation angle about
1236/// Z axis and is done first, theta is the rotation about new X and is done
1237/// second, psi is the rotation angle about new Z and is done third. All angles are in
1238/// degrees.
1239
1241{
1242 Double_t degrad = TMath::Pi()/180.;
1243 Double_t sinphi = TMath::Sin(degrad*phi);
1244 Double_t cosphi = TMath::Cos(degrad*phi);
1245 Double_t sinthe = TMath::Sin(degrad*theta);
1246 Double_t costhe = TMath::Cos(degrad*theta);
1247 Double_t sinpsi = TMath::Sin(degrad*psi);
1248 Double_t cospsi = TMath::Cos(degrad*psi);
1249
1250 fRotationMatrix[0] = cospsi*cosphi - costhe*sinphi*sinpsi;
1251 fRotationMatrix[1] = -sinpsi*cosphi - costhe*sinphi*cospsi;
1252 fRotationMatrix[2] = sinthe*sinphi;
1253 fRotationMatrix[3] = cospsi*sinphi + costhe*cosphi*sinpsi;
1254 fRotationMatrix[4] = -sinpsi*sinphi + costhe*cosphi*cospsi;
1255 fRotationMatrix[5] = -sinthe*cosphi;
1256 fRotationMatrix[6] = sinpsi*sinthe;
1257 fRotationMatrix[7] = cospsi*sinthe;
1258 fRotationMatrix[8] = costhe;
1259
1260 if (!IsValid()) Error("SetAngles", "invalid rotation (Euler angles : phi=%f theta=%f psi=%f)",phi,theta,psi);
1261 CheckMatrix();
1262}
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// Set matrix elements in the GEANT3 way
1266
1268 Double_t theta3, Double_t phi3)
1269{
1270 Double_t degrad = TMath::Pi()/180.;
1271 fRotationMatrix[0] = TMath::Cos(degrad*phi1)*TMath::Sin(degrad*theta1);
1272 fRotationMatrix[3] = TMath::Sin(degrad*phi1)*TMath::Sin(degrad*theta1);
1273 fRotationMatrix[6] = TMath::Cos(degrad*theta1);
1274 fRotationMatrix[1] = TMath::Cos(degrad*phi2)*TMath::Sin(degrad*theta2);
1275 fRotationMatrix[4] = TMath::Sin(degrad*phi2)*TMath::Sin(degrad*theta2);
1276 fRotationMatrix[7] = TMath::Cos(degrad*theta2);
1277 fRotationMatrix[2] = TMath::Cos(degrad*phi3)*TMath::Sin(degrad*theta3);
1278 fRotationMatrix[5] = TMath::Sin(degrad*phi3)*TMath::Sin(degrad*theta3);
1279 fRotationMatrix[8] = TMath::Cos(degrad*theta3);
1280 // do the trick to eliminate most of the floating point errors
1281 for (Int_t i=0; i<9; i++) {
1282 if (TMath::Abs(fRotationMatrix[i])<1E-15) fRotationMatrix[i] = 0;
1283 if (TMath::Abs(fRotationMatrix[i]-1)<1E-15) fRotationMatrix[i] = 1;
1284 if (TMath::Abs(fRotationMatrix[i]+1)<1E-15) fRotationMatrix[i] = -1;
1285 }
1286 if (!IsValid()) Error("SetAngles", "invalid rotation (G3 angles, th1=%f phi1=%f, th2=%f ph2=%f, th3=%f phi3=%f)",
1287 theta1,phi1,theta2,phi2,theta3,phi3);
1288 CheckMatrix();
1289}
1290
1291////////////////////////////////////////////////////////////////////////////////
1292/// Retrieve rotation angles
1293
1294void TGeoRotation::GetAngles(Double_t &theta1, Double_t &phi1, Double_t &theta2, Double_t &phi2,
1295 Double_t &theta3, Double_t &phi3) const
1296{
1297 Double_t raddeg = 180./TMath::Pi();
1298 theta1 = raddeg*TMath::ACos(fRotationMatrix[6]);
1299 theta2 = raddeg*TMath::ACos(fRotationMatrix[7]);
1300 theta3 = raddeg*TMath::ACos(fRotationMatrix[8]);
1301 if (TMath::Abs(fRotationMatrix[0])<1E-6 && TMath::Abs(fRotationMatrix[3])<1E-6) phi1=0.;
1302 else phi1 = raddeg*TMath::ATan2(fRotationMatrix[3],fRotationMatrix[0]);
1303 if (phi1<0) phi1+=360.;
1304 if (TMath::Abs(fRotationMatrix[1])<1E-6 && TMath::Abs(fRotationMatrix[4])<1E-6) phi2=0.;
1305 else phi2 = raddeg*TMath::ATan2(fRotationMatrix[4],fRotationMatrix[1]);
1306 if (phi2<0) phi2+=360.;
1307 if (TMath::Abs(fRotationMatrix[2])<1E-6 && TMath::Abs(fRotationMatrix[5])<1E-6) phi3=0.;
1308 else phi3 = raddeg*TMath::ATan2(fRotationMatrix[5],fRotationMatrix[2]);
1309 if (phi3<0) phi3+=360.;
1310}
1311
1312////////////////////////////////////////////////////////////////////////////////
1313/// Retrieve Euler angles.
1314
1316{
1317 const Double_t *m = fRotationMatrix;
1318 // Check if theta is 0 or 180.
1319 if (TMath::Abs(1.-TMath::Abs(m[8]))<1.e-9) {
1320 theta = TMath::ACos(m[8])*TMath::RadToDeg();
1321 phi = TMath::ATan2(-m[8]*m[1],m[0])*TMath::RadToDeg();
1322 psi = 0.; // convention, phi+psi matters
1323 return;
1324 }
1325 // sin(theta) != 0
1326 phi = TMath::ATan2(m[2],-m[5]);
1327 Double_t sphi = TMath::Sin(phi);
1328 if (TMath::Abs(sphi)<1.e-9) theta = -TMath::ASin(m[5]/TMath::Cos(phi))*TMath::RadToDeg();
1329 else theta = TMath::ASin(m[2]/sphi)*TMath::RadToDeg();
1330 phi *= TMath::RadToDeg();
1331 psi = TMath::ATan2(m[6],m[7])*TMath::RadToDeg();
1332}
1333
1334////////////////////////////////////////////////////////////////////////////////
1335/// computes determinant of the rotation matrix
1336
1338{
1339 Double_t
1346 return det;
1347}
1348
1349////////////////////////////////////////////////////////////////////////////////
1350/// performes an orthogonality check and finds if the matrix is a reflection
1351/// Warning("CheckMatrix", "orthogonality check not performed yet");
1352
1354{
1355 if (Determinant() < 0) SetBit(kGeoReflection);
1357 if (TMath::Abs(dd) < 1.E-12) ResetBit(kGeoRotation);
1358 else SetBit(kGeoRotation);
1359}
1360
1361////////////////////////////////////////////////////////////////////////////////
1362/// Get the inverse rotation matrix (which is simply the transpose)
1363
1365{
1366 if (!invmat) {
1367 Error("GetInverse", "no place to store the inverse matrix");
1368 return;
1369 }
1370 for (Int_t i=0; i<3; i++) {
1371 for (Int_t j=0; j<3; j++) {
1372 invmat[3*i+j] = fRotationMatrix[3*j+i];
1373 }
1374 }
1375}
1376
1377////////////////////////////////////////////////////////////////////////////////
1378/// Multiply this rotation with the one specified by ROT.
1379/// - after=TRUE (default): THIS*ROT
1380/// - after=FALSE : ROT*THIS
1381
1383{
1384 const Double_t *matleft, *matright;
1386 Double_t newmat[9] = {0};
1387 if (after) {
1388 matleft = &fRotationMatrix[0];
1389 matright = rot->GetRotationMatrix();
1390 } else {
1391 matleft = rot->GetRotationMatrix();
1392 matright = &fRotationMatrix[0];
1393 }
1394 for (Int_t i=0; i<3; i++) {
1395 for (Int_t j=0; j<3; j++) {
1396 for (Int_t k=0; k<3; k++) {
1397 newmat[3*i+j] += matleft[3*i+k] * matright[3*k+j];
1398 }
1399 }
1400 }
1401 memcpy(&fRotationMatrix[0], &newmat[0], kN9);
1402}
1403
1404/** \class TGeoScale
1405\ingroup Geometry_classes
1406Class describing scale transformations. A scale is an
1407array of 3 doubles (sx, sy, sz) multiplying elements 0, 5 and 10
1408of the homogenous matrix. A scale is normalized : sx*sy*sz = 1
1409*/
1410
1412
1413////////////////////////////////////////////////////////////////////////////////
1414/// default constructor
1415
1417{
1419 for (Int_t i=0; i<3; i++) fScale[i] = 1.;
1420}
1421
1422////////////////////////////////////////////////////////////////////////////////
1423/// Copy constructor
1424
1426 :TGeoMatrix(other)
1427{
1428 SetScale(other);
1429}
1430
1431////////////////////////////////////////////////////////////////////////////////
1432/// Ctor. based on a general matrix
1433
1435 :TGeoMatrix(other)
1436{
1439 SetScale(other);
1440}
1441
1442////////////////////////////////////////////////////////////////////////////////
1443/// default constructor
1444
1446 :TGeoMatrix("")
1447{
1449 SetScale(sx, sy, sz);
1450}
1451
1452////////////////////////////////////////////////////////////////////////////////
1453/// default constructor
1454
1457{
1459 SetScale(sx, sy, sz);
1460}
1461
1462////////////////////////////////////////////////////////////////////////////////
1463/// destructor
1464
1466{
1467}
1468
1469////////////////////////////////////////////////////////////////////////////////
1470/// Assignment from a general matrix
1471
1473{
1474 if (&matrix == this) return *this;
1475 Bool_t registered = TestBit(kGeoRegistered);
1476 TNamed::operator=(matrix);
1477 SetScale(matrix);
1478 SetBit(kGeoRegistered,registered);
1481 return *this;
1482}
1483
1484////////////////////////////////////////////////////////////////////////////////
1485/// Scale composition
1486
1488{
1489 const Double_t *scl = right.GetScale();
1490 fScale[0] *= scl[0];
1491 fScale[1] *= scl[1];
1492 fScale[2] *= scl[2];
1493 SetBit(kGeoReflection, fScale[0] * fScale[1] * fScale[2] < 0);
1494 if (!IsScale()) SetBit(kGeoScale, right.IsScale());
1495 return *this;
1496}
1497
1499{
1500 TGeoScale s = *this;
1501 s *= right;
1502 return s;
1503}
1504
1506{
1507 TGeoHMatrix t = *this;
1508 t *= right;
1509 return t;
1510}
1511
1512////////////////////////////////////////////////////////////////////////////////
1513/// Is-equal operator
1514
1516{
1517 if (&other == this) return kTRUE;
1518 const Double_t *scl = GetScale();
1519 const Double_t *oscl = other.GetScale();
1520 for (auto i=0; i<3; i++)
1521 if (TMath::Abs(scl[i]-oscl[i])>1.E-10) return kFALSE;
1522 return kTRUE;
1523}
1524
1525////////////////////////////////////////////////////////////////////////////////
1526/// Return a temporary inverse of this.
1527
1529{
1530 TGeoHMatrix h;
1531 h = *this;
1533 Double_t scale[3];
1534 scale[0] = 1./fScale[0];
1535 scale[1] = 1./fScale[1];
1536 scale[2] = 1./fScale[2];
1537 h.SetScale(scale);
1538 return h;
1539}
1540
1541////////////////////////////////////////////////////////////////////////////////
1542/// scale setter
1543
1545{
1546 if (TMath::Abs(sx*sy*sz) < 1.E-10) {
1547 Error("SetScale", "Invalid scale %f, %f, %f for transformation %s",sx,sy,sx,GetName());
1548 return;
1549 }
1550 fScale[0] = sx;
1551 fScale[1] = sy;
1552 fScale[2] = sz;
1553 if (sx*sy*sz<0) SetBit(kGeoReflection);
1555}
1556
1557////////////////////////////////////////////////////////////////////////////////
1558/// Set scale from other transformation
1559
1561{
1562 SetBit(kGeoScale, other.IsScale());
1564 memcpy(fScale, other.GetScale(), kN3);
1565}
1566
1567////////////////////////////////////////////////////////////////////////////////
1568/// Convert a local point to the master frame.
1569
1570void TGeoScale::LocalToMaster(const Double_t *local, Double_t *master) const
1571{
1572 master[0] = local[0]*fScale[0];
1573 master[1] = local[1]*fScale[1];
1574 master[2] = local[2]*fScale[2];
1575}
1576
1577////////////////////////////////////////////////////////////////////////////////
1578/// Convert the local distance along unit vector DIR to master frame. If DIR
1579/// is not specified perform a conversion such as the returned distance is the
1580/// the minimum for all possible directions.
1581
1583{
1584 Double_t scale;
1585 if (!dir) {
1586 scale = TMath::Abs(fScale[0]);
1587 if (TMath::Abs(fScale[1])<scale) scale = TMath::Abs(fScale[1]);
1588 if (TMath::Abs(fScale[2])<scale) scale = TMath::Abs(fScale[2]);
1589 } else {
1590 scale = fScale[0]*fScale[0]*dir[0]*dir[0] +
1591 fScale[1]*fScale[1]*dir[1]*dir[1] +
1592 fScale[2]*fScale[2]*dir[2]*dir[2];
1593 scale = TMath::Sqrt(scale);
1594 }
1595 return scale*dist;
1596}
1597
1598////////////////////////////////////////////////////////////////////////////////
1599/// Make a clone of this matrix.
1600
1602{
1603 TGeoMatrix *matrix = new TGeoScale(*this);
1604 return matrix;
1605}
1606
1607////////////////////////////////////////////////////////////////////////////////
1608/// Convert a global point to local frame.
1609
1610void TGeoScale::MasterToLocal(const Double_t *master, Double_t *local) const
1611{
1612 local[0] = master[0]/fScale[0];
1613 local[1] = master[1]/fScale[1];
1614 local[2] = master[2]/fScale[2];
1615}
1616
1617////////////////////////////////////////////////////////////////////////////////
1618/// Convert the distance along unit vector DIR to local frame. If DIR
1619/// is not specified perform a conversion such as the returned distance is the
1620/// the minimum for all possible directions.
1621
1623{
1624 Double_t scale;
1625 if (!dir) {
1626 scale = TMath::Abs(fScale[0]);
1627 if (TMath::Abs(fScale[1])>scale) scale = TMath::Abs(fScale[1]);
1628 if (TMath::Abs(fScale[2])>scale) scale = TMath::Abs(fScale[2]);
1629 scale = 1./scale;
1630 } else {
1631 scale = (dir[0]*dir[0])/(fScale[0]*fScale[0]) +
1632 (dir[1]*dir[1])/(fScale[1]*fScale[1]) +
1633 (dir[2]*dir[2])/(fScale[2]*fScale[2]);
1634 scale = TMath::Sqrt(scale);
1635 }
1636 return scale*dist;
1637}
1638
1639/** \class TGeoCombiTrans
1640\ingroup Geometry_classes
1641Class describing rotation + translation. Most frequently used in the description
1642of TGeoNode 's
1643*/
1644
1646
1647////////////////////////////////////////////////////////////////////////////////
1648/// dummy ctor
1649
1651{
1652 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
1653 fRotation = 0;
1654}
1655
1656////////////////////////////////////////////////////////////////////////////////
1657/// Copy ctor from generic matrix.
1658
1660 :TGeoMatrix(other)
1661{
1663 if (other.IsTranslation()) {
1665 memcpy(fTranslation,other.GetTranslation(),kN3);
1666 } else {
1667 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
1668 }
1669 if (other.IsRotation()) {
1672 fRotation = new TGeoRotation(other);
1673 }
1674 else fRotation = 0;
1675}
1676
1677////////////////////////////////////////////////////////////////////////////////
1678/// Constructor from a translation and a rotation.
1679
1681{
1682 if (tr.IsTranslation()) {
1684 const Double_t *trans = tr.GetTranslation();
1685 memcpy(fTranslation, trans, kN3);
1686 } else {
1687 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
1688 }
1689 if (rot.IsRotation()) {
1692 fRotation = new TGeoRotation(rot);
1694 }
1695 else fRotation = 0;
1696}
1697
1698////////////////////////////////////////////////////////////////////////////////
1699/// Named ctor.
1700
1703{
1704 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
1705 fRotation = 0;
1706}
1707
1708////////////////////////////////////////////////////////////////////////////////
1709/// Constructor from a translation specified by X,Y,Z and a pointer to a rotation. The rotation will not be owned by this.
1710
1712 :TGeoMatrix("")
1713{
1714 SetTranslation(dx, dy, dz);
1715 fRotation = 0;
1716 SetRotation(rot);
1717}
1718
1719////////////////////////////////////////////////////////////////////////////////
1720/// Named ctor
1721
1724{
1725 SetTranslation(dx, dy, dz);
1726 fRotation = 0;
1727 SetRotation(rot);
1728}
1729
1730////////////////////////////////////////////////////////////////////////////////
1731/// Assignment operator with generic matrix.
1732
1734{
1735 if (&matrix == this) return *this;
1736 Bool_t registered = TestBit(kGeoRegistered);
1737 Clear();
1738 TNamed::operator=(matrix);
1739
1740 if (matrix.IsTranslation()) {
1741 memcpy(fTranslation,matrix.GetTranslation(),kN3);
1742 }
1743 if (matrix.IsRotation()) {
1744 if (!fRotation) {
1745 fRotation = new TGeoRotation();
1747 } else {
1748 if (!TestBit(kGeoMatrixOwned)) {
1749 fRotation = new TGeoRotation();
1751 }
1752 }
1756 } else {
1759 fRotation = 0;
1760 }
1761 SetBit(kGeoRegistered,registered);
1763 return *this;
1764}
1765
1766////////////////////////////////////////////////////////////////////////////////
1767/// Is-equal operator
1768
1770{
1771 if (&other == this) return kTRUE;
1772 const Double_t *tr = GetTranslation();
1773 const Double_t *otr = other.GetTranslation();
1774 for (auto i=0; i<3; i++) if (TMath::Abs(tr[i]-otr[i])>1.E-10) return kFALSE;
1775 const Double_t *rot = GetRotationMatrix();
1776 const Double_t *orot = other.GetRotationMatrix();
1777 for (auto i=0; i<9; i++) if (TMath::Abs(rot[i]-orot[i])>1.E-10) return kFALSE;
1778 return kTRUE;
1779}
1780
1781////////////////////////////////////////////////////////////////////////////////
1782/// Composition
1783
1785{
1786 Multiply(&right);
1787 return *this;
1788}
1789
1791{
1792 TGeoHMatrix h = *this;
1793 h *= right;
1794 return h;
1795}
1796
1797////////////////////////////////////////////////////////////////////////////////
1798/// destructor
1799
1801{
1802 if (fRotation) {
1804 }
1805}
1806
1807////////////////////////////////////////////////////////////////////////////////
1808/// Reset translation/rotation to identity
1809
1811{
1812 if (IsTranslation()) {
1814 memset(fTranslation, 0, kN3);
1815 }
1816 if (fRotation) {
1817 if (TestBit(kGeoMatrixOwned)) delete fRotation;
1818 fRotation = 0;
1819 }
1823}
1824
1825////////////////////////////////////////////////////////////////////////////////
1826/// Return a temporary inverse of this.
1827
1829{
1830 TGeoHMatrix h;
1831 h = *this;
1833 Bool_t is_tr = IsTranslation();
1834 Bool_t is_rot = IsRotation();
1835 Double_t tr[3];
1836 Double_t newrot[9];
1837 const Double_t *rot = GetRotationMatrix();
1838 tr[0] = -fTranslation[0]*rot[0] - fTranslation[1]*rot[3] - fTranslation[2]*rot[6];
1839 tr[1] = -fTranslation[0]*rot[1] - fTranslation[1]*rot[4] - fTranslation[2]*rot[7];
1840 tr[2] = -fTranslation[0]*rot[2] - fTranslation[1]*rot[5] - fTranslation[2]*rot[8];
1841 h.SetTranslation(tr);
1842 newrot[0] = rot[0];
1843 newrot[1] = rot[3];
1844 newrot[2] = rot[6];
1845 newrot[3] = rot[1];
1846 newrot[4] = rot[4];
1847 newrot[5] = rot[7];
1848 newrot[6] = rot[2];
1849 newrot[7] = rot[5];
1850 newrot[8] = rot[8];
1851 h.SetRotation(newrot);
1852 h.SetBit(kGeoTranslation,is_tr);
1853 h.SetBit(kGeoRotation,is_rot);
1854 return h;
1855}
1856
1857////////////////////////////////////////////////////////////////////////////////
1858/// Make a clone of this matrix.
1859
1861{
1862 TGeoMatrix *matrix = new TGeoCombiTrans(*this);
1863 return matrix;
1864}
1865
1866////////////////////////////////////////////////////////////////////////////////
1867/// multiply to the right with an other transformation
1868/// if right is identity matrix, just return
1869
1871{
1872 if (right->IsIdentity()) return;
1873 TGeoHMatrix h = *this;
1874 h.Multiply(right);
1875 operator=(h);
1876}
1877
1878////////////////////////////////////////////////////////////////////////////////
1879/// Register the matrix in the current manager, which will become the owner.
1880
1882{
1885}
1886
1887////////////////////////////////////////////////////////////////////////////////
1888/// Rotate about X axis with angle expressed in degrees.
1889
1891{
1892 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
1894 else fRotation = new TGeoRotation();
1896 }
1898 const Double_t *rot = fRotation->GetRotationMatrix();
1900 Double_t c = TMath::Cos(phi);
1901 Double_t s = TMath::Sin(phi);
1902 Double_t v[9];
1903 v[0] = rot[0];
1904 v[1] = rot[1];
1905 v[2] = rot[2];
1906 v[3] = c*rot[3]-s*rot[6];
1907 v[4] = c*rot[4]-s*rot[7];
1908 v[5] = c*rot[5]-s*rot[8];
1909 v[6] = s*rot[3]+c*rot[6];
1910 v[7] = s*rot[4]+c*rot[7];
1911 v[8] = s*rot[5]+c*rot[8];
1914 if (!IsTranslation()) return;
1915 v[0] = fTranslation[0];
1916 v[1] = c*fTranslation[1]-s*fTranslation[2];
1917 v[2] = s*fTranslation[1]+c*fTranslation[2];
1918 memcpy(fTranslation,v,kN3);
1919}
1920
1921////////////////////////////////////////////////////////////////////////////////
1922/// Rotate about Y axis with angle expressed in degrees.
1923
1925{
1926 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
1928 else fRotation = new TGeoRotation();
1930 }
1932 const Double_t *rot = fRotation->GetRotationMatrix();
1934 Double_t c = TMath::Cos(phi);
1935 Double_t s = TMath::Sin(phi);
1936 Double_t v[9];
1937 v[0] = c*rot[0]+s*rot[6];
1938 v[1] = c*rot[1]+s*rot[7];
1939 v[2] = c*rot[2]+s*rot[8];
1940 v[3] = rot[3];
1941 v[4] = rot[4];
1942 v[5] = rot[5];
1943 v[6] = -s*rot[0]+c*rot[6];
1944 v[7] = -s*rot[1]+c*rot[7];
1945 v[8] = -s*rot[2]+c*rot[8];
1948 if (!IsTranslation()) return;
1949 v[0] = c*fTranslation[0]+s*fTranslation[2];
1950 v[1] = fTranslation[1];
1951 v[2] = -s*fTranslation[0]+c*fTranslation[2];
1952 memcpy(fTranslation,v,kN3);
1953}
1954
1955////////////////////////////////////////////////////////////////////////////////
1956/// Rotate about Z axis with angle expressed in degrees.
1957
1959{
1960 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
1962 else fRotation = new TGeoRotation();
1964 }
1966 const Double_t *rot = fRotation->GetRotationMatrix();
1968 Double_t c = TMath::Cos(phi);
1969 Double_t s = TMath::Sin(phi);
1970 Double_t v[9];
1971 v[0] = c*rot[0]-s*rot[3];
1972 v[1] = c*rot[1]-s*rot[4];
1973 v[2] = c*rot[2]-s*rot[5];
1974 v[3] = s*rot[0]+c*rot[3];
1975 v[4] = s*rot[1]+c*rot[4];
1976 v[5] = s*rot[2]+c*rot[5];
1977 v[6] = rot[6];
1978 v[7] = rot[7];
1979 v[8] = rot[8];
1982 if (!IsTranslation()) return;
1983 v[0] = c*fTranslation[0]-s*fTranslation[1];
1984 v[1] = s*fTranslation[0]+c*fTranslation[1];
1985 v[2] = fTranslation[2];
1986 memcpy(fTranslation,v,kN3);
1987}
1988
1989////////////////////////////////////////////////////////////////////////////////
1990/// Multiply by a reflection respect to YZ.
1991
1993{
1994 if (leftside && !rotonly) fTranslation[0] = -fTranslation[0];
1995 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
1997 else fRotation = new TGeoRotation();
1999 }
2001 fRotation->ReflectX(leftside);
2003}
2004
2005////////////////////////////////////////////////////////////////////////////////
2006/// Multiply by a reflection respect to ZX.
2007
2009{
2010 if (leftside && !rotonly) fTranslation[1] = -fTranslation[1];
2011 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
2013 else fRotation = new TGeoRotation();
2015 }
2017 fRotation->ReflectY(leftside);
2019}
2020
2021////////////////////////////////////////////////////////////////////////////////
2022/// Multiply by a reflection respect to XY.
2023
2025{
2026 if (leftside && !rotonly) fTranslation[2] = -fTranslation[2];
2027 if (!fRotation || !TestBit(kGeoMatrixOwned)) {
2029 else fRotation = new TGeoRotation();
2031 }
2033 fRotation->ReflectZ(leftside);
2035}
2036
2037////////////////////////////////////////////////////////////////////////////////
2038/// Save a primitive as a C++ statement(s) on output stream "out".
2039
2040void TGeoCombiTrans::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
2041{
2042 if (TestBit(kGeoSavePrimitive)) return;
2043 out << " // Combi transformation: " << GetName() << std::endl;
2044 out << " dx = " << fTranslation[0] << ";" << std::endl;
2045 out << " dy = " << fTranslation[1] << ";" << std::endl;
2046 out << " dz = " << fTranslation[2] << ";" << std::endl;
2047 if (fRotation && fRotation->IsRotation()) {
2049 out << " auto " << GetPointerName() << " = new TGeoCombiTrans(\"" << GetName() << "\", dx, dy, dz, "
2050 << fRotation->GetPointerName() << ");" << std::endl;
2051 } else {
2052 out << " auto " << GetPointerName() << " = new TGeoCombiTrans(\"" << GetName() << "\");" << std::endl;
2053 out << " " << GetPointerName() << "->SetTranslation(dx, dy, dz);" << std::endl;
2054 }
2056}
2057
2058////////////////////////////////////////////////////////////////////////////////
2059/// Assign a foreign rotation to the combi. The rotation is NOT owned by this.
2060
2062{
2064 fRotation = 0;
2068 if (!rot) return;
2069 if (!rot->IsRotation()) return;
2070
2073 TGeoRotation *rr = (TGeoRotation*)rot;
2074 fRotation = rr;
2075}
2076
2077////////////////////////////////////////////////////////////////////////////////
2078/// Copy the rotation from another one.
2079
2081{
2083 fRotation = 0;
2084 if (!rot.IsRotation()) {
2088 return;
2089 }
2090
2093 fRotation = new TGeoRotation(rot);
2095}
2096
2097////////////////////////////////////////////////////////////////////////////////
2098/// copy the translation component
2099
2101{
2102 if (tr.IsTranslation()) {
2104 const Double_t *trans = tr.GetTranslation();
2105 memcpy(fTranslation, trans, kN3);
2106 } else {
2107 if (!IsTranslation()) return;
2108 memset(fTranslation, 0, kN3);
2110 }
2111}
2112
2113////////////////////////////////////////////////////////////////////////////////
2114/// set the translation component
2115
2117{
2118 fTranslation[0] = dx;
2119 fTranslation[1] = dy;
2120 fTranslation[2] = dz;
2123}
2124
2125////////////////////////////////////////////////////////////////////////////////
2126/// set the translation component
2127
2129{
2130 fTranslation[0] = vect[0];
2131 fTranslation[1] = vect[1];
2132 fTranslation[2] = vect[2];
2135}
2136
2137////////////////////////////////////////////////////////////////////////////////
2138/// get the rotation array
2139
2141{
2142 if (!fRotation) return kIdentityMatrix;
2143 return fRotation->GetRotationMatrix();
2144}
2145
2146/** \class TGeoGenTrans
2147\ingroup Geometry_classes
2148Most general transformation, holding a translation, a rotation and a scale
2149*/
2150
2152
2153////////////////////////////////////////////////////////////////////////////////
2154/// dummy ctor
2155
2157{
2159 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
2160 for (Int_t j=0; j<3; j++) fScale[j] = 1.0;
2161 fRotation = 0;
2162}
2163
2164////////////////////////////////////////////////////////////////////////////////
2165/// constructor
2166
2169{
2171 for (Int_t i=0; i<3; i++) fTranslation[i] = 0.0;
2172 for (Int_t j=0; j<3; j++) fScale[j] = 1.0;
2173 fRotation = 0;
2174}
2175
2176////////////////////////////////////////////////////////////////////////////////
2177/// constructor
2178
2180 Double_t sx, Double_t sy, Double_t sz, TGeoRotation *rot)
2181 :TGeoCombiTrans("")
2182{
2184 SetTranslation(dx, dy, dz);
2185 SetScale(sx, sy, sz);
2186 SetRotation(rot);
2187}
2188
2189////////////////////////////////////////////////////////////////////////////////
2190/// constructor
2191
2193 Double_t sx, Double_t sy, Double_t sz, TGeoRotation *rot)
2195{
2197 SetTranslation(dx, dy, dz);
2198 SetScale(sx, sy, sz);
2199 SetRotation(rot);
2200}
2201
2202////////////////////////////////////////////////////////////////////////////////
2203/// destructor
2204
2206{
2207}
2208
2209////////////////////////////////////////////////////////////////////////////////
2210/// clear the fields of this transformation
2211
2213{
2214 memset(&fTranslation[0], 0, kN3);
2215 memset(&fScale[0], 0, kN3);
2216 if (fRotation) fRotation->Clear();
2217}
2218
2219////////////////////////////////////////////////////////////////////////////////
2220/// set the scale
2221
2223{
2224 if (sx<1.E-5 || sy<1.E-5 || sz<1.E-5) {
2225 Error("ctor", "Invalid scale");
2226 return;
2227 }
2228 fScale[0] = sx;
2229 fScale[1] = sy;
2230 fScale[2] = sz;
2231}
2232
2233////////////////////////////////////////////////////////////////////////////////
2234/// Return a temporary inverse of this.
2235
2237{
2238 TGeoHMatrix h = *this;
2240 return h;
2241}
2242
2243////////////////////////////////////////////////////////////////////////////////
2244/// A scale transformation should be normalized by sx*sy*sz factor
2245
2247{
2248 Double_t normfactor = fScale[0]*fScale[1]*fScale[2];
2249 if (normfactor <= 1E-5) return kFALSE;
2250 for (Int_t i=0; i<3; i++)
2251 fScale[i] /= normfactor;
2252 return kTRUE;
2253}
2254
2255/** \class TGeoIdentity
2256\ingroup Geometry_classes
2257An identity transformation. It holds no data member
2258and returns pointers to static null translation and identity
2259transformations for rotation and scale
2260*/
2261
2263
2264////////////////////////////////////////////////////////////////////////////////
2265/// dummy ctor
2266
2268{
2269 if (!gGeoIdentity) gGeoIdentity = this;
2271}
2272
2273////////////////////////////////////////////////////////////////////////////////
2274/// constructor
2275
2278{
2279 if (!gGeoIdentity) gGeoIdentity = this;
2281}
2282
2283////////////////////////////////////////////////////////////////////////////////
2284/// Return a temporary inverse of this.
2285
2287{
2289 return h;
2290}
2291
2292/** \class TGeoHMatrix
2293\ingroup Geometry_classes
2294
2295Matrix class used for computing global transformations
2296Should NOT be used for node definition. An instance of this class
2297is generally used to pile-up local transformations starting from
2298the top level physical node, down to the current node.
2299*/
2300
2302
2303////////////////////////////////////////////////////////////////////////////////
2304/// dummy ctor
2305
2307{
2308 memset(&fTranslation[0], 0, kN3);
2310 memcpy(fScale,kUnitScale,kN3);
2311}
2312
2313////////////////////////////////////////////////////////////////////////////////
2314/// constructor
2315
2318{
2319 memset(&fTranslation[0], 0, kN3);
2321 memcpy(fScale,kUnitScale,kN3);
2322}
2323
2324////////////////////////////////////////////////////////////////////////////////
2325/// assignment
2326
2328 :TGeoMatrix(matrix)
2329{
2330 memset(&fTranslation[0], 0, kN3);
2332 memcpy(fScale,kUnitScale,kN3);
2333 if (matrix.IsIdentity()) return;
2334 if (matrix.IsTranslation())
2336 if (matrix.IsRotation())
2337 memcpy(fRotationMatrix,matrix.GetRotationMatrix(),kN9);
2338 if (matrix.IsScale())
2339 memcpy(fScale,matrix.GetScale(),kN3);
2340}
2341
2342////////////////////////////////////////////////////////////////////////////////
2343/// destructor
2344
2346{
2347}
2348
2349////////////////////////////////////////////////////////////////////////////////
2350/// assignment
2351
2353{
2354 return TGeoHMatrix::operator=(*matrix);
2355}
2356
2357////////////////////////////////////////////////////////////////////////////////
2358/// assignment
2359
2361{
2362 if (&matrix == this) return *this;
2363 Clear();
2364 Bool_t registered = TestBit(kGeoRegistered);
2365 TNamed::operator=(matrix);
2366 if (matrix.IsIdentity()) return *this;
2367 if (matrix.IsTranslation())
2368 memcpy(fTranslation,matrix.GetTranslation(),kN3);
2369 if (matrix.IsRotation())
2370 memcpy(fRotationMatrix,matrix.GetRotationMatrix(),kN9);
2371 if (matrix.IsScale())
2372 memcpy(fScale,matrix.GetScale(),kN3);
2373 SetBit(kGeoRegistered,registered);
2374 return *this;
2375}
2376
2377////////////////////////////////////////////////////////////////////////////////
2378/// Composition
2379
2381{
2382 Multiply(&right);
2383 return *this;
2384}
2385
2387{
2388 TGeoHMatrix h = *this;
2389 h *= right;
2390 return h;
2391}
2392
2393////////////////////////////////////////////////////////////////////////////////
2394/// Is-equal operator
2395
2397{
2398 if (&other == this) return kTRUE;
2399 const Double_t *tr = GetTranslation();
2400 const Double_t *otr = other.GetTranslation();
2401 for (auto i=0; i<3; i++) if (TMath::Abs(tr[i]-otr[i])>1.E-10) return kFALSE;
2402 const Double_t *rot = GetRotationMatrix();
2403 const Double_t *orot = other.GetRotationMatrix();
2404 for (auto i=0; i<9; i++) if (TMath::Abs(rot[i]-orot[i])>1.E-10) return kFALSE;
2405 const Double_t *scl = GetScale();
2406 const Double_t *oscl = other.GetScale();
2407 for (auto i=0; i<3; i++) if (TMath::Abs(scl[i]-oscl[i])>1.E-10) return kFALSE;
2408 return kTRUE;
2409}
2410
2411////////////////////////////////////////////////////////////////////////////////
2412/// Fast copy method.
2413
2415{
2417 SetBit(kGeoRotation, other->IsRotation());
2419 memcpy(fTranslation,other->GetTranslation(),kN3);
2420 memcpy(fRotationMatrix,other->GetRotationMatrix(),kN9);
2421}
2422
2423////////////////////////////////////////////////////////////////////////////////
2424/// clear the data for this matrix
2425
2427{
2429 if (IsIdentity()) return;
2435 memcpy(fScale,kUnitScale,kN3);
2436}
2437
2438////////////////////////////////////////////////////////////////////////////////
2439/// Make a clone of this matrix.
2440
2442{
2443 TGeoMatrix *matrix = new TGeoHMatrix(*this);
2444 return matrix;
2445}
2446
2447////////////////////////////////////////////////////////////////////////////////
2448/// Perform a rotation about Z having the sine/cosine of the rotation angle.
2449
2451{
2452 fRotationMatrix[0] = sincos[1];
2453 fRotationMatrix[1] = -sincos[0];
2454 fRotationMatrix[3] = sincos[0];
2455 fRotationMatrix[4] = sincos[1];
2457}
2458
2459////////////////////////////////////////////////////////////////////////////////
2460/// Return a temporary inverse of this.
2461
2463{
2464 TGeoHMatrix h;
2465 h = *this;
2467 if (IsTranslation()) {
2468 Double_t tr[3];
2472 h.SetTranslation(tr);
2473 }
2474 if (IsRotation()) {
2475 Double_t newrot[9];
2476 newrot[0] = fRotationMatrix[0];
2477 newrot[1] = fRotationMatrix[3];
2478 newrot[2] = fRotationMatrix[6];
2479 newrot[3] = fRotationMatrix[1];
2480 newrot[4] = fRotationMatrix[4];
2481 newrot[5] = fRotationMatrix[7];
2482 newrot[6] = fRotationMatrix[2];
2483 newrot[7] = fRotationMatrix[5];
2484 newrot[8] = fRotationMatrix[8];
2485 h.SetRotation(newrot);
2486 }
2487 if (IsScale()) {
2488 Double_t sc[3];
2489 sc[0] = 1./fScale[0];
2490 sc[1] = 1./fScale[1];
2491 sc[2] = 1./fScale[2];
2492 h.SetScale(sc);
2493 }
2494 return h;
2495}
2496
2497////////////////////////////////////////////////////////////////////////////////
2498/// computes determinant of the rotation matrix
2499
2501{
2502 Double_t
2509 return det;
2510}
2511
2512////////////////////////////////////////////////////////////////////////////////
2513/// multiply to the right with an other transformation
2514/// if right is identity matrix, just return
2515
2517{
2518 if (right->IsIdentity()) return;
2519 const Double_t *r_tra = right->GetTranslation();
2520 const Double_t *r_rot = right->GetRotationMatrix();
2521 const Double_t *r_scl = right->GetScale();
2522 if (IsIdentity()) {
2523 if (right->IsRotation()) {
2525 memcpy(fRotationMatrix,r_rot,kN9);
2527 }
2528 if (right->IsScale()) {
2530 memcpy(fScale,r_scl,kN3);
2531 }
2532 if (right->IsTranslation()) {
2534 memcpy(fTranslation,r_tra,kN3);
2535 }
2536 return;
2537 }
2538 Int_t i, j;
2539 Double_t new_rot[9];
2540
2541 if (right->IsRotation()) {
2544 }
2545 if (right->IsScale()) SetBit(kGeoScale);
2546 if (right->IsTranslation()) SetBit(kGeoTranslation);
2547
2548 // new translation
2549 if (IsTranslation()) {
2550 for (i=0; i<3; i++) {
2551 fTranslation[i] += fRotationMatrix[3*i]*r_tra[0]
2552 + fRotationMatrix[3*i+1]*r_tra[1]
2553 + fRotationMatrix[3*i+2]*r_tra[2];
2554 }
2555 }
2556 if (IsRotation()) {
2557 // new rotation
2558 for (i=0; i<3; i++) {
2559 for (j=0; j<3; j++) {
2560 new_rot[3*i+j] = fRotationMatrix[3*i]*r_rot[j] +
2561 fRotationMatrix[3*i+1]*r_rot[3+j] +
2562 fRotationMatrix[3*i+2]*r_rot[6+j];
2563 }
2564 }
2565 memcpy(fRotationMatrix,new_rot,kN9);
2566 }
2567 // new scale
2568 if (IsScale()) {
2569 for (i=0; i<3; i++) fScale[i] *= r_scl[i];
2570 }
2571}
2572
2573////////////////////////////////////////////////////////////////////////////////
2574/// multiply to the left with an other transformation
2575/// if right is identity matrix, just return
2576
2578{
2579 if (left == gGeoIdentity) return;
2580 const Double_t *l_tra = left->GetTranslation();
2581 const Double_t *l_rot = left->GetRotationMatrix();
2582 const Double_t *l_scl = left->GetScale();
2583 if (IsIdentity()) {
2584 if (left->IsRotation()) {
2587 memcpy(fRotationMatrix,l_rot,kN9);
2588 }
2589 if (left->IsScale()) {
2591 memcpy(fScale,l_scl,kN3);
2592 }
2593 if (left->IsTranslation()) {
2595 memcpy(fTranslation,l_tra,kN3);
2596 }
2597 return;
2598 }
2599 Int_t i, j;
2600 Double_t new_tra[3];
2601 Double_t new_rot[9];
2602
2603 if (left->IsRotation()) {
2606 }
2607 if (left->IsScale()) SetBit(kGeoScale);
2608 if (left->IsTranslation()) SetBit(kGeoTranslation);
2609
2610 // new translation
2611 if (IsTranslation()) {
2612 for (i=0; i<3; i++) {
2613 new_tra[i] = l_tra[i]
2614 + l_rot[3*i]* fTranslation[0]
2615 + l_rot[3*i+1]*fTranslation[1]
2616 + l_rot[3*i+2]*fTranslation[2];
2617 }
2618 memcpy(fTranslation,new_tra,kN3);
2619 }
2620 if (IsRotation()) {
2621 // new rotation
2622 for (i=0; i<3; i++) {
2623 for (j=0; j<3; j++) {
2624 new_rot[3*i+j] = l_rot[3*i]*fRotationMatrix[j] +
2625 l_rot[3*i+1]*fRotationMatrix[3+j] +
2626 l_rot[3*i+2]*fRotationMatrix[6+j];
2627 }
2628 }
2629 memcpy(fRotationMatrix,new_rot,kN9);
2630 }
2631 // new scale
2632 if (IsScale()) {
2633 for (i=0; i<3; i++) fScale[i] *= l_scl[i];
2634 }
2635}
2636
2637////////////////////////////////////////////////////////////////////////////////
2638/// Rotate about X axis with angle expressed in degrees.
2639
2641{
2644 Double_t c = TMath::Cos(phi);
2645 Double_t s = TMath::Sin(phi);
2646 Double_t v[9];
2647 v[0] = fRotationMatrix[0];
2648 v[1] = fRotationMatrix[1];
2649 v[2] = fRotationMatrix[2];
2650 v[3] = c*fRotationMatrix[3]-s*fRotationMatrix[6];
2651 v[4] = c*fRotationMatrix[4]-s*fRotationMatrix[7];
2652 v[5] = c*fRotationMatrix[5]-s*fRotationMatrix[8];
2653 v[6] = s*fRotationMatrix[3]+c*fRotationMatrix[6];
2654 v[7] = s*fRotationMatrix[4]+c*fRotationMatrix[7];
2655 v[8] = s*fRotationMatrix[5]+c*fRotationMatrix[8];
2656 memcpy(fRotationMatrix, v, kN9);
2657
2658 v[0] = fTranslation[0];
2659 v[1] = c*fTranslation[1]-s*fTranslation[2];
2660 v[2] = s*fTranslation[1]+c*fTranslation[2];
2661 memcpy(fTranslation,v,kN3);
2662}
2663
2664////////////////////////////////////////////////////////////////////////////////
2665/// Rotate about Y axis with angle expressed in degrees.
2666
2668{
2671 Double_t c = TMath::Cos(phi);
2672 Double_t s = TMath::Sin(phi);
2673 Double_t v[9];
2674 v[0] = c*fRotationMatrix[0]+s*fRotationMatrix[6];
2675 v[1] = c*fRotationMatrix[1]+s*fRotationMatrix[7];
2676 v[2] = c*fRotationMatrix[2]+s*fRotationMatrix[8];
2677 v[3] = fRotationMatrix[3];
2678 v[4] = fRotationMatrix[4];
2679 v[5] = fRotationMatrix[5];
2680 v[6] = -s*fRotationMatrix[0]+c*fRotationMatrix[6];
2681 v[7] = -s*fRotationMatrix[1]+c*fRotationMatrix[7];
2682 v[8] = -s*fRotationMatrix[2]+c*fRotationMatrix[8];
2683 memcpy(fRotationMatrix, v, kN9);
2684
2685 v[0] = c*fTranslation[0]+s*fTranslation[2];
2686 v[1] = fTranslation[1];
2687 v[2] = -s*fTranslation[0]+c*fTranslation[2];
2688 memcpy(fTranslation,v,kN3);
2689}
2690
2691////////////////////////////////////////////////////////////////////////////////
2692/// Rotate about Z axis with angle expressed in degrees.
2693
2695{
2698 Double_t c = TMath::Cos(phi);
2699 Double_t s = TMath::Sin(phi);
2700 Double_t v[9];
2701 v[0] = c*fRotationMatrix[0]-s*fRotationMatrix[3];
2702 v[1] = c*fRotationMatrix[1]-s*fRotationMatrix[4];
2703 v[2] = c*fRotationMatrix[2]-s*fRotationMatrix[5];
2704 v[3] = s*fRotationMatrix[0]+c*fRotationMatrix[3];
2705 v[4] = s*fRotationMatrix[1]+c*fRotationMatrix[4];
2706 v[5] = s*fRotationMatrix[2]+c*fRotationMatrix[5];
2707 v[6] = fRotationMatrix[6];
2708 v[7] = fRotationMatrix[7];
2709 v[8] = fRotationMatrix[8];
2710 memcpy(&fRotationMatrix[0],v,kN9);
2711
2712 v[0] = c*fTranslation[0]-s*fTranslation[1];
2713 v[1] = s*fTranslation[0]+c*fTranslation[1];
2714 v[2] = fTranslation[2];
2715 memcpy(fTranslation,v,kN3);
2716}
2717
2718////////////////////////////////////////////////////////////////////////////////
2719/// Multiply by a reflection respect to YZ.
2720
2721void TGeoHMatrix::ReflectX(Bool_t leftside, Bool_t rotonly)
2722{
2723 if (leftside && !rotonly) fTranslation[0] = -fTranslation[0];
2724 if (leftside) {
2728 } else {
2732 }
2735}
2736
2737////////////////////////////////////////////////////////////////////////////////
2738/// Multiply by a reflection respect to ZX.
2739
2740void TGeoHMatrix::ReflectY(Bool_t leftside, Bool_t rotonly)
2741{
2742 if (leftside && !rotonly) fTranslation[1] = -fTranslation[1];
2743 if (leftside) {
2747 } else {
2751 }
2754}
2755
2756////////////////////////////////////////////////////////////////////////////////
2757/// Multiply by a reflection respect to XY.
2758
2759void TGeoHMatrix::ReflectZ(Bool_t leftside, Bool_t rotonly)
2760{
2761 if (leftside && !rotonly) fTranslation[2] = -fTranslation[2];
2762 if (leftside) {
2766 } else {
2770 }
2773}
2774
2775////////////////////////////////////////////////////////////////////////////////
2776/// Save a primitive as a C++ statement(s) on output stream "out".
2777
2778void TGeoHMatrix::SavePrimitive(std::ostream &out, Option_t * /*option*/ /*= ""*/)
2779{
2780 if (TestBit(kGeoSavePrimitive)) return;
2781 const Double_t *tr = fTranslation;
2782 const Double_t *rot = fRotationMatrix;
2783 out << " // HMatrix: " << GetName() << std::endl;
2784 out << " tr[0] = " << tr[0] << "; " << "tr[1] = " << tr[1] << "; " << "tr[2] = " << tr[2] << ";" << std::endl;
2785 out << " rot[0] =" << rot[0] << "; " << "rot[1] = " << rot[1] << "; " << "rot[2] = " << rot[2] << ";" << std::endl;
2786 out << " rot[3] =" << rot[3] << "; " << "rot[4] = " << rot[4] << "; " << "rot[5] = " << rot[5] << ";" << std::endl;
2787 out << " rot[6] =" << rot[6] << "; " << "rot[7] = " << rot[7] << "; " << "rot[8] = " << rot[8] << ";" << std::endl;
2788 const char *name = GetPointerName();
2789 out << " auto " << name << " = new TGeoHMatrix(\"" << GetName() << "\");" << std::endl;
2790 out << " " << name << "->SetTranslation(tr);" << std::endl;
2791 out << " " << name << "->SetRotation(rot);" << std::endl;
2792 if (IsTranslation()) out << " " << name << "->SetBit(TGeoMatrix::kGeoTranslation);" << std::endl;
2793 if (IsRotation()) out << " " << name << "->SetBit(TGeoMatrix::kGeoRotation);" << std::endl;
2794 if (IsReflection()) out << " " << name << "->SetBit(TGeoMatrix::kGeoReflection);" << std::endl;
2796}
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t option
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 GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t TPoint TPoint angle
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
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeoManager * gGeoManager
const Int_t kN3
TGeoIdentity * gGeoIdentity
const Int_t kN9
const Double_t kUnitScale[3]
Definition TGeoMatrix.h:30
const Double_t kIdentityMatrix[3 *3]
Definition TGeoMatrix.h:26
R__EXTERN TGeoIdentity * gGeoIdentity
Definition TGeoMatrix.h:478
const Double_t kNullVector[3]
Definition TGeoMatrix.h:24
Class describing rotation + translation.
Definition TGeoMatrix.h:292
void Multiply(const TGeoMatrix *right)
multiply to the right with an other transformation if right is identity matrix, just return
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
TGeoCombiTrans & operator*=(const TGeoMatrix &other)
Composition.
Bool_t operator==(const TGeoMatrix &other) const
Is-equal operator.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
virtual const Double_t * GetTranslation() const
Definition TGeoMatrix.h:336
TGeoCombiTrans()
dummy ctor
virtual TGeoMatrix * MakeClone() const
Make a clone of this matrix.
virtual ~TGeoCombiTrans()
destructor
virtual const Double_t * GetRotationMatrix() const
get the rotation array
virtual void RotateZ(Double_t angle)
Rotate about Z axis with angle expressed in degrees.
Double_t fTranslation[3]
Definition TGeoMatrix.h:294
virtual void RegisterYourself()
Register the matrix in the current manager, which will become the owner.
TGeoRotation * fRotation
Definition TGeoMatrix.h:295
void SetTranslation(const TGeoTranslation &tr)
copy the translation component
void SetRotation(const TGeoRotation &other)
Copy the rotation from another one.
virtual void ReflectX(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to YZ.
TGeoCombiTrans & operator=(const TGeoCombiTrans &other)
Definition TGeoMatrix.h:305
virtual void ReflectY(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to ZX.
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
TGeoCombiTrans operator*(const TGeoMatrix &other) const
void Clear(Option_t *option="")
Reset translation/rotation to identity.
virtual void RotateY(Double_t angle)
Rotate about Y axis with angle expressed in degrees.
virtual void RotateX(Double_t angle)
Rotate about X axis with angle expressed in degrees.
Most general transformation, holding a translation, a rotation and a scale.
Definition TGeoMatrix.h:351
virtual ~TGeoGenTrans()
destructor
Double_t fScale[3]
Definition TGeoMatrix.h:353
void Clear(Option_t *option="")
clear the fields of this transformation
Bool_t Normalize()
A scale transformation should be normalized by sx*sy*sz factor.
TGeoGenTrans()
dummy ctor
void SetScale(Double_t sx, Double_t sy, Double_t sz)
set the scale
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:421
TGeoHMatrix & operator*=(const TGeoMatrix &other)
Composition.
TGeoHMatrix()
dummy ctor
virtual void RotateX(Double_t angle)
Rotate about X axis with angle expressed in degrees.
virtual void RotateZ(Double_t angle)
Rotate about Z axis with angle expressed in degrees.
void MultiplyLeft(const TGeoMatrix *left)
multiply to the left with an other transformation if right is identity matrix, just return
virtual ~TGeoHMatrix()
destructor
Double_t Determinant() const
computes determinant of the rotation matrix
virtual void ReflectY(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to ZX.
void Clear(Option_t *option="")
clear the data for this matrix
void CopyFrom(const TGeoMatrix *other)
Fast copy method.
void FastRotZ(const Double_t *sincos)
Perform a rotation about Z having the sine/cosine of the rotation angle.
Double_t fTranslation[3]
Definition TGeoMatrix.h:423
Bool_t operator==(const TGeoMatrix &other) const
Is-equal operator.
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
virtual const Double_t * GetScale() const
Definition TGeoMatrix.h:469
Double_t fRotationMatrix[9]
Definition TGeoMatrix.h:424
virtual const Double_t * GetTranslation() const
Definition TGeoMatrix.h:467
virtual TGeoMatrix * MakeClone() const
Make a clone of this matrix.
virtual void RotateY(Double_t angle)
Rotate about Y axis with angle expressed in degrees.
TGeoHMatrix & operator=(const TGeoHMatrix &other)
Definition TGeoMatrix.h:434
void Multiply(const TGeoMatrix *right)
multiply to the right with an other transformation if right is identity matrix, just return
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
TGeoHMatrix operator*(const TGeoMatrix &other) const
void SetTranslation(const Double_t *vect)
Definition TGeoMatrix.h:462
Double_t fScale[3]
Definition TGeoMatrix.h:425
virtual const Double_t * GetRotationMatrix() const
Definition TGeoMatrix.h:468
virtual void ReflectX(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to YZ.
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
An identity transformation.
Definition TGeoMatrix.h:384
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
TGeoIdentity()
dummy ctor
TObjArray * GetListOfMatrices() const
void RegisterMatrix(const TGeoMatrix *matrix)
Register a matrix to the list of matrices.
void BombTranslation(const Double_t *tr, Double_t *bombtr)
Get the new 'bombed' translation vector according current exploded view mode.
void UnbombTranslation(const Double_t *tr, Double_t *bombtr)
Get the new 'unbombed' translation vector according current exploded view mode.
Bool_t IsCleaning() const
Geometrical transformation package.
Definition TGeoMatrix.h:41
virtual void LocalToMasterVect(const Double_t *local, Double_t *master) const
convert a vector by multiplying its column vector (x, y, z, 1) to matrix inverse
Bool_t IsScale() const
Definition TGeoMatrix.h:70
void SetDefaultName()
If no name was supplied in the ctor, the type of transformation is checked.
Bool_t IsGeneral() const
Definition TGeoMatrix.h:75
@ kGeoSavePrimitive
Definition TGeoMatrix.h:51
@ kGeoTranslation
Definition TGeoMatrix.h:46
@ kGeoMatrixOwned
Definition TGeoMatrix.h:52
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void MasterToLocalVect(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
virtual const Double_t * GetTranslation() const =0
Bool_t IsTranslation() const
Definition TGeoMatrix.h:67
Bool_t IsReflection() const
Definition TGeoMatrix.h:69
Bool_t IsRotation() const
Definition TGeoMatrix.h:68
virtual void LocalToMasterBomb(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
virtual void RegisterYourself()
Register the matrix in the current manager, which will become the owner.
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
virtual void MasterToLocalBomb(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
Bool_t IsRotAboutZ() const
Returns true if no rotation or the rotation is about Z axis.
void GetHomogenousMatrix(Double_t *hmat) const
The homogenous matrix associated with the transformation is used for piling up's and visualization.
TGeoMatrix()
dummy constructor
Bool_t IsOwned() const
Definition TGeoMatrix.h:72
virtual const Double_t * GetScale() const =0
static void Normalize(Double_t *vect)
Normalize a vector.
void Print(Option_t *option="") const
print the matrix in 4x4 format
Bool_t IsIdentity() const
Definition TGeoMatrix.h:66
const char * GetPointerName() const
Provide a pointer name containing uid.
Bool_t IsCombi() const
Definition TGeoMatrix.h:73
Bool_t IsRegistered() const
Definition TGeoMatrix.h:77
virtual ~TGeoMatrix()
Destructor.
Bool_t IsShared() const
Definition TGeoMatrix.h:71
virtual const Double_t * GetRotationMatrix() const =0
virtual Int_t GetByteCount() const
Get total size in bytes of this.
virtual void ReflectY(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to ZX.
virtual void ReflectX(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to YZ.
Class describing rotations.
Definition TGeoMatrix.h:175
virtual void RotateY(Double_t angle)
Rotate about Y axis of the master frame with angle expressed in degrees.
virtual void ReflectX(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to YZ.
TGeoRotation()
Default constructor.
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
void SetAngles(Double_t phi, Double_t theta, Double_t psi)
Set matrix elements according to Euler angles.
void Clear(Option_t *option="")
reset data members
virtual const Double_t * GetRotationMatrix() const
Definition TGeoMatrix.h:230
void MultiplyBy(const TGeoRotation *rot, Bool_t after=kTRUE)
Multiply this rotation with the one specified by ROT.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
Bool_t operator==(const TGeoRotation &other) const
Is-equal operator.
void SetMatrix(const Double_t *rot)
Definition TGeoMatrix.h:225
virtual TGeoMatrix * MakeClone() const
Make a clone of this matrix.
virtual void ReflectY(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to ZX.
TGeoRotation & operator*=(const TGeoRotation &other)
Composition.
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void RotateX(Double_t angle)
Rotate about X axis of the master frame with angle expressed in degrees.
void CheckMatrix()
performes an orthogonality check and finds if the matrix is a reflection Warning("CheckMatrix",...
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
Double_t GetPhiRotation(Bool_t fixX=kFALSE) const
Returns rotation angle about Z axis in degrees.
void FastRotZ(const Double_t *sincos)
Perform a rotation about Z having the sine/cosine of the rotation angle.
void GetInverse(Double_t *invmat) const
Get the inverse rotation matrix (which is simply the transpose)
Double_t Determinant() const
computes determinant of the rotation matrix
void GetAngles(Double_t &theta1, Double_t &phi1, Double_t &theta2, Double_t &phi2, Double_t &theta3, Double_t &phi3) const
Retrieve rotation angles.
Bool_t IsValid() const
Perform orthogonality test for rotation.
virtual void RotateZ(Double_t angle)
Rotate about Z axis of the master frame with angle expressed in degrees.
TGeoRotation operator*(const TGeoRotation &other) const
Double_t fRotationMatrix[3 *3]
Definition TGeoMatrix.h:177
void SetRotation(const TGeoMatrix &other)
Copy rotation elements from other rotation matrix.
TGeoRotation & operator=(const TGeoRotation &other)
Definition TGeoMatrix.h:191
Class describing scale transformations.
Definition TGeoMatrix.h:245
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
Convert a global point to local frame.
TGeoScale()
default constructor
TGeoScale & operator=(const TGeoScale &other)
Definition TGeoMatrix.h:256
Bool_t operator==(const TGeoScale &other) const
Is-equal operator.
void SetScale(Double_t sx, Double_t sy, Double_t sz)
scale setter
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
virtual ~TGeoScale()
destructor
Double_t fScale[3]
Definition TGeoMatrix.h:247
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
Convert a local point to the master frame.
virtual const Double_t * GetScale() const
Definition TGeoMatrix.h:279
TGeoScale operator*(const TGeoScale &other) const
TGeoScale & operator*=(const TGeoScale &other)
Scale composition.
virtual TGeoMatrix * MakeClone() const
Make a clone of this matrix.
Class describing translations.
Definition TGeoMatrix.h:122
virtual const Double_t * GetTranslation() const
Definition TGeoMatrix.h:160
Bool_t operator==(const TGeoTranslation &other) const
Is-equal operator.
TGeoTranslation & operator*=(const TGeoTranslation &other)
Translation composition.
virtual void MasterToLocalBomb(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
virtual void RotateX(Double_t angle)
Rotate about X axis of the master frame with angle expressed in degrees.
TGeoTranslation operator*(const TGeoTranslation &right) const
void Add(const TGeoTranslation *other)
Adding a translation to this one.
virtual void LocalToMasterBomb(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
virtual void LocalToMasterVect(const Double_t *local, Double_t *master) const
convert a vector to MARS
virtual void RotateZ(Double_t angle)
Rotate about Z axis of the master frame with angle expressed in degrees.
virtual void MasterToLocalVect(const Double_t *master, Double_t *local) const
convert a vector from MARS to local
virtual void RotateY(Double_t angle)
Rotate about Y axis of the master frame with angle expressed in degrees.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save a primitive as a C++ statement(s) on output stream "out".
void SetTranslation(Double_t dx, Double_t dy, Double_t dz)
Set translation components.
TGeoTranslation & operator=(const TGeoTranslation &other)
Definition TGeoMatrix.h:133
Double_t fTranslation[3]
Definition TGeoMatrix.h:124
virtual void LocalToMaster(const Double_t *local, Double_t *master) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix inverse
TGeoHMatrix Inverse() const
Return a temporary inverse of this.
virtual TGeoMatrix * MakeClone() const
Make a clone of this matrix.
TGeoTranslation()
Default constructor.
void Subtract(const TGeoTranslation *other)
Subtracting a translation from this one.
virtual void MasterToLocal(const Double_t *master, Double_t *local) const
convert a point by multiplying its column vector (x, y, z, 1) to matrix
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
TObject * Remove(TObject *obj) override
Remove object from array.
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual UInt_t GetUniqueID() const
Return the unique object id.
Definition TObject.cxx:457
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:956
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:774
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:970
void ResetBit(UInt_t f)
Definition TObject.h:200
Basic string class.
Definition TString.h:139
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition TString.cxx:2334
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:630
Double_t ASin(Double_t)
Returns the principal value of the arc sine of x, expressed in radians.
Definition TMath.h:622
Double_t ATan2(Double_t y, Double_t x)
Returns the principal value of the arc tangent of y/x, expressed in radians.
Definition TMath.h:644
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:660
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:592
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:586
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:72
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
auto * th3
Definition textalign.C:22
auto * th2
Definition textalign.C:18
auto * th1
Definition textalign.C:14
TMarker m
Definition textangle.C:8