Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeometry.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Rene Brun 22/09/95
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#include "TROOT.h"
13#include "TBuffer.h"
14#include "THashList.h"
15#include "TObjArray.h"
16#include "TGeometry.h"
17#include "TNode.h"
18#include "TMaterial.h"
19#include "TBrowser.h"
20
22
24
25/** \class TGeometry
26\ingroup g3d
27TGeometry description.
28
29The Geometry class describes the geometry of a detector.
30The current implementation supports the GEANT3 style description.
31A special program provided in the ROOT utilities (toroot) can be used
32to automatically translate a GEANT detector geometry into a ROOT geometry.
33
34a Geometry object is entered into the list of geometries into the
35ROOT main object (see TROOT description) when the TGeometry
36constructor is invoked.
37Several geometries may coexist in memory.
38/
39A Geometry object consist of the following linked lists:
40
41 - the TMaterial list (material definition only).
42 - the TRotmatrix list (Rotation matrices definition only).
43 - the TShape list (volume definition only).
44 - the TNode list assembling all detector elements.
45
46Only the Build and Draw functions for a geometry are currently supported.
47
48The conversion program from Geant to Root has been added in the list
49of utilities in utils directory.(see g2root)
50The executable module of g2root can be found in $ROOTSYS/bin/g2root.
51
52To use this conversion program, type the shell command:
53
54~~~ {.cpp}
55 g2root geant_rzfile macro_name
56~~~
57
58for example
59
60~~~ {.cpp}
61 g2root na49.geom na49.C
62~~~
63
64will convert the GEANT RZ file na49.geom into a ROOT macro na49.C
65
66To generate the Geometry structure within Root, do:
67
68~~~ {.cpp}
69 Root > .x na49.C
70 Root > na49.Draw()
71 Root > wh.x3d() (this invokes the 3-d Root viewer)
72 Root > TFile gna49("na49.root","NEW") //open a new root file
73 Root > na49.Write() //Write the na49 geometry structure
74 Root > gna49.Write() //Write all keys (in this case only one)
75~~~
76
77Note: all keys are also written on closing of the file, gna49.Close or
78when the program exits, Root closes all open files correctly.
79Once this file has been written, in a subsequent session, simply do:
80
81~~~ {.cpp}
82 Root > TFile gna49("na49.root")
83 Root > na49.Draw()
84~~~
85
86The figure below shows the geometry above using the x3d viewer.
87This x3d viewer is invoked by selecting "View x3d" in the View menu
88of a canvas (See example of this tool bar in TCanvas).
89
90\image html g3d_na49.png
91*/
92
93////////////////////////////////////////////////////////////////////////////////
94/// Geometry default constructor.
95
97{
98 fMaterials = new THashList(100,3);
99 fMatrices = new THashList(100,3);
100 fShapes = new THashList(500,3);
101 fNodes = new TList;
102 fCurrentNode = 0;
104 fMatrixPointer = 0;
105 fShapePointer = 0;
106 gGeometry = this;
107 fBomb = 1;
108 fMatrix = 0;
109 fX=fY=fZ =0.0;
110 fGeomLevel =0;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Geometry normal constructor.
116
117TGeometry::TGeometry(const char *name,const char *title ) : TNamed (name, title)
118{
119 fMaterials = new THashList(1000,3);
120 fMatrices = new THashList(1000,3);
121 fShapes = new THashList(5000,3);
122 fNodes = new TList;
123 fCurrentNode = 0;
125 fMatrixPointer = 0;
126 fShapePointer = 0;
127 gGeometry = this;
128 fBomb = 1;
129 fMatrix = 0;
130 fX=fY=fZ =0.0;
131 gROOT->GetListOfGeometries()->Add(this);
132 fGeomLevel =0;
134}
135
136////////////////////////////////////////////////////////////////////////////////
137/// copy constructor
138
140 TNamed(geo),
141 fMaterials(geo.fMaterials),
142 fMatrices(geo.fMatrices),
143 fShapes(geo.fShapes),
144 fNodes(geo.fNodes),
145 fMatrix(geo.fMatrix),
146 fCurrentNode(geo.fCurrentNode),
147 fMaterialPointer(geo.fMaterialPointer),
148 fMatrixPointer(geo.fMatrixPointer),
149 fShapePointer(geo.fShapePointer),
150 fBomb(geo.fBomb),
151 fGeomLevel(geo.fGeomLevel),
152 fX(geo.fX),
153 fY(geo.fY),
154 fZ(geo.fZ)
155{
156 for(Int_t i=0; i<kMAXLEVELS; i++) {
157 for(Int_t j=0; j<kVectorSize; j++)
158 fTranslation[i][j]=geo.fTranslation[i][j];
159 for(Int_t j=0; j<kMatrixSize; j++)
160 fRotMatrix[i][j]=geo.fRotMatrix[i][j];
162 }
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// assignment operator
167
169{
170 if(this!=&geo) {
174 fShapes=geo.fShapes;
175 fNodes=geo.fNodes;
176 fMatrix=geo.fMatrix;
181 fBomb=geo.fBomb;
183 fX=geo.fX;
184 fY=geo.fY;
185 fZ=geo.fZ;
186 for(Int_t i=0; i<kMAXLEVELS; i++) {
187 for(Int_t j=0; j<kVectorSize; j++)
188 fTranslation[i][j]=geo.fTranslation[i][j];
189 for(Int_t j=0; j<kMatrixSize; j++)
190 fRotMatrix[i][j]=geo.fRotMatrix[i][j];
192 }
193 }
194 return *this;
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// Geometry default destructor.
199
201{
202 if (!fMaterials) return;
204 fMatrices->Delete();
205 fShapes->Delete();
206 fNodes->Delete();
207 delete fMaterials;
208 delete fMatrices;
209 delete fShapes;
210 delete fNodes;
211 delete [] fMaterialPointer;
212 delete [] fMatrixPointer;
213 delete [] fShapePointer;
214 fMaterials = 0;
215 fMatrices = 0;
216 fShapes = 0;
217 fNodes = 0;
219 fMatrixPointer = 0;
220 fShapePointer = 0;
221
222 if (gGeometry == this) {
223 gGeometry = (TGeometry*) gROOT->GetListOfGeometries()->First();
224 if (gGeometry == this)
225 gGeometry = (TGeometry*) gROOT->GetListOfGeometries()->After(gGeometry);
226 }
227 gROOT->GetListOfGeometries()->Remove(this);
228}
229
230////////////////////////////////////////////////////////////////////////////////
231/// Browse.
232
234{
235 if( b ) {
236 b->Add( fMaterials, "Materials" );
237 b->Add( fMatrices, "Rotation Matrices" );
238 b->Add( fShapes, "Shapes" );
239 b->Add( fNodes, "Nodes" );
240 }
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Change Current Geometry to this.
245
246void TGeometry::cd(const char *)
247{
248 gGeometry = this;
249}
250
251////////////////////////////////////////////////////////////////////////////////
252/// Draw this Geometry.
253
255{
256 TNode *node1 = (TNode*)fNodes->First();
257 if (node1) node1->Draw(option);
258
259}
260
261////////////////////////////////////////////////////////////////////////////////
262/// Find object in a geometry node, material, etc
263
265{
266 Error("FindObject","Not yet implemented");
267 return 0;
268}
269
270////////////////////////////////////////////////////////////////////////////////
271/// Search object identified by name in the geometry tree
272
274{
276 if (loc) return loc->At(0);
277 return 0;
278}
279
280////////////////////////////////////////////////////////////////////////////////
281/// Static function called by TROOT to search name in the geometry.
282/// Returns a TObjArray containing a pointer to the found object
283/// and a pointer to the container where the object was found.
284
286{
287 static TObjArray *locs = 0;
288 if (!locs) locs = new TObjArray(2);
289 TObjArray &loc = *locs;
290 loc[0] = 0;
291 loc[1] = 0;
292
293 if (!gGeometry) return &loc;
294
295 TObject *temp;
296 TObject *where;
297
299 where = gGeometry->GetListOfMaterials();
300
301 if (!temp) {
303 where = gGeometry->GetListOfShapes();
304 }
305 if (!temp) {
307 where = gGeometry->GetListOfMatrices();
308 }
309 if (!temp) {
310 temp = gGeometry->GetNode(name);
311 where = gGeometry;
312 }
313 loc[0] = temp;
314 loc[1] = where;
315
316 return &loc;
317}
318
319////////////////////////////////////////////////////////////////////////////////
320/// Return pointer to Material with name.
321
323{
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// Return pointer to Material with number.
329
331{
332 TMaterial *mat;
333 if (number < 0 || number >= fMaterials->GetSize()) return 0;
334 if (fMaterialPointer) return fMaterialPointer[number];
335 TIter next(fMaterials);
336 while ((mat = (TMaterial*) next())) {
337 if (mat->GetNumber() == number) return mat;
338 }
339 return 0;
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Return pointer to node with name in the geometry tree.
344
345TNode *TGeometry::GetNode(const char *name) const
346{
347 TNode *node= (TNode*)GetListOfNodes()->First();
348 if (!node) return 0;
349 if (node->TestBit(kNotDeleted)) return node->GetNode(name);
350 return 0;
351}
352
353////////////////////////////////////////////////////////////////////////////////
354/// Return pointer to RotMatrix with name.
355
357{
359}
360
361////////////////////////////////////////////////////////////////////////////////
362/// Return pointer to RotMatrix with number.
363
365{
366 TRotMatrix *matrix;
367 if (number < 0 || number >= fMatrices->GetSize()) return 0;
368 if (fMatrixPointer) return fMatrixPointer[number];
369 TIter next(fMatrices);
370 while ((matrix = (TRotMatrix*) next())) {
371 if (matrix->GetNumber() == number) return matrix;
372 }
373 return 0;
374}
375
376////////////////////////////////////////////////////////////////////////////////
377/// Return pointer to Shape with name.
378
379TShape *TGeometry::GetShape(const char *name) const
380{
381 return (TShape*)fShapes->FindObject(name);
382}
383
384////////////////////////////////////////////////////////////////////////////////
385/// Return pointer to Shape with number.
386
388{
389 TShape *shape;
390 if (number < 0 || number >= fShapes->GetSize()) return 0;
391 if (fShapePointer) return fShapePointer[number];
392 TIter next(fShapes);
393 while ((shape = (TShape*) next())) {
394 if (shape->GetNumber() == number) return shape;
395 }
396 return 0;
397}
398
399////////////////////////////////////////////////////////////////////////////////
400/// Convert one point from local system to master reference system.
401///
402/// Note that before invoking this function, the global rotation matrix
403/// and translation vector for this node must have been computed.
404/// This is automatically done by the Paint functions.
405/// Otherwise TNode::UpdateMatrix should be called before.
406
408{
409 if (GeomLevel()) {
410 Double_t x,y,z;
411 Double_t bomb = GetBomb();
412 Double_t *matrix = &fRotMatrix[GeomLevel()][0];
413 x = bomb*fX
414 + local[0]*matrix[0]
415 + local[1]*matrix[3]
416 + local[2]*matrix[6];
417
418 y = bomb*fY
419 + local[0]*matrix[1]
420 + local[1]*matrix[4]
421 + local[2]*matrix[7];
422
423 z = bomb*fZ
424 + local[0]*matrix[2]
425 + local[1]*matrix[5]
426 + local[2]*matrix[8];
427 master[0] = x; master[1] = y; master[2] = z;
428 }
429 else
430 for (Int_t i=0;i<3;i++) master[i] = local[i];
431}
432
433////////////////////////////////////////////////////////////////////////////////
434/// Convert one point from local system to master reference system.
435///
436/// Note that before invoking this function, the global rotation matrix
437/// and translation vector for this node must have been computed.
438/// This is automatically done by the Paint functions.
439/// Otherwise TNode::UpdateMatrix should be called before.
440
442{
443 if (GeomLevel()) {
444 Float_t x,y,z;
445 Float_t bomb = GetBomb();
446
447 Double_t *matrix = &fRotMatrix[GeomLevel()][0];
448
449 x = bomb*fX
450 + local[0]*matrix[0]
451 + local[1]*matrix[3]
452 + local[2]*matrix[6];
453
454 y = bomb*fY
455 + local[0]*matrix[1]
456 + local[1]*matrix[4]
457 + local[2]*matrix[7];
458
459 z = bomb*fZ
460 + local[0]*matrix[2]
461 + local[1]*matrix[5]
462 + local[2]*matrix[8];
463
464 master[0] = x; master[1] = y; master[2] = z;
465 }
466 else
467 for (Int_t i=0;i<3;i++) master[i] = local[i];
468}
469
470////////////////////////////////////////////////////////////////////////////////
471/// List this geometry.
472
473void TGeometry::ls(Option_t *option) const
474{
475 TString opt = option;
476 opt.ToLower();
477 if (opt.Contains("m")) {
478 Printf("=================List of Materials================");
479 fMaterials->ls(option);
480 }
481 if (opt.Contains("r")) {
482 Printf("=================List of RotationMatrices================");
483 fMatrices->ls(option);
484 }
485 if (opt.Contains("s")) {
486 Printf("=================List of Shapes==========================");
487 fShapes->ls(option);
488 }
489 if (opt.Contains("n")) {
490 Printf("=================List of Nodes===========================");
491 fNodes->ls(option);
492 }
493}
494
495////////////////////////////////////////////////////////////////////////////////
496/// Convert one point from master system to local reference system.
497///
498/// Note that before invoking this function, the global rotation matrix
499/// and translation vector for this node must have been computed.
500/// This is automatically done by the Paint functions.
501/// Otherwise TNode::UpdateMatrix should be called before.
502
504{
505 if (GeomLevel()) {
506 Double_t x,y,z;
507 Double_t bomb = GetBomb();
508 Double_t *matrix = &fRotMatrix[GeomLevel()][0];
509
510 Double_t xms = master[0] - bomb*fX;
511 Double_t yms = master[1] - bomb*fY;
512 Double_t zms = master[2] - bomb*fZ;
513
514 x = xms*matrix[0] + yms*matrix[1] + zms*matrix[2];
515 y = xms*matrix[3] + yms*matrix[4] + zms*matrix[5];
516 z = xms*matrix[6] + yms*matrix[7] + zms*matrix[8];
517
518 local[0] = x; local[1] = y; local[2] = z;
519 }
520 else
521 memcpy(local,master,sizeof(Double_t)* kVectorSize);
522}
523
524////////////////////////////////////////////////////////////////////////////////
525/// Convert one point from master system to local reference system.
526///
527/// Note that before invoking this function, the global rotation matrix
528/// and translation vector for this node must have been computed.
529/// This is automatically done by the Paint functions.
530/// Otherwise TNode::UpdateMatrix should be called before.
531
533{
534 if (GeomLevel()) {
535 Float_t x,y,z;
536 Float_t bomb = GetBomb();
537
538 Double_t *matrix = &fRotMatrix[GeomLevel()][0];
539
540 Double_t xms = master[0] - bomb*fX;
541 Double_t yms = master[1] - bomb*fY;
542 Double_t zms = master[2] - bomb*fZ;
543
544 x = xms*matrix[0] + yms*matrix[1] + zms*matrix[2];
545 y = xms*matrix[3] + yms*matrix[4] + zms*matrix[5];
546 z = xms*matrix[6] + yms*matrix[7] + zms*matrix[8];
547
548 local[0] = x; local[1] = y; local[2] = z;
549 }
550 else
551 memcpy(local,master,sizeof(Float_t)* kVectorSize);
552}
553
554////////////////////////////////////////////////////////////////////////////////
555/// Add a node to the current node in this geometry.
556
557void TGeometry::Node(const char *name, const char *title, const char *shapename, Double_t x, Double_t y, Double_t z, const char *matrixname, Option_t *option)
558{
559 new TNode(name,title,shapename,x,y,z,matrixname,option);
560}
561
562////////////////////////////////////////////////////////////////////////////////
563/// Recursively remove object from a Geometry list.
564
566{
567 if (fNodes) fNodes->RecursiveRemove(obj);
568}
569
570////////////////////////////////////////////////////////////////////////////////
571/// Stream a class object.
572
573void TGeometry::Streamer(TBuffer &b)
574{
575 if (b.IsReading()) {
576 UInt_t R__s, R__c;
577 Version_t R__v = b.ReadVersion(&R__s, &R__c);
578 if (R__v > 1) {
579 b.ReadClassBuffer(TGeometry::Class(), this, R__v, R__s, R__c);
580 } else {
581 //====process old versions before automatic schema evolution
582 TNamed::Streamer(b);
583 fMaterials->Streamer(b);
584 fMatrices->Streamer(b);
585 fShapes->Streamer(b);
586 fNodes->Streamer(b);
587 b >> fBomb;
588 b.CheckByteCount(R__s, R__c, TGeometry::IsA());
589 //====end of old versions
590 }
591 // Build direct access pointers to individual materials,matrices and shapes
592 Int_t i;
593 TMaterial *onemat;
594 TRotMatrix *onematrix;
595 TShape *oneshape;
596 Int_t nmat = fMaterials->GetSize();
597 if (nmat) fMaterialPointer = new TMaterial* [nmat];
598 TIter nextmat(fMaterials);
599 i = 0;
600 while ((onemat = (TMaterial*) nextmat())) {
601 fMaterialPointer[i] = onemat;
602 i++;
603 }
604
605 Int_t nrot = fMatrices->GetSize();
606 if (nrot) fMatrixPointer = new TRotMatrix* [nrot];
607 TIter nextmatrix(fMatrices);
608 i = 0;
609 while ((onematrix = (TRotMatrix*) nextmatrix())) {
610 fMatrixPointer[i] = onematrix;
611 i++;
612 }
613
614 Int_t nsha = fShapes->GetSize();
615 if (nsha) fShapePointer = new TShape* [nsha];
616 TIter nextshape(fShapes);
617 i = 0;
618 while ((oneshape = (TShape*) nextshape())) {
619 fShapePointer[i] = oneshape;
620 i++;
621 }
622
623 gROOT->GetListOfGeometries()->Add(this);
624
626 } else {
627 b.WriteClassBuffer(TGeometry::Class(),this);
628 }
629}
630
631////////////////////////////////////////////////////////////////////////////////
632/// Update global rotation matrix/translation vector for this node
633/// this function must be called before invoking Local2Master
634
636{
637 TNode *nodes[kMAXLEVELS];
638 for (Int_t i=0;i<kVectorSize;i++) fTranslation[0][i] = 0;
639 for (Int_t i=0;i<kMatrixSize;i++) fRotMatrix[0][i] = 0;
640 fRotMatrix[0][0] = 1; fRotMatrix[0][4] = 1; fRotMatrix[0][8] = 1;
641
642 fGeomLevel = 0;
643 //build array of parent nodes
644 while (node) {
645 nodes[fGeomLevel] = node;
646 node = node->GetParent();
647 fGeomLevel++;
648 }
649 fGeomLevel--;
650 Int_t saveGeomLevel = fGeomLevel;
651 //Update matrices in the hierarchy
652 for (fGeomLevel=1;fGeomLevel<=saveGeomLevel;fGeomLevel++) {
653 node = nodes[fGeomLevel-1];
654 UpdateTempMatrix(node->GetX(),node->GetY(),node->GetZ(),node->GetMatrix());
655 }
656}
657
658////////////////////////////////////////////////////////////////////////////////
659/// Update temp matrix.
660
662{
663 Double_t *matrix = 0;
664 Bool_t isReflection = kFALSE;
665 if (rotMatrix && rotMatrix->GetType()) {
666 matrix = rotMatrix->GetMatrix();
667 isReflection = rotMatrix->IsReflection();
668 }
669 UpdateTempMatrix( x,y,z, matrix,isReflection);
670}
671
672////////////////////////////////////////////////////////////////////////////////
673/// Update temp matrix.
674
676{
677 Int_t i=GeomLevel();
678 if (i) {
679 if(matrix) {
680 UpdateTempMatrix(&(fTranslation[i-1][0]),&fRotMatrix[i-1][0]
681 ,x,y,z,matrix
682 ,&fTranslation[i][0],&fRotMatrix[i][0]);
683 fX = fTranslation[i][0];
684 fY = fTranslation[i][1];
685 fZ = fTranslation[i][2];
686 fIsReflection[i] = fIsReflection[i-1] ^ isReflection;
687 } else {
688 fX = fTranslation[i][0] = fTranslation[i-1][0] + x;
689 fY = fTranslation[i][1] = fTranslation[i-1][1] + y;
690 fZ = fTranslation[i][2] = fTranslation[i-1][2] + z;
691 }
692 } else {
693 fX=fY=fZ=0;
695 for (i=0;i<kVectorSize;i++) fTranslation[0][i] = 0;
696 for (i=0;i<kMatrixSize;i++) fRotMatrix[0][i] = 0;
697 fRotMatrix[0][0] = 1; fRotMatrix[0][4] = 1; fRotMatrix[0][8] = 1;
698 }
699}
700
701////////////////////////////////////////////////////////////////////////////////
702/// Compute new translation vector and global matrix.
703///
704/// - dx old translation vector
705/// - rmat old global matrix
706/// - x,y,z offset of new local system with respect to mother
707/// - dxnew new translation vector
708/// - rmatnew new global rotation matrix
709
711 , Double_t x, Double_t y, Double_t z, Double_t *matrix
712 , Double_t *dxnew, Double_t *rmatnew)
713{
714 dxnew[0] = dx[0] + x*rmat[0] + y*rmat[3] + z*rmat[6];
715 dxnew[1] = dx[1] + x*rmat[1] + y*rmat[4] + z*rmat[7];
716 dxnew[2] = dx[2] + x*rmat[2] + y*rmat[5] + z*rmat[8];
717
718 rmatnew[0] = rmat[0]*matrix[0] + rmat[3]*matrix[1] + rmat[6]*matrix[2];
719 rmatnew[1] = rmat[1]*matrix[0] + rmat[4]*matrix[1] + rmat[7]*matrix[2];
720 rmatnew[2] = rmat[2]*matrix[0] + rmat[5]*matrix[1] + rmat[8]*matrix[2];
721 rmatnew[3] = rmat[0]*matrix[3] + rmat[3]*matrix[4] + rmat[6]*matrix[5];
722 rmatnew[4] = rmat[1]*matrix[3] + rmat[4]*matrix[4] + rmat[7]*matrix[5];
723 rmatnew[5] = rmat[2]*matrix[3] + rmat[5]*matrix[4] + rmat[8]*matrix[5];
724 rmatnew[6] = rmat[0]*matrix[6] + rmat[3]*matrix[7] + rmat[6]*matrix[8];
725 rmatnew[7] = rmat[1]*matrix[6] + rmat[4]*matrix[7] + rmat[7]*matrix[8];
726 rmatnew[8] = rmat[2]*matrix[6] + rmat[5]*matrix[7] + rmat[8]*matrix[8];
727}
#define b(i)
Definition RSha256.hxx:100
int Int_t
Definition RtypesCore.h:45
short Version_t
Definition RtypesCore.h:65
unsigned int UInt_t
Definition RtypesCore.h:46
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
TGeometry * gGeometry
Definition TGeometry.cxx:21
const Int_t kVectorSize
Definition TGeometry.h:28
const Int_t kMAXLEVELS
Definition TGeometry.h:27
const Int_t kMatrixSize
Definition TGeometry.h:29
R__EXTERN TGeometry * gGeometry
Definition TGeometry.h:158
#define gROOT
Definition TROOT.h:406
void Printf(const char *fmt,...)
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual void ls(Option_t *option="") const
List (ls) all objects in this collection.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
TGeometry description.
Definition TGeometry.h:39
TShape * GetShapeByNumber(Int_t number) const
Return pointer to Shape with number.
Double_t fZ
The global translation of the current node.
Definition TGeometry.h:55
Double_t fRotMatrix[kMAXLEVELS][kMatrixSize]
Definition TGeometry.h:57
Float_t GetBomb() const
Definition TGeometry.h:73
TMaterial * GetMaterialByNumber(Int_t number) const
Return pointer to Material with number.
TGeometry & operator=(const TGeometry &)
assignment operator
virtual void Local2Master(Double_t *local, Double_t *master)
Convert one point from local system to master reference system.
Int_t GeomLevel() const
Definition TGeometry.h:74
virtual void UpdateTempMatrix(Double_t x=0, Double_t y=0, Double_t z=0, TRotMatrix *matrix=0)
Update temp matrix.
TRotMatrix * fMatrix
Definition TGeometry.h:46
TMaterial ** fMaterialPointer
Pointer to current node.
Definition TGeometry.h:48
TMaterial * GetMaterial(const char *name) const
Return pointer to Material with name.
virtual TObject * FindObject(const char *name) const
Search object identified by name in the geometry tree.
virtual void Master2Local(Double_t *master, Double_t *local)
Convert one point from master system to local reference system.
TGeometry()
Geometry default constructor.
Definition TGeometry.cxx:96
THashList * fMaterials
Definition TGeometry.h:42
TRotMatrix ** fMatrixPointer
Pointers to materials.
Definition TGeometry.h:49
TShape ** fShapePointer
Pointers to rotation matrices.
Definition TGeometry.h:50
THashList * GetListOfShapes() const
Definition TGeometry.h:75
TRotMatrix * GetRotMatrixByNumber(Int_t number) const
Return pointer to RotMatrix with number.
TList * fNodes
Definition TGeometry.h:45
virtual ~TGeometry()
Geometry default destructor.
TNode * fCurrentNode
Pointers to current rotation matrices.
Definition TGeometry.h:47
Int_t fGeomLevel
Definition TGeometry.h:52
Bool_t fIsReflection[kMAXLEVELS]
Definition TGeometry.h:58
THashList * fMatrices
Definition TGeometry.h:43
Double_t fY
Definition TGeometry.h:54
virtual void RecursiveRemove(TObject *obj)
Recursively remove object from a Geometry list.
THashList * fShapes
Definition TGeometry.h:44
THashList * GetListOfMatrices() const
Definition TGeometry.h:78
THashList * GetListOfMaterials() const
Definition TGeometry.h:77
Float_t fBomb
Pointers to shapes.
Definition TGeometry.h:51
static TObjArray * Get(const char *name)
Static function called by TROOT to search name in the geometry.
TRotMatrix * GetRotMatrix(const char *name) const
Return pointer to RotMatrix with name.
virtual void Browse(TBrowser *b)
Browse.
virtual void Draw(Option_t *option="")
Draw this Geometry.
virtual void Node(const char *name, const char *title, const char *shapename, Double_t x=0, Double_t y=0, Double_t z=0, const char *matrixname="", Option_t *option="")
Add a node to the current node in this geometry.
TShape * GetShape(const char *name) const
Return pointer to Shape with name.
virtual void UpdateMatrix(TNode *node)
Update global rotation matrix/translation vector for this node this function must be called before in...
virtual void cd(const char *path=0)
Change Current Geometry to this.
TList * GetListOfNodes() const
Definition TGeometry.h:76
Double_t fX
Definition TGeometry.h:53
TNode * GetNode(const char *name) const
Return pointer to node with name in the geometry tree.
virtual void ls(Option_t *option="rsn2") const
List this geometry.
Double_t fTranslation[kMAXLEVELS][kVectorSize]
Definition TGeometry.h:56
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition THashList.h:34
TObject * FindObject(const char *name) const
Find object using its name.
void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
A doubly linked list.
Definition TList.h:44
virtual void RecursiveRemove(TObject *obj)
Remove object from this collection and recursively remove the object from all other objects (and coll...
Definition TList.cxx:764
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:470
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:659
Manages a detector material.
Definition TMaterial.h:28
virtual Int_t GetNumber() const
Definition TMaterial.h:42
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition TNamed.cxx:51
TNode description.
Definition TNode.h:33
virtual void Draw(Option_t *option="")
Draw Referenced node with current parameters.
Definition TNode.cxx:324
virtual TNode * GetParent() const
Definition TNode.h:70
virtual Double_t GetY() const
Definition TNode.h:74
virtual Double_t GetX() const
Definition TNode.h:73
virtual TRotMatrix * GetMatrix() const
Definition TNode.h:66
virtual Double_t GetZ() const
Definition TNode.h:75
virtual TNode * GetNode(const char *name) const
Return pointer to node with name in the node tree.
Definition TNode.cxx:379
An array of TObjects.
Definition TObjArray.h:37
TObject * At(Int_t idx) const
Definition TObjArray.h:166
Mother of all ROOT objects.
Definition TObject.h:37
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:187
@ kNotDeleted
object has not been deleted
Definition TObject.h:78
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:893
Manages a detector rotation matrix.
Definition TRotMatrix.h:28
virtual Int_t GetType() const
Definition TRotMatrix.h:56
virtual Int_t GetNumber() const
Definition TRotMatrix.h:55
virtual Double_t * GetMatrix()
Definition TRotMatrix.h:54
virtual Bool_t IsReflection() const
Definition TRotMatrix.h:61
This is the base class for all geometry shapes.
Definition TShape.h:35
virtual Int_t GetNumber() const
Definition TShape.h:57
Basic string class.
Definition TString.h:136
void ToLower()
Change string to lower-case.
Definition TString.cxx:1145
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17