Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGeoManager.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 TGeoManager
13\ingroup Geometry_classes
14
15The manager class for any TGeo geometry. Provides user
16interface for geometry creation, navigation, state querying,
17visualization, IO, geometry checking and other utilities.
18
19## General architecture
20
21 The ROOT geometry package is a tool designed for building, browsing,
22tracking and visualizing a detector geometry. The code is independent from
23other external MC for simulation, therefore it does not contain any
24constraints related to physics. However, the package defines a number of
25hooks for tracking, such as media, materials, magnetic field or track state flags,
26in order to allow interfacing to tracking MC's. The final goal is to be
27able to use the same geometry for several purposes, such as tracking,
28reconstruction or visualization, taking advantage of the ROOT features
29related to bookkeeping, I/O, histogramming, browsing and GUI's.
30
31 The geometrical modeler is the most important component of the package and
32it provides answers to the basic questions like "Where am I ?" or "How far
33from the next boundary ?", but also to more complex ones like "How far from
34the closest surface ?" or "Which is the next crossing along a helix ?".
35
36 The architecture of the modeler is a combination between a GEANT-like
37containment scheme and a normal CSG binary tree at the level of shapes. An
38important common feature of all detector geometry descriptions is the
39mother-daughter concept. This is the most natural approach when tracking
40is concerned and imposes a set of constraints to the way geometry is defined.
41Constructive solid geometry composition is used only in order to create more
42complex shapes from an existing set of primitives through boolean operations.
43This feature is not implemented yet but in future full definition of boolean
44expressions will be supported.
45
46 Practically every geometry defined in GEANT style can be mapped by the modeler.
47The basic components used for building the logical hierarchy of the geometry
48are called "volumes" and "nodes". Volumes (sometimes called "solids") are fully
49defined geometrical objects having a given shape and medium and possibly
50containing a list of nodes. Nodes represent just positioned instances of volumes
51inside a container volume and they are not directly defined by user. They are
52automatically created as a result of adding one volume inside other or dividing
53a volume. The geometrical transformation hold by nodes is always defined with
54respect to their mother (relative positioning). Reflection matrices are allowed.
55All volumes have to be fully aware of their containees when the geometry is
56closed. They will build additional structures (voxels) in order to fasten-up
57the search algorithms. Finally, nodes can be regarded as bidirectional links
58between containers and containees objects.
59
60 The structure defined in this way is a graph structure since volumes are
61replicable (same volume can become daughter node of several other volumes),
62every volume becoming a branch in this graph. Any volume in the logical graph
63can become the actual top volume at run time (see TGeoManager::SetTopVolume()).
64All functionalities of the modeler will behave in this case as if only the
65corresponding branch starting from this volume is the registered geometry.
66
67\image html geom_graf.jpg
68
69 A given volume can be positioned several times in the geometry. A volume
70can be divided according default or user-defined patterns, creating automatically
71the list of division nodes inside. The elementary volumes created during the
72dividing process follow the same scheme as usual volumes, therefore it is possible
73to position further geometrical structures inside or to divide them further more
74(see TGeoVolume::Divide()).
75
76 The primitive shapes supported by the package are basically the GEANT3
77shapes (see class TGeoShape), arbitrary wedges with eight vertices on two parallel
78planes. All basic primitives inherits from class TGeoBBox since the bounding box
79of a solid is essential for the tracking algorithms. They also implement the
80virtual methods defined in the virtual class TGeoShape (point and segment
81classification). User-defined primitives can be directly plugged into the modeler
82provided that they override these methods. Composite shapes will be soon supported
83by the modeler. In order to build a TGeoCompositeShape, one will have to define
84first the primitive components. The object that handle boolean
85operations among components is called TGeoBoolCombinator and it has to be
86constructed providing a string boolean expression between the components names.
87
88
89## Example for building a simple geometry
90
91Begin_Macro(source)
92../../../tutorials/geom/rootgeom.C
93End_Macro
94
95## TGeoManager - the manager class for the geometry package.
96
97 TGeoManager class is embedding all the API needed for building and tracking
98a geometry. It defines a global pointer (gGeoManager) in order to be fully
99accessible from external code. The mechanism of handling multiple geometries
100at the same time will be soon implemented.
101
102 TGeoManager is the owner of all geometry objects defined in a session,
103therefore users must not try to control their deletion. It contains lists of
104media, materials, transformations, shapes and volumes. Logical nodes (positioned
105volumes) are created and destroyed by the TGeoVolume class. Physical
106nodes and their global transformations are subjected to a caching mechanism
107due to the sometimes very large memory requirements of logical graph expansion.
108The caching mechanism is triggered by the total number of physical instances
109of volumes and the cache manager is a client of TGeoManager. The manager class
110also controls the painter client. This is linked with ROOT graphical libraries
111loaded on demand in order to control visualization actions.
112
113## Rules for building a valid geometry
114
115 A given geometry can be built in various ways, but there are mandatory steps
116that have to be followed in order to be validated by the modeler. There are
117general rules : volumes needs media and shapes in order to be created,
118both container and containee volumes must be created before linking them together,
119and the relative transformation matrix must be provided. All branches must
120have an upper link point otherwise they will not be considered as part of the
121geometry. Visibility or tracking properties of volumes can be provided both
122at build time or after geometry is closed, but global visualization settings
123(see TGeoPainter class) should not be provided at build time, otherwise the
124drawing package will be loaded. There is also a list of specific rules :
125positioned daughters should not extrude their mother or intersect with sisters
126unless this is specified (see TGeoVolume::AddNodeOverlap()), the top volume
127(containing all geometry tree) must be specified before closing the geometry
128and must not be positioned - it represents the global reference frame. After
129building the full geometry tree, the geometry must be closed
130(see TGeoManager::CloseGeometry()). Voxelization can be redone per volume after
131this process.
132
133
134 Below is the general scheme of the manager class.
135
136\image html geom_mgr.jpg
137
138## An interactive session
139
140 Provided that a geometry was successfully built and closed (for instance the
141previous example $ROOTSYS/tutorials/geom/rootgeom.C ), the manager class will register
142itself to ROOT and the logical/physical structures will become immediately browsable.
143The ROOT browser will display starting from the geometry folder : the list of
144transformations and media, the top volume and the top logical node. These last
145two can be fully expanded, any intermediate volume/node in the browser being subject
146of direct access context menu operations (right mouse button click). All user
147utilities of classes TGeoManager, TGeoVolume and TGeoNode can be called via the
148context menu.
149
150\image html geom_browser.jpg
151
152### Drawing the geometry
153
154 Any logical volume can be drawn via TGeoVolume::Draw() member function.
155This can be directly accessed from the context menu of the volume object
156directly from the browser.
157 There are several drawing options that can be set with
158TGeoManager::SetVisOption(Int_t opt) method :
159
160#### opt=0
161 only the content of the volume is drawn, N levels down (default N=3).
162 This is the default behavior. The number of levels to be drawn can be changed
163 via TGeoManager::SetVisLevel(Int_t level) method.
164
165\image html geom_frame0.jpg
166
167#### opt=1
168 the final leaves (e.g. daughters with no containment) of the branch
169 starting from volume are drawn down to the current number of levels.
170 WARNING : This mode is memory consuming
171 depending of the size of geometry, so drawing from top level within this mode
172 should be handled with care for expensive geometries. In future there will be
173 a limitation on the maximum number of nodes to be visualized.
174
175\image html geom_frame1.jpg
176
177#### opt=2
178 only the clicked volume is visualized. This is automatically set by
179 TGeoVolume::DrawOnly() method
180
181#### opt=3 - only a given path is visualized. This is automatically set by
182 TGeoVolume::DrawPath(const char *path) method
183
184 The current view can be exploded in cartesian, cylindrical or spherical
185coordinates :
186 TGeoManager::SetExplodedView(Int_t opt). Options may be :
187- 0 - default (no bombing)
188- 1 - cartesian coordinates. The bomb factor on each axis can be set with
189 TGeoManager::SetBombX(Double_t bomb) and corresponding Y and Z.
190- 2 - bomb in cylindrical coordinates. Only the bomb factors on Z and R
191 are considered
192 \image html geom_frameexp.jpg
193
194- 3 - bomb in radial spherical coordinate : TGeoManager::SetBombR()
195
196Volumes themselves support different visualization settings :
197 - TGeoVolume::SetVisibility() : set volume visibility.
198 - TGeoVolume::VisibleDaughters() : set daughters visibility.
199All these actions automatically updates the current view if any.
200
201### Checking the geometry
202
203 Several checking methods are accessible from the volume context menu. They
204generally apply only to the visible parts of the drawn geometry in order to
205ease geometry checking, and their implementation is in the TGeoChecker class
206from the painting package.
207
208#### Checking a given point.
209 Can be called from TGeoManager::CheckPoint(Double_t x, Double_t y, Double_t z).
210This method is drawing the daughters of the volume containing the point one
211level down, printing the path to the deepest physical node holding this point.
212It also computes the closest distance to any boundary. The point will be drawn
213in red, as well as a sphere having this closest distance as radius. In case a
214non-zero distance is given by the user as fifth argument of CheckPoint, this
215distance will be used as radius of the safety sphere.
216
217\image html geom_checkpoint.jpg
218
219#### Shooting random points.
220 Can be called from TGeoVolume::RandomPoints() (context menu function) and
221it will draw this volume with current visualization settings. Random points
222are generated in the bounding box of the top drawn volume. The points are
223classified and drawn with the color of their deepest container. Only points
224in visible nodes will be drawn.
225
226\image html geom_random1.jpg
227
228
229#### Raytracing.
230 Can be called from TGeoVolume::RandomRays() (context menu of volumes) and
231will shoot rays from a given point in the local reference frame with random
232directions. The intersections with displayed nodes will appear as segments
233having the color of the touched node. Drawn geometry will be then made invisible
234in order to enhance rays.
235
236\image html geom_random2.jpg
237*/
238
239#include <cstdlib>
240#include <iostream>
241#include <fstream>
242
243#include "TROOT.h"
244#include "TGeoManager.h"
245#include "TStyle.h"
246#include "TVirtualPad.h"
247#include "TBrowser.h"
248#include "TFile.h"
249#include "TKey.h"
250#include "THashList.h"
251#include "TClass.h"
252#include "ThreadLocalStorage.h"
253#include "TBufferText.h"
254
255#include "TGeoVoxelFinder.h"
256#include "TGeoElement.h"
257#include "TGeoMaterial.h"
258#include "TGeoMedium.h"
259#include "TGeoMatrix.h"
260#include "TGeoNode.h"
261#include "TGeoPhysicalNode.h"
262#include "TGeoPara.h"
263#include "TGeoParaboloid.h"
264#include "TGeoTube.h"
265#include "TGeoEltu.h"
266#include "TGeoHype.h"
267#include "TGeoCone.h"
268#include "TGeoSphere.h"
269#include "TGeoArb8.h"
270#include "TGeoPgon.h"
271#include "TGeoTrd1.h"
272#include "TGeoTrd2.h"
273#include "TGeoTorus.h"
274#include "TGeoXtru.h"
275#include "TGeoCompositeShape.h"
276#include "TGeoBoolNode.h"
277#include "TGeoBuilder.h"
278#include "TVirtualGeoPainter.h"
279#include "TPluginManager.h"
280#include "TVirtualGeoTrack.h"
281#include "TQObject.h"
282#include "TMath.h"
283#include "TEnv.h"
284#include "TGeoParallelWorld.h"
285#include "TGeoRegion.h"
286#include "TGDMLMatrix.h"
287#include "TGeoOpticalSurface.h"
288
289// statics and globals
290
292
294
295std::mutex TGeoManager::fgMutex;
307
308////////////////////////////////////////////////////////////////////////////////
309/// Default constructor.
310
312{
313 if (!fgThreadId)
317 fTmin = 0.;
318 fTmax = 999.;
319 fPhiCut = kFALSE;
320 fPhimin = 0;
321 fPhimax = 360;
326 fClosed = kFALSE;
328 fBits = nullptr;
329 fCurrentNavigator = nullptr;
330 fMaterials = nullptr;
331 fHashPNE = nullptr;
332 fArrayPNE = nullptr;
333 fMatrices = nullptr;
334 fNodes = nullptr;
335 fOverlaps = nullptr;
336 fRegions = nullptr;
337 fNNodes = 0;
338 fMaxVisNodes = 10000;
339 fVolumes = nullptr;
340 fPhysicalNodes = nullptr;
341 fShapes = nullptr;
342 fGVolumes = nullptr;
343 fGShapes = nullptr;
344 fTracks = nullptr;
345 fMedia = nullptr;
346 fNtracks = 0;
347 fNpdg = 0;
348 fPdgNames = nullptr;
349 fGDMLMatrices = nullptr;
350 fOpticalSurfaces = nullptr;
351 fSkinSurfaces = nullptr;
352 fBorderSurfaces = nullptr;
353 memset(fPdgId, 0, 1024 * sizeof(Int_t));
354 // TObjArray *fNavigators; //! list of navigators
355 fCurrentTrack = nullptr;
356 fCurrentVolume = nullptr;
357 fTopVolume = nullptr;
358 fTopNode = nullptr;
359 fMasterVolume = nullptr;
360 fPainter = nullptr;
363 fVisDensity = 0.;
364 fVisLevel = 3;
365 fVisOption = 1;
366 fExplodedView = 0;
367 fNsegments = 20;
368 fNLevel = 0;
369 fUniqueVolumes = nullptr;
370 fClippingShape = nullptr;
373 fGLMatrix = nullptr;
374 fPaintVolume = nullptr;
375 fUserPaintVolume = nullptr;
376 fElementTable = nullptr;
377 fHashVolumes = nullptr;
378 fHashGVolumes = nullptr;
379 fSizePNEId = 0;
380 fNPNEId = 0;
381 fKeyPNEId = nullptr;
382 fValuePNEId = nullptr;
384 fRaytraceMode = 0;
385 fMaxThreads = 0;
387 fParallelWorld = nullptr;
389 } else {
390 Init();
392 gGeoIdentity = new TGeoIdentity("Identity");
394 }
395}
396
397////////////////////////////////////////////////////////////////////////////////
398/// Constructor.
399
400TGeoManager::TGeoManager(const char *name, const char *title) : TNamed(name, title)
401{
402 if (!gROOT->GetListOfGeometries()->FindObject(this))
403 gROOT->GetListOfGeometries()->Add(this);
404 if (!gROOT->GetListOfBrowsables()->FindObject(this))
405 gROOT->GetListOfBrowsables()->Add(this);
406 Init();
407 gGeoIdentity = new TGeoIdentity("Identity");
409 if (fgVerboseLevel > 0)
410 Info("TGeoManager", "Geometry %s, %s created", GetName(), GetTitle());
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// Initialize manager class.
415
417{
418 if (gGeoManager) {
419 Warning("Init", "Deleting previous geometry: %s/%s", gGeoManager->GetName(), gGeoManager->GetTitle());
420 delete gGeoManager;
421 if (fgLock)
422 Fatal("Init", "New geometry created while the old one locked !!!");
423 }
424
425 gGeoManager = this;
426 if (!fgThreadId)
429 fTmin = 0.;
430 fTmax = 999.;
431 fPhiCut = kFALSE;
432 fPhimin = 0;
433 fPhimax = 360;
438 fClosed = kFALSE;
440 fBits = new UChar_t[50000]; // max 25000 nodes per volume
441 fCurrentNavigator = nullptr;
442 fHashPNE = new THashList(256, 3);
443 fArrayPNE = nullptr;
444 fMaterials = new THashList(200, 3);
445 fMatrices = new TObjArray(256);
446 fNodes = new TObjArray(30);
447 fOverlaps = new TObjArray(256);
448 fRegions = new TObjArray(256);
449 fNNodes = 0;
450 fMaxVisNodes = 10000;
451 fVolumes = new TObjArray(256);
452 fPhysicalNodes = new TObjArray(256);
453 fShapes = new TObjArray(256);
454 fGVolumes = new TObjArray(256);
455 fGShapes = new TObjArray(256);
456 fTracks = new TObjArray(256);
457 fMedia = new THashList(200, 3);
458 fNtracks = 0;
459 fNpdg = 0;
460 fPdgNames = nullptr;
461 fGDMLMatrices = new TObjArray();
463 fSkinSurfaces = new TObjArray();
465 memset(fPdgId, 0, 1024 * sizeof(Int_t));
466 fCurrentTrack = nullptr;
467 fCurrentVolume = nullptr;
468 fTopVolume = nullptr;
469 fTopNode = nullptr;
470 fMasterVolume = nullptr;
471 fPainter = nullptr;
474 fVisDensity = 0.;
475 fVisLevel = 3;
476 fVisOption = 1;
477 fExplodedView = 0;
478 fNsegments = 20;
479 fNLevel = 0;
480 fUniqueVolumes = new TObjArray(256);
481 fClippingShape = nullptr;
484 fGLMatrix = new TGeoHMatrix();
485 fPaintVolume = nullptr;
486 fUserPaintVolume = nullptr;
487 fElementTable = nullptr;
488 fHashVolumes = nullptr;
489 fHashGVolumes = nullptr;
490 fSizePNEId = 0;
491 fNPNEId = 0;
492 fKeyPNEId = nullptr;
493 fValuePNEId = nullptr;
495 fRaytraceMode = 0;
496 fMaxThreads = 0;
498 fParallelWorld = nullptr;
500}
501
502////////////////////////////////////////////////////////////////////////////////
503/// Destructor
504
506{
507 if (gGeoManager != this)
508 gGeoManager = this;
510
511 if (gROOT->GetListOfFiles()) { // in case this function is called from TROOT destructor
512 gROOT->GetListOfGeometries()->Remove(this);
513 gROOT->GetListOfBrowsables()->Remove(this);
514 }
515 // TSeqCollection *brlist = gROOT->GetListOfBrowsers();
516 // TIter next(brlist);
517 // TBrowser *browser = 0;
518 // while ((browser=(TBrowser*)next())) browser->RecursiveRemove(this);
521 delete TGeoBuilder::Instance(this);
522 if (fBits)
523 delete[] fBits;
526 if (fOverlaps) {
527 fOverlaps->Delete();
529 }
530 if (fRegions) {
531 fRegions->Delete();
533 }
534 if (fMaterials) {
537 }
539 if (fMedia) {
540 fMedia->Delete();
542 }
543 if (fHashVolumes) {
544 fHashVolumes->Clear("nodelete");
546 }
547 if (fHashGVolumes) {
548 fHashGVolumes->Clear("nodelete");
550 }
551 if (fHashPNE) {
552 fHashPNE->Delete();
554 }
555 if (fArrayPNE) {
556 delete fArrayPNE;
557 }
558 if (fVolumes) {
559 fVolumes->Delete();
561 }
562 if (fShapes) {
563 fShapes->Delete();
565 }
566 if (fPhysicalNodes) {
569 }
570 if (fMatrices) {
571 fMatrices->Delete();
573 }
574 if (fTracks) {
575 fTracks->Delete();
577 }
579 if (fPdgNames) {
580 fPdgNames->Delete();
582 }
583 if (fGDMLMatrices) {
586 }
587 if (fOpticalSurfaces) {
590 }
591 if (fSkinSurfaces) {
594 }
595 if (fBorderSurfaces) {
598 }
600 CleanGarbage();
603 if (fSizePNEId) {
604 delete[] fKeyPNEId;
605 delete[] fValuePNEId;
606 }
607 delete fParallelWorld;
609 gGeoIdentity = nullptr;
610 gGeoManager = nullptr;
611}
612
613////////////////////////////////////////////////////////////////////////////////
614/// Add a material to the list. Returns index of the material in list.
615
617{
618 return TGeoBuilder::Instance(this)->AddMaterial((TGeoMaterial *)material);
619}
620
621////////////////////////////////////////////////////////////////////////////////
622/// Add an illegal overlap/extrusion to the list.
623
625{
627 fOverlaps->Add((TObject *)ovlp);
628 return size;
629}
630
631////////////////////////////////////////////////////////////////////////////////
632/// Add a new region of volumes.
634{
636 fRegions->Add(region);
637 return size;
638}
639
640////////////////////////////////////////////////////////////////////////////////
641/// Add a user-defined property. Returns true if added, false if existing.
642
644{
645 auto pos = fProperties.insert(ConstPropMap_t::value_type(property, value));
646 if (!pos.second) {
647 Warning("AddProperty", "Property \"%s\" already exists with value %g", property, (pos.first)->second);
648 return false;
649 }
650 return true;
651}
652
653////////////////////////////////////////////////////////////////////////////////
654/// Get a user-defined property
655
657{
658 auto pos = fProperties.find(property);
659 if (pos == fProperties.end()) {
660 if (error)
661 *error = kTRUE;
662 return 0.;
663 }
664 if (error)
665 *error = kFALSE;
666 return pos->second;
667}
668
669////////////////////////////////////////////////////////////////////////////////
670/// Get a user-defined property from a given index
671
673{
674 // This is a quite inefficient way to access map elements, but needed for the GDML writer to
675 if (i >= fProperties.size()) {
676 if (error)
677 *error = kTRUE;
678 return 0.;
679 }
680 size_t pos = 0;
681 auto it = fProperties.begin();
682 while (pos < i) {
683 ++it;
684 ++pos;
685 }
686 if (error)
687 *error = kFALSE;
688 name = (*it).first;
689 return (*it).second;
690}
691
692////////////////////////////////////////////////////////////////////////////////
693/// Add a matrix to the list. Returns index of the matrix in list.
694
696{
697 return TGeoBuilder::Instance(this)->AddTransformation((TGeoMatrix *)matrix);
698}
699
700////////////////////////////////////////////////////////////////////////////////
701/// Add a shape to the list. Returns index of the shape in list.
702
704{
705 return TGeoBuilder::Instance(this)->AddShape((TGeoShape *)shape);
706}
707
708////////////////////////////////////////////////////////////////////////////////
709/// Add a track to the list of tracks. Use this for primaries only. For secondaries,
710/// add them to the parent track. The method create objects that are registered
711/// to the analysis manager but have to be cleaned-up by the user via ClearTracks().
712
714{
716 fTracks->AddAtAndExpand(GetGeomPainter()->AddTrack(id, pdgcode, particle), fNtracks++);
717 return index;
718}
719
720////////////////////////////////////////////////////////////////////////////////
721/// Add a track to the list of tracks
722
724{
727 return index;
728}
729
730////////////////////////////////////////////////////////////////////////////////
731/// Makes a primary track but do not attach it to the list of tracks. The track
732/// can be attached as daughter to another one with TVirtualGeoTrack::AddTrack
733
735{
736 TVirtualGeoTrack *track = GetGeomPainter()->AddTrack(id, pdgcode, particle);
737 return track;
738}
739
740////////////////////////////////////////////////////////////////////////////////
741/// Add a volume to the list. Returns index of the volume in list.
742
744{
745 if (!volume) {
746 Error("AddVolume", "invalid volume");
747 return -1;
748 }
750 if (!uid)
751 uid++;
752 if (!fCurrentVolume) {
753 fCurrentVolume = volume;
754 fUniqueVolumes->AddAtAndExpand(volume, uid);
755 } else {
756 if (!strcmp(volume->GetName(), fCurrentVolume->GetName())) {
757 uid = fCurrentVolume->GetNumber();
758 } else {
759 fCurrentVolume = volume;
760 Int_t olduid = GetUID(volume->GetName());
761 if (olduid < 0) {
762 fUniqueVolumes->AddAtAndExpand(volume, uid);
763 } else {
764 uid = olduid;
765 }
766 }
767 }
768 volume->SetNumber(uid);
769 if (!fHashVolumes) {
770 fHashVolumes = new THashList(256);
771 fHashGVolumes = new THashList(256);
772 }
773 TObjArray *list = fVolumes;
774 if (!volume->GetShape() || volume->IsRunTime() || volume->IsVolumeMulti()) {
775 list = fGVolumes;
776 fHashGVolumes->Add(volume);
777 } else {
778 fHashVolumes->Add(volume);
779 }
780 Int_t index = list->GetEntriesFast();
781 list->AddAtAndExpand(volume, index);
782 return uid;
783}
784
785////////////////////////////////////////////////////////////////////////////////
786/// Add a navigator in the list of navigators. If it is the first one make it
787/// current navigator.
788
790{
791 if (fMultiThread) {
793 fgMutex.lock();
794 }
795 std::thread::id threadId = std::this_thread::get_id();
796 NavigatorsMap_t::const_iterator it = fNavigators.find(threadId);
797 TGeoNavigatorArray *array = nullptr;
798 if (it != fNavigators.end())
799 array = it->second;
800 else {
801 array = new TGeoNavigatorArray(this);
802 fNavigators.insert(NavigatorsMap_t::value_type(threadId, array));
803 }
804 TGeoNavigator *nav = array->AddNavigator();
805 if (fClosed)
806 nav->GetCache()->BuildInfoBranch();
807 if (fMultiThread)
808 fgMutex.unlock();
809 return nav;
810}
811
812////////////////////////////////////////////////////////////////////////////////
813/// Returns current navigator for the calling thread.
814
816{
817 TTHREAD_TLS(TGeoNavigator *) tnav = nullptr;
818 if (!fMultiThread)
819 return fCurrentNavigator;
820 TGeoNavigator *nav = tnav; // TTHREAD_TLS_GET(TGeoNavigator*,tnav);
821 if (nav)
822 return nav;
823 std::thread::id threadId = std::this_thread::get_id();
824 NavigatorsMap_t::const_iterator it = fNavigators.find(threadId);
825 if (it == fNavigators.end())
826 return nullptr;
827 TGeoNavigatorArray *array = it->second;
828 nav = array->GetCurrentNavigator();
829 tnav = nav; // TTHREAD_TLS_SET(TGeoNavigator*,tnav,nav);
830 return nav;
831}
832
833////////////////////////////////////////////////////////////////////////////////
834/// Get list of navigators for the calling thread.
835
837{
838 std::thread::id threadId = std::this_thread::get_id();
839 NavigatorsMap_t::const_iterator it = fNavigators.find(threadId);
840 if (it == fNavigators.end())
841 return nullptr;
842 TGeoNavigatorArray *array = it->second;
843 return array;
844}
845
846////////////////////////////////////////////////////////////////////////////////
847/// Switch to another existing navigator for the calling thread.
848
850{
851 std::thread::id threadId = std::this_thread::get_id();
852 NavigatorsMap_t::const_iterator it = fNavigators.find(threadId);
853 if (it == fNavigators.end()) {
854 Error("SetCurrentNavigator", "No navigator defined for this thread\n");
855 std::cout << " thread id: " << threadId << std::endl;
856 return kFALSE;
857 }
858 TGeoNavigatorArray *array = it->second;
860 if (!nav) {
861 Error("SetCurrentNavigator", "Navigator %d not existing for this thread\n", index);
862 std::cout << " thread id: " << threadId << std::endl;
863 return kFALSE;
864 }
865 if (!fMultiThread)
866 fCurrentNavigator = nav;
867 return kTRUE;
868}
869
870////////////////////////////////////////////////////////////////////////////////
871/// Set the lock for navigators.
872
874{
875 fgLockNavigators = flag;
876}
877
878////////////////////////////////////////////////////////////////////////////////
879/// Clear all navigators.
880
882{
883 if (fMultiThread)
884 fgMutex.lock();
885 TGeoNavigatorArray *arr = nullptr;
886 for (NavigatorsMap_t::iterator it = fNavigators.begin(); it != fNavigators.end(); ++it) {
887 arr = (*it).second;
888 if (arr)
889 delete arr;
890 }
891 fNavigators.clear();
892 if (fMultiThread)
893 fgMutex.unlock();
894}
895
896////////////////////////////////////////////////////////////////////////////////
897/// Clear a single navigator.
898
900{
901 if (fMultiThread)
902 fgMutex.lock();
903 for (NavigatorsMap_t::iterator it = fNavigators.begin(); it != fNavigators.end(); ++it) {
904 TGeoNavigatorArray *arr = (*it).second;
905 if (arr) {
906 if ((TGeoNavigator *)arr->Remove((TObject *)nav)) {
907 delete nav;
908 if (!arr->GetEntries())
909 fNavigators.erase(it);
910 if (fMultiThread)
911 fgMutex.unlock();
912 return;
913 }
914 }
915 }
916 Error("Remove navigator", "Navigator %p not found", nav);
917 if (fMultiThread)
918 fgMutex.unlock();
919}
920
921////////////////////////////////////////////////////////////////////////////////
922/// Set maximum number of threads for navigation.
923
925{
926 if (!fClosed) {
927 Error("SetMaxThreads", "Cannot set maximum number of threads before closing the geometry");
928 return;
929 }
930 if (!fMultiThread) {
932 std::thread::id threadId = std::this_thread::get_id();
933 NavigatorsMap_t::const_iterator it = fNavigators.find(threadId);
934 if (it != fNavigators.end()) {
935 TGeoNavigatorArray *array = it->second;
936 fNavigators.erase(it);
937 fNavigators.insert(NavigatorsMap_t::value_type(threadId, array));
938 }
939 }
940 if (fMaxThreads) {
943 }
944 fMaxThreads = nthreads + 1;
945 if (fMaxThreads > 0) {
948 }
949}
950
951////////////////////////////////////////////////////////////////////////////////
952
954{
955 if (!fMaxThreads)
956 return;
957 fgMutex.lock();
958 TIter next(fVolumes);
959 TGeoVolume *vol;
960 while ((vol = (TGeoVolume *)next()))
961 vol->ClearThreadData();
962 fgMutex.unlock();
963}
964
965////////////////////////////////////////////////////////////////////////////////
966/// Create thread private data for all geometry objects.
967
969{
970 if (!fMaxThreads)
971 return;
972 fgMutex.lock();
973 TIter next(fVolumes);
974 TGeoVolume *vol;
975 while ((vol = (TGeoVolume *)next()))
977 fgMutex.unlock();
978}
979
980////////////////////////////////////////////////////////////////////////////////
981/// Clear the current map of threads. This will be filled again by the calling
982/// threads via ThreadId calls.
983
985{
987 return;
988 fgMutex.lock();
989 if (!fgThreadId->empty())
990 fgThreadId->clear();
991 fgNumThreads = 0;
992 fgMutex.unlock();
993}
994
995////////////////////////////////////////////////////////////////////////////////
996/// Translates the current thread id to an ordinal number. This can be used to
997/// manage data which is specific for a given thread.
998
1000{
1001 TTHREAD_TLS(Int_t) tid = -1;
1002 Int_t ttid = tid; // TTHREAD_TLS_GET(Int_t,tid);
1003 if (ttid > -1)
1004 return ttid;
1006 return 0;
1007 std::thread::id threadId = std::this_thread::get_id();
1008 TGeoManager::ThreadsMapIt_t it = fgThreadId->find(threadId);
1009 if (it != fgThreadId->end())
1010 return it->second;
1011 // Map needs to be updated.
1012 fgMutex.lock();
1013 (*fgThreadId)[threadId] = fgNumThreads;
1014 tid = fgNumThreads; // TTHREAD_TLS_SET(Int_t,tid,fgNumThreads);
1015 ttid = fgNumThreads++;
1016 fgMutex.unlock();
1017 return ttid;
1018}
1019
1020////////////////////////////////////////////////////////////////////////////////
1021/// Describe how to browse this object.
1022
1024{
1025 if (!b)
1026 return;
1027 if (fMaterials)
1028 b->Add(fMaterials, "Materials");
1029 if (fMedia)
1030 b->Add(fMedia, "Media");
1031 if (fMatrices)
1032 b->Add(fMatrices, "Local transformations");
1033 if (fOverlaps)
1034 b->Add(fOverlaps, "Illegal overlaps");
1035 if (fTracks)
1036 b->Add(fTracks, "Tracks");
1037 if (fMasterVolume)
1038 b->Add(fMasterVolume, "Master Volume", fMasterVolume->IsVisible());
1039 if (fTopVolume)
1040 b->Add(fTopVolume, "Top Volume", fTopVolume->IsVisible());
1041 if (fTopNode)
1042 b->Add(fTopNode);
1043 TString browserImp(gEnv->GetValue("Browser.Name", "TRootBrowserLite"));
1044 TQObject::Connect(browserImp.Data(), "Checked(TObject*,Bool_t)", "TGeoManager", this,
1045 "SetVisibility(TObject*,Bool_t)");
1046}
1047
1048////////////////////////////////////////////////////////////////////////////////
1049/// Append a pad for this geometry.
1050
1052{
1053 AppendPad("");
1055}
1056
1057////////////////////////////////////////////////////////////////////////////////
1058/// Set visibility for a volume.
1059
1061{
1062 if (obj->IsA() == TGeoVolume::Class()) {
1063 TGeoVolume *vol = (TGeoVolume *)obj;
1064 vol->SetVisibility(vis);
1065 } else {
1066 if (obj->InheritsFrom(TGeoNode::Class())) {
1067 TGeoNode *node = (TGeoNode *)obj;
1068 node->SetVisibility(vis);
1069 } else
1070 return;
1071 }
1073}
1074
1075////////////////////////////////////////////////////////////////////////////////
1076/// Get the new 'bombed' translation vector according current exploded view mode.
1077
1079{
1080 if (fPainter)
1081 fPainter->BombTranslation(tr, bombtr);
1082 return;
1083}
1084
1085////////////////////////////////////////////////////////////////////////////////
1086/// Get the new 'unbombed' translation vector according current exploded view mode.
1087
1089{
1090 if (fPainter)
1091 fPainter->UnbombTranslation(tr, bombtr);
1092 return;
1093}
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Backup the current state without affecting the cache stack.
1097
1099{
1101}
1102
1103////////////////////////////////////////////////////////////////////////////////
1104/// Restore a backed-up state without affecting the cache stack.
1105
1107{
1109}
1110
1111////////////////////////////////////////////////////////////////////////////////
1112/// Register a matrix to the list of matrices. It will be cleaned-up at the
1113/// destruction TGeoManager.
1114
1116{
1117 return TGeoBuilder::Instance(this)->RegisterMatrix((TGeoMatrix *)matrix);
1118}
1119
1120////////////////////////////////////////////////////////////////////////////////
1121/// Replaces all occurrences of VORIG with VNEW in the geometry tree. The volume VORIG
1122/// is not replaced from the list of volumes, but all node referencing it will reference
1123/// VNEW instead. Returns number of occurrences changed.
1124
1126{
1127 Int_t nref = 0;
1128 if (!vorig || !vnew)
1129 return nref;
1130 TGeoMedium *morig = vorig->GetMedium();
1131 Bool_t checkmed = kFALSE;
1132 if (morig)
1133 checkmed = kTRUE;
1134 TGeoMedium *mnew = vnew->GetMedium();
1135 // Try to limit the damage produced by incorrect usage.
1136 if (!mnew && !vnew->IsAssembly()) {
1137 Error("ReplaceVolume", "Replacement volume %s has no medium and it is not an assembly", vnew->GetName());
1138 return nref;
1139 }
1140 if (mnew && checkmed) {
1141 if (mnew->GetId() != morig->GetId())
1142 Warning("ReplaceVolume", "Replacement volume %s has different medium than original volume %s", vnew->GetName(),
1143 vorig->GetName());
1144 checkmed = kFALSE;
1145 }
1146
1147 // Medium checking now performed only if replacement is an assembly and old volume a real one.
1148 // Check result is dependent on positioning.
1149 Int_t nvol = fVolumes->GetEntriesFast();
1150 Int_t i, j, nd;
1151 Int_t ierr = 0;
1152 TGeoVolume *vol;
1153 TGeoNode *node;
1154 TGeoVoxelFinder *voxels;
1155 for (i = 0; i < nvol; i++) {
1156 vol = (TGeoVolume *)fVolumes->At(i);
1157 if (!vol)
1158 continue;
1159 if (vol == vorig || vol == vnew)
1160 continue;
1161 nd = vol->GetNdaughters();
1162 for (j = 0; j < nd; j++) {
1163 node = vol->GetNode(j);
1164 if (node->GetVolume() == vorig) {
1165 if (checkmed) {
1166 mnew = node->GetMotherVolume()->GetMedium();
1167 if (mnew && mnew->GetId() != morig->GetId())
1168 ierr++;
1169 }
1170 nref++;
1171 if (node->IsOverlapping()) {
1172 node->SetOverlapping(kFALSE);
1173 Info("ReplaceVolume", "%s replaced with assembly and declared NON-OVERLAPPING!", node->GetName());
1174 }
1175 node->SetVolume(vnew);
1176 voxels = node->GetMotherVolume()->GetVoxels();
1177 if (voxels)
1178 voxels->SetNeedRebuild();
1179 } else {
1180 if (node->GetMotherVolume() == vorig) {
1181 nref++;
1182 node->SetMotherVolume(vnew);
1183 if (node->IsOverlapping()) {
1184 node->SetOverlapping(kFALSE);
1185 Info("ReplaceVolume", "%s inside substitute assembly %s declared NON-OVERLAPPING!", node->GetName(),
1186 vnew->GetName());
1187 }
1188 }
1189 }
1190 }
1191 }
1192 if (ierr)
1193 Warning("ReplaceVolume",
1194 "Volumes should not be replaced with assemblies if they are positioned in containers having a different "
1195 "medium ID.\n %i occurrences for assembly replacing volume %s",
1196 ierr, vorig->GetName());
1197 return nref;
1198}
1199
1200////////////////////////////////////////////////////////////////////////////////
1201/// Transform all volumes named VNAME to assemblies. The volumes must be virtual.
1202
1204{
1205 TGeoVolume *toTransform = FindVolumeFast(vname);
1206 if (!toTransform) {
1207 Warning("TransformVolumeToAssembly", "Volume %s not found", vname);
1208 return 0;
1209 }
1210 Int_t index = fVolumes->IndexOf(toTransform);
1211 Int_t count = 0;
1212 Int_t indmax = fVolumes->GetEntries();
1213 Bool_t replace = kTRUE;
1214 TGeoVolume *transformed;
1215 while (index < indmax) {
1216 if (replace) {
1217 replace = kFALSE;
1218 transformed = TGeoVolumeAssembly::MakeAssemblyFromVolume(toTransform);
1219 if (transformed) {
1220 ReplaceVolume(toTransform, transformed);
1221 count++;
1222 } else {
1223 if (toTransform->IsAssembly())
1224 Warning("TransformVolumeToAssembly", "Volume %s already assembly", toTransform->GetName());
1225 if (!toTransform->GetNdaughters())
1226 Warning("TransformVolumeToAssembly", "Volume %s has no daughters, cannot transform",
1227 toTransform->GetName());
1228 if (toTransform->IsVolumeMulti())
1229 Warning("TransformVolumeToAssembly", "Volume %s divided, cannot transform", toTransform->GetName());
1230 }
1231 }
1232 index++;
1233 if (index >= indmax)
1234 return count;
1235 toTransform = (TGeoVolume *)fVolumes->At(index);
1236 if (!strcmp(toTransform->GetName(), vname))
1237 replace = kTRUE;
1238 }
1239 return count;
1240}
1241
1242////////////////////////////////////////////////////////////////////////////////
1243/// Create a new volume by dividing an existing one (GEANT3 like)
1244///
1245/// Divides MOTHER into NDIV divisions called NAME
1246/// along axis IAXIS starting at coordinate value START
1247/// and having size STEP. The created volumes will have tracking
1248/// media ID=NUMED (if NUMED=0 -> same media as MOTHER)
1249/// The behavior of the division operation can be triggered using OPTION :
1250///
1251/// OPTION (case insensitive) :
1252/// - N - divide all range in NDIV cells (same effect as STEP<=0) (GSDVN in G3)
1253/// - NX - divide range starting with START in NDIV cells (GSDVN2 in G3)
1254/// - S - divide all range with given STEP. NDIV is computed and divisions will be centered
1255/// in full range (same effect as NDIV<=0) (GSDVS, GSDVT in G3)
1256/// - SX - same as DVS, but from START position. (GSDVS2, GSDVT2 in G3)
1257
1258TGeoVolume *TGeoManager::Division(const char *name, const char *mother, Int_t iaxis, Int_t ndiv, Double_t start,
1259 Double_t step, Int_t numed, Option_t *option)
1260{
1261 return TGeoBuilder::Instance(this)->Division(name, mother, iaxis, ndiv, start, step, numed, option);
1262}
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// Create rotation matrix named 'mat<index>'.
1266///
1267/// - index rotation matrix number
1268/// - theta1 polar angle for axis X
1269/// - phi1 azimuthal angle for axis X
1270/// - theta2 polar angle for axis Y
1271/// - phi2 azimuthal angle for axis Y
1272/// - theta3 polar angle for axis Z
1273/// - phi3 azimuthal angle for axis Z
1274///
1275
1277 Double_t phi3)
1278{
1279 TGeoBuilder::Instance(this)->Matrix(index, theta1, phi1, theta2, phi2, theta3, phi3);
1280}
1281
1282////////////////////////////////////////////////////////////////////////////////
1283/// Create material with given A, Z and density, having an unique id.
1284
1286 Double_t intlen)
1287{
1288 return TGeoBuilder::Instance(this)->Material(name, a, z, dens, uid, radlen, intlen);
1289}
1290
1291////////////////////////////////////////////////////////////////////////////////
1292/// Create mixture OR COMPOUND IMAT as composed by THE BASIC nelem
1293/// materials defined by arrays A,Z and WMAT, having an unique id.
1294
1296TGeoManager::Mixture(const char *name, Float_t *a, Float_t *z, Double_t dens, Int_t nelem, Float_t *wmat, Int_t uid)
1297{
1298 return TGeoBuilder::Instance(this)->Mixture(name, a, z, dens, nelem, wmat, uid);
1299}
1300
1301////////////////////////////////////////////////////////////////////////////////
1302/// Create mixture OR COMPOUND IMAT as composed by THE BASIC nelem
1303/// materials defined by arrays A,Z and WMAT, having an unique id.
1304
1306TGeoManager::Mixture(const char *name, Double_t *a, Double_t *z, Double_t dens, Int_t nelem, Double_t *wmat, Int_t uid)
1307{
1308 return TGeoBuilder::Instance(this)->Mixture(name, a, z, dens, nelem, wmat, uid);
1309}
1310
1311////////////////////////////////////////////////////////////////////////////////
1312/// Create tracking medium
1313///
1314/// - numed tracking medium number assigned
1315/// - name tracking medium name
1316/// - nmat material number
1317/// - isvol sensitive volume flag
1318/// - ifield magnetic field
1319/// - fieldm max. field value (kilogauss)
1320/// - tmaxfd max. angle due to field (deg/step)
1321/// - stemax max. step allowed
1322/// - deemax max. fraction of energy lost in a step
1323/// - epsil tracking precision (cm)
1324/// - stmin min. step due to continuous processes (cm)
1325///
1326/// - ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1327/// - ifield = 1 if tracking performed with g3rkuta; ifield = 2 if tracking
1328/// performed with g3helix; ifield = 3 if tracking performed with g3helx3.
1329///
1330
1331TGeoMedium *TGeoManager::Medium(const char *name, Int_t numed, Int_t nmat, Int_t isvol, Int_t ifield, Double_t fieldm,
1332 Double_t tmaxfd, Double_t stemax, Double_t deemax, Double_t epsil, Double_t stmin)
1333{
1334 return TGeoBuilder::Instance(this)->Medium(name, numed, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax, epsil,
1335 stmin);
1336}
1337
1338////////////////////////////////////////////////////////////////////////////////
1339/// Create a node called `<name_nr>` pointing to the volume called `<name>`
1340/// as daughter of the volume called `<mother>` (gspos). The relative matrix is
1341/// made of : a translation (x,y,z) and a rotation matrix named `<matIROT>`.
1342/// In case npar>0, create the volume to be positioned in mother, according
1343/// its actual parameters (gsposp).
1344/// - NAME Volume name
1345/// - NUMBER Copy number of the volume
1346/// - MOTHER Mother volume name
1347/// - X X coord. of the volume in mother ref. sys.
1348/// - Y Y coord. of the volume in mother ref. sys.
1349/// - Z Z coord. of the volume in mother ref. sys.
1350/// - IROT Rotation matrix number w.r.t. mother ref. sys.
1351/// - ISONLY ONLY/MANY flag
1352
1353void TGeoManager::Node(const char *name, Int_t nr, const char *mother, Double_t x, Double_t y, Double_t z, Int_t irot,
1354 Bool_t isOnly, Float_t *upar, Int_t npar)
1355{
1356 TGeoBuilder::Instance(this)->Node(name, nr, mother, x, y, z, irot, isOnly, upar, npar);
1357}
1358
1359////////////////////////////////////////////////////////////////////////////////
1360/// Create a node called `<name_nr>` pointing to the volume called `<name>`
1361/// as daughter of the volume called `<mother>` (gspos). The relative matrix is
1362/// made of : a translation (x,y,z) and a rotation matrix named `<matIROT>`.
1363/// In case npar>0, create the volume to be positioned in mother, according
1364/// its actual parameters (gsposp).
1365/// - NAME Volume name
1366/// - NUMBER Copy number of the volume
1367/// - MOTHER Mother volume name
1368/// - X X coord. of the volume in mother ref. sys.
1369/// - Y Y coord. of the volume in mother ref. sys.
1370/// - Z Z coord. of the volume in mother ref. sys.
1371/// - IROT Rotation matrix number w.r.t. mother ref. sys.
1372/// - ISONLY ONLY/MANY flag
1373
1374void TGeoManager::Node(const char *name, Int_t nr, const char *mother, Double_t x, Double_t y, Double_t z, Int_t irot,
1375 Bool_t isOnly, Double_t *upar, Int_t npar)
1376{
1377 TGeoBuilder::Instance(this)->Node(name, nr, mother, x, y, z, irot, isOnly, upar, npar);
1378}
1379
1380////////////////////////////////////////////////////////////////////////////////
1381/// Create a volume in GEANT3 style.
1382/// - NAME Volume name
1383/// - SHAPE Volume type
1384/// - NMED Tracking medium number
1385/// - NPAR Number of shape parameters
1386/// - UPAR Vector containing shape parameters
1387
1388TGeoVolume *TGeoManager::Volume(const char *name, const char *shape, Int_t nmed, Float_t *upar, Int_t npar)
1389{
1390 return TGeoBuilder::Instance(this)->Volume(name, shape, nmed, upar, npar);
1391}
1392
1393////////////////////////////////////////////////////////////////////////////////
1394/// Create a volume in GEANT3 style.
1395/// - NAME Volume name
1396/// - SHAPE Volume type
1397/// - NMED Tracking medium number
1398/// - NPAR Number of shape parameters
1399/// - UPAR Vector containing shape parameters
1400
1401TGeoVolume *TGeoManager::Volume(const char *name, const char *shape, Int_t nmed, Double_t *upar, Int_t npar)
1402{
1403 return TGeoBuilder::Instance(this)->Volume(name, shape, nmed, upar, npar);
1404}
1405
1406////////////////////////////////////////////////////////////////////////////////
1407/// Assigns uid's for all materials,media and matrices.
1408
1410{
1411 Int_t index = 1;
1412 TIter next(fMaterials);
1413 TGeoMaterial *mater;
1414 while ((mater = (TGeoMaterial *)next())) {
1415 mater->SetUniqueID(index++);
1417 }
1418 index = 1;
1419 TIter next1(fMedia);
1420 TGeoMedium *med;
1421 while ((med = (TGeoMedium *)next1())) {
1422 med->SetUniqueID(index++);
1424 }
1425 index = 1;
1426 TIter next2(fShapes);
1427 TGeoShape *shape;
1428 while ((shape = (TGeoShape *)next2())) {
1429 shape->SetUniqueID(index++);
1430 if (shape->IsComposite())
1431 ((TGeoCompositeShape *)shape)->GetBoolNode()->RegisterMatrices();
1432 }
1433
1434 TIter next3(fMatrices);
1435 TGeoMatrix *matrix;
1436 while ((matrix = (TGeoMatrix *)next3())) {
1437 matrix->RegisterYourself();
1438 }
1439 TIter next4(fMatrices);
1440 index = 1;
1441 while ((matrix = (TGeoMatrix *)next4())) {
1442 matrix->SetUniqueID(index++);
1444 }
1445 TIter next5(fVolumes);
1446 TGeoVolume *vol;
1447 while ((vol = (TGeoVolume *)next5()))
1448 vol->UnmarkSaved();
1449}
1450
1451////////////////////////////////////////////////////////////////////////////////
1452/// Reset all attributes to default ones. Default attributes for visualization
1453/// are those defined before closing the geometry.
1454
1456{
1457 if (gPad)
1458 delete gPad;
1459 gPad = nullptr;
1460 SetVisOption(0);
1461 SetVisLevel(3);
1462 SetExplodedView(0);
1464 if (!gStyle)
1465 return;
1466 TIter next(fVolumes);
1467 TGeoVolume *vol = nullptr;
1468 while ((vol = (TGeoVolume *)next())) {
1469 if (!vol->IsVisTouched())
1470 continue;
1471 vol->SetVisTouched(kFALSE);
1472 }
1473}
1474////////////////////////////////////////////////////////////////////////////////
1475/// Closing geometry implies checking the geometry validity, fixing shapes
1476/// with negative parameters (run-time shapes)building the cache manager,
1477/// voxelizing all volumes, counting the total number of physical nodes and
1478/// registering the manager class to the browser.
1479
1481{
1482 if (fClosed) {
1483 Warning("CloseGeometry", "geometry already closed");
1484 return;
1485 }
1486 if (!fMasterVolume) {
1487 Error("CloseGeometry", "you MUST call SetTopVolume() first !");
1488 return;
1489 }
1490 if (!gROOT->GetListOfGeometries()->FindObject(this))
1491 gROOT->GetListOfGeometries()->Add(this);
1492 if (!gROOT->GetListOfBrowsables()->FindObject(this))
1493 gROOT->GetListOfBrowsables()->Add(this);
1494 // TSeqCollection *brlist = gROOT->GetListOfBrowsers();
1495 // TIter next(brlist);
1496 // TBrowser *browser = 0;
1497 // while ((browser=(TBrowser*)next())) browser->Refresh();
1498 TString opt(option);
1499 opt.ToLower();
1500 // Bool_t dummy = opt.Contains("d");
1501 Bool_t nodeid = opt.Contains("i");
1502 // Create a geometry navigator if not present
1503 TGeoNavigator *nav = nullptr;
1504 Int_t nnavigators = 0;
1505 // Check if the geometry is streamed from file
1506 if (fIsGeomReading) {
1507 if (fgVerboseLevel > 0)
1508 Info("CloseGeometry", "Geometry loaded from file...");
1510 if (!fElementTable)
1512 if (!fTopNode) {
1513 if (!fMasterVolume) {
1514 Error("CloseGeometry", "Master volume not streamed");
1515 return;
1516 }
1518 if (fStreamVoxels && fgVerboseLevel > 0)
1519 Info("CloseGeometry", "Voxelization retrieved from file");
1520 }
1521 // Create a geometry navigator if not present
1522 if (!GetCurrentNavigator())
1524 nnavigators = GetListOfNavigators()->GetEntriesFast();
1525 if (!opt.Contains("nv")) {
1526 Voxelize("ALL");
1527 }
1528 CountLevels();
1529 for (Int_t i = 0; i < nnavigators; i++) {
1530 nav = (TGeoNavigator *)GetListOfNavigators()->At(i);
1531 nav->GetCache()->BuildInfoBranch();
1532 if (nodeid)
1533 nav->GetCache()->BuildIdArray();
1534 }
1535 if (!fHashVolumes) {
1536 Int_t nvol = fVolumes->GetEntriesFast();
1537 Int_t ngvol = fGVolumes->GetEntriesFast();
1538 fHashVolumes = new THashList(nvol + 1);
1539 fHashGVolumes = new THashList(ngvol + 1);
1540 Int_t i;
1541 for (i = 0; i < ngvol; i++)
1543 for (i = 0; i < nvol; i++)
1545 }
1546 fClosed = kTRUE;
1547 if (fParallelWorld) {
1548 if (fgVerboseLevel > 0)
1549 Info("CloseGeometry", "Recreating parallel world %s ...", fParallelWorld->GetName());
1551 }
1552
1553 if (fgVerboseLevel > 0)
1554 Info("CloseGeometry", "%i nodes/ %i volume UID's in %s", fNNodes, fUniqueVolumes->GetEntriesFast() - 1,
1555 GetTitle());
1556 if (fgVerboseLevel > 0)
1557 Info("CloseGeometry", "----------------modeler ready----------------");
1558 return;
1559 }
1560
1561 // Create a geometry navigator if not present
1562 if (!GetCurrentNavigator())
1564 nnavigators = GetListOfNavigators()->GetEntriesFast();
1566 CheckGeometry();
1567 if (fgVerboseLevel > 0)
1568 Info("CloseGeometry", "Counting nodes...");
1569 fNNodes = CountNodes();
1570 fNLevel = fMasterVolume->CountNodes(1, 3) + 1;
1571 if (fNLevel < 30)
1572 fNLevel = 100;
1573
1574 // BuildIdArray();
1575 // avoid voxelization if requested to speed up geometry startup
1576 if (!opt.Contains("nv")) {
1577 Voxelize("ALL");
1578 } else {
1579 TGeoVolume *vol;
1580 TIter next(fVolumes);
1581 while ((vol = (TGeoVolume *)next())) {
1582 vol->SortNodes();
1583 }
1584 }
1585 if (fgVerboseLevel > 0)
1586 Info("CloseGeometry", "Building cache...");
1587 CountLevels();
1588 for (Int_t i = 0; i < nnavigators; i++) {
1589 nav = (TGeoNavigator *)GetListOfNavigators()->At(i);
1590 nav->GetCache()->BuildInfoBranch();
1591 if (nodeid)
1592 nav->GetCache()->BuildIdArray();
1593 }
1594 fClosed = kTRUE;
1595 if (fgVerboseLevel > 0) {
1596 Info("CloseGeometry", "%i nodes/ %i volume UID's in %s", fNNodes, fUniqueVolumes->GetEntriesFast() - 1,
1597 GetTitle());
1598 Info("CloseGeometry", "----------------modeler ready----------------");
1599 }
1600}
1601
1602////////////////////////////////////////////////////////////////////////////////
1603/// Clear the list of overlaps.
1604
1606{
1607 if (fOverlaps) {
1608 fOverlaps->Delete();
1609 delete fOverlaps;
1610 }
1611 fOverlaps = new TObjArray();
1612}
1613
1614////////////////////////////////////////////////////////////////////////////////
1615/// Remove a shape from the list of shapes.
1616
1618{
1619 if (fShapes->FindObject(shape))
1620 fShapes->Remove((TGeoShape *)shape);
1621 delete shape;
1622}
1623
1624////////////////////////////////////////////////////////////////////////////////
1625/// Clean temporary volumes and shapes from garbage collection.
1626
1628{
1629 if (!fGVolumes && !fGShapes)
1630 return;
1631 Int_t i, nentries;
1632 if (fGVolumes) {
1634 TGeoVolume *vol = nullptr;
1635 for (i = 0; i < nentries; i++) {
1636 vol = (TGeoVolume *)fGVolumes->At(i);
1637 if (vol)
1638 vol->SetFinder(nullptr);
1639 }
1640 fGVolumes->Delete();
1641 delete fGVolumes;
1642 fGVolumes = nullptr;
1643 }
1644 if (fGShapes) {
1645 fGShapes->Delete();
1646 delete fGShapes;
1647 fGShapes = nullptr;
1648 }
1649}
1650
1651////////////////////////////////////////////////////////////////////////////////
1652/// Change current path to point to the node having this id.
1653/// Node id has to be in range : 0 to fNNodes-1 (no check for performance reasons)
1654
1656{
1657 GetCurrentNavigator()->CdNode(nodeid);
1658}
1659
1660////////////////////////////////////////////////////////////////////////////////
1661/// Get the unique ID of the current node.
1662
1664{
1666}
1667
1668////////////////////////////////////////////////////////////////////////////////
1669/// Make top level node the current node. Updates the cache accordingly.
1670/// Determine the overlapping state of current node.
1671
1673{
1675}
1676
1677////////////////////////////////////////////////////////////////////////////////
1678/// Go one level up in geometry. Updates cache accordingly.
1679/// Determine the overlapping state of current node.
1680
1682{
1684}
1685
1686////////////////////////////////////////////////////////////////////////////////
1687/// Make a daughter of current node current. Can be called only with a valid
1688/// daughter index (no check). Updates cache accordingly.
1689
1691{
1693}
1694
1695////////////////////////////////////////////////////////////////////////////////
1696/// Do a cd to the node found next by FindNextBoundary
1697
1699{
1701}
1702
1703////////////////////////////////////////////////////////////////////////////////
1704/// Browse the tree of nodes starting from fTopNode according to pathname.
1705/// Changes the path accordingly.
1706
1707Bool_t TGeoManager::cd(const char *path)
1708{
1709 return GetCurrentNavigator()->cd(path);
1710}
1711
1712////////////////////////////////////////////////////////////////////////////////
1713/// Check if a geometry path is valid without changing the state of the current navigator.
1714
1715Bool_t TGeoManager::CheckPath(const char *path) const
1716{
1717 return GetCurrentNavigator()->CheckPath(path);
1718}
1719
1720////////////////////////////////////////////////////////////////////////////////
1721/// Convert all reflections in geometry to normal rotations + reflected shapes.
1722
1724{
1725 if (!fTopNode)
1726 return;
1727 if (fgVerboseLevel > 0)
1728 Info("ConvertReflections", "Converting reflections in: %s - %s ...", GetName(), GetTitle());
1730 TGeoNode *node;
1731 TGeoNodeMatrix *nodematrix;
1732 TGeoMatrix *matrix, *mclone;
1733 TGeoVolume *reflected;
1734 while ((node = next())) {
1735 matrix = node->GetMatrix();
1736 if (matrix->IsReflection()) {
1737 // printf("%s before\n", node->GetName());
1738 // matrix->Print();
1739 mclone = new TGeoCombiTrans(*matrix);
1740 mclone->RegisterYourself();
1741 // Reflect just the rotation component
1742 mclone->ReflectZ(kFALSE, kTRUE);
1743 nodematrix = (TGeoNodeMatrix *)node;
1744 nodematrix->SetMatrix(mclone);
1745 // printf("%s after\n", node->GetName());
1746 // node->GetMatrix()->Print();
1747 reflected = node->GetVolume()->MakeReflectedVolume();
1748 node->SetVolume(reflected);
1749 }
1750 }
1751 if (fgVerboseLevel > 0)
1752 Info("ConvertReflections", "Done");
1753}
1754
1755////////////////////////////////////////////////////////////////////////////////
1756/// Count maximum number of nodes per volume, maximum depth and maximum
1757/// number of xtru vertices.
1758
1760{
1761 if (!fTopNode) {
1762 Error("CountLevels", "Top node not defined.");
1763 return;
1764 }
1766 Bool_t fixrefs = fIsGeomReading && (fMasterVolume->GetRefCount() == 1);
1767 if (fMasterVolume->GetRefCount() > 1)
1769 if (fgVerboseLevel > 1 && fixrefs)
1770 Info("CountLevels", "Fixing volume reference counts");
1771 TGeoNode *node;
1772 Int_t maxlevel = 1;
1773 Int_t maxnodes = fTopVolume->GetNdaughters();
1774 Int_t maxvertices = 1;
1775 while ((node = next())) {
1776 if (fixrefs) {
1777 node->GetVolume()->Grab();
1778 for (Int_t ibit = 10; ibit < 14; ibit++) {
1779 node->SetBit(BIT(ibit + 4), node->TestBit(BIT(ibit)));
1780 // node->ResetBit(BIT(ibit)); // cannot overwrite old crap for reproducibility
1781 }
1782 }
1783 if (node->GetNdaughters() > maxnodes)
1784 maxnodes = node->GetNdaughters();
1785 if (next.GetLevel() > maxlevel)
1786 maxlevel = next.GetLevel();
1787 if (node->GetVolume()->GetShape()->IsA() == TGeoXtru::Class()) {
1788 TGeoXtru *xtru = (TGeoXtru *)node->GetVolume()->GetShape();
1789 if (xtru->GetNvert() > maxvertices)
1790 maxvertices = xtru->GetNvert();
1791 }
1792 }
1793 fgMaxLevel = maxlevel;
1794 fgMaxDaughters = maxnodes;
1795 fgMaxXtruVert = maxvertices;
1796 if (fgVerboseLevel > 0)
1797 Info("CountLevels", "max level = %d, max placements = %d", fgMaxLevel, fgMaxDaughters);
1798}
1799
1800////////////////////////////////////////////////////////////////////////////////
1801/// Count the total number of nodes starting from a volume, nlevels down.
1802
1804{
1805 TGeoVolume *top;
1806 if (!vol) {
1807 top = fTopVolume;
1808 } else {
1809 top = (TGeoVolume *)vol;
1810 }
1811 Int_t count = top->CountNodes(nlevels, option);
1812 return count;
1813}
1814
1815////////////////////////////////////////////////////////////////////////////////
1816/// Set default angles for a given view.
1817
1819{
1820 if (fPainter)
1822}
1823
1824////////////////////////////////////////////////////////////////////////////////
1825/// Draw current point in the same view.
1826
1828{
1829 if (fPainter)
1830 fPainter->DrawCurrentPoint(color);
1831}
1832
1833////////////////////////////////////////////////////////////////////////////////
1834/// Draw animation of tracks
1835
1837{
1840 if (tmin < 0 || tmin >= tmax || nframes < 1)
1841 return;
1843 box[0] = box[1] = box[2] = 0;
1844 box[3] = box[4] = box[5] = 100;
1845 Double_t dt = (tmax - tmin) / Double_t(nframes);
1846 Double_t delt = 2E-9;
1847 Double_t t = tmin;
1848 Int_t i, j;
1849 TString opt(option);
1850 Bool_t save = kFALSE, geomanim = kFALSE;
1851 TString fname;
1852 if (opt.Contains("/S"))
1853 save = kTRUE;
1854
1855 if (opt.Contains("/G"))
1856 geomanim = kTRUE;
1857 SetTminTmax(0, 0);
1858 DrawTracks(opt.Data());
1859 Double_t start[6] = {0, 0, 0, 0, 0, 0};
1860 Double_t end[6] = {0, 0, 0, 0, 0, 0};
1861 Double_t dd[6] = {0, 0, 0, 0, 0, 0};
1862 Double_t dlat = 0, dlong = 0, dpsi = 0;
1863 if (geomanim) {
1864 fPainter->EstimateCameraMove(tmin + 5 * dt, tmin + 15 * dt, start, end);
1865 for (i = 0; i < 3; i++) {
1866 start[i + 3] = 20 + 1.3 * start[i + 3];
1867 end[i + 3] = 20 + 0.9 * end[i + 3];
1868 }
1869 for (i = 0; i < 6; i++) {
1870 dd[i] = (end[i] - start[i]) / 10.;
1871 }
1872 memcpy(box, start, 6 * sizeof(Double_t));
1873 fPainter->GetViewAngles(dlong, dlat, dpsi);
1874 dlong = (-206 - dlong) / Double_t(nframes);
1875 dlat = (126 - dlat) / Double_t(nframes);
1876 dpsi = (75 - dpsi) / Double_t(nframes);
1878 }
1879
1880 for (i = 0; i < nframes; i++) {
1881 if (t - delt < 0)
1882 SetTminTmax(t - delt, t);
1883 else
1884 gGeoManager->SetTminTmax(t - delt, t);
1885 if (geomanim) {
1886 for (j = 0; j < 6; j++)
1887 box[j] += dd[j];
1888 fPainter->GrabFocus(1, dlong, dlat, dpsi);
1889 } else {
1890 ModifiedPad();
1891 }
1892 if (save) {
1893 fname = TString::Format("anim%04d.gif", i);
1894 gPad->Print(fname);
1895 }
1896 t += dt;
1897 }
1899}
1900
1901////////////////////////////////////////////////////////////////////////////////
1902/// Draw tracks over the geometry, according to option. By default, only
1903/// primaries are drawn. See TGeoTrack::Draw() for additional options.
1904
1906{
1907 TVirtualGeoTrack *track;
1908 // SetVisLevel(1);
1909 // SetVisOption(1);
1911 for (Int_t i = 0; i < fNtracks; i++) {
1912 track = GetTrack(i);
1913 if (track)
1914 track->Draw(option);
1915 }
1917 ModifiedPad();
1918}
1919
1920////////////////////////////////////////////////////////////////////////////////
1921/// Draw current path
1922
1923void TGeoManager::DrawPath(const char *path, Option_t *option)
1924{
1925 if (!fTopVolume)
1926 return;
1928 GetGeomPainter()->DrawPath(path, option);
1929}
1930
1931////////////////////////////////////////////////////////////////////////////////
1932/// Draw random points in the bounding box of a volume.
1933
1935{
1936 GetGeomPainter()->RandomPoints((TGeoVolume *)vol, npoints, option);
1937}
1938
1939////////////////////////////////////////////////////////////////////////////////
1940/// Check time of finding "Where am I" for n points.
1941
1943{
1944 GetGeomPainter()->Test(npoints, option);
1945}
1946
1947////////////////////////////////////////////////////////////////////////////////
1948/// Geometry overlap checker based on sampling.
1949
1950void TGeoManager::TestOverlaps(const char *path)
1951{
1953}
1954
1955////////////////////////////////////////////////////////////////////////////////
1956/// Fill volume names of current branch into an array.
1957
1959{
1961}
1962
1963////////////////////////////////////////////////////////////////////////////////
1964/// Get name for given pdg code;
1965
1966const char *TGeoManager::GetPdgName(Int_t pdg) const
1967{
1968 static char defaultname[5] = {"XXX"};
1969 if (!fPdgNames || !pdg)
1970 return defaultname;
1971 for (Int_t i = 0; i < fNpdg; i++) {
1972 if (fPdgId[i] == pdg)
1973 return fPdgNames->At(i)->GetName();
1974 }
1975 return defaultname;
1976}
1977
1978////////////////////////////////////////////////////////////////////////////////
1979/// Set a name for a particle having a given pdg.
1980
1981void TGeoManager::SetPdgName(Int_t pdg, const char *name)
1982{
1983 if (!pdg)
1984 return;
1985 if (!fPdgNames) {
1986 fPdgNames = new TObjArray(1024);
1987 }
1988 if (!strcmp(name, GetPdgName(pdg)))
1989 return;
1990 // store pdg name
1991 if (fNpdg > 1023) {
1992 Warning("SetPdgName", "No more than 256 different pdg codes allowed");
1993 return;
1994 }
1995 fPdgId[fNpdg] = pdg;
1996 TNamed *pdgname = new TNamed(name, "");
1997 fPdgNames->AddAtAndExpand(pdgname, fNpdg++);
1998}
1999
2000////////////////////////////////////////////////////////////////////////////////
2001/// Get GDML matrix with a given name;
2002
2004{
2006}
2007
2008////////////////////////////////////////////////////////////////////////////////
2009/// Add GDML matrix;
2011{
2012 if (GetGDMLMatrix(mat->GetName())) {
2013 Error("AddGDMLMatrix", "Matrix %s already added to manager", mat->GetName());
2014 return;
2015 }
2016 fGDMLMatrices->Add(mat);
2017}
2018
2019////////////////////////////////////////////////////////////////////////////////
2020/// Get optical surface with a given name;
2021
2023{
2025}
2026
2027////////////////////////////////////////////////////////////////////////////////
2028/// Add optical surface;
2030{
2031 if (GetOpticalSurface(optsurf->GetName())) {
2032 Error("AddOpticalSurface", "Surface %s already added to manager", optsurf->GetName());
2033 return;
2034 }
2035 fOpticalSurfaces->Add(optsurf);
2036}
2037
2038////////////////////////////////////////////////////////////////////////////////
2039/// Get skin surface with a given name;
2040
2042{
2044}
2045
2046////////////////////////////////////////////////////////////////////////////////
2047/// Add skin surface;
2049{
2050 if (GetSkinSurface(surf->GetName())) {
2051 Error("AddSkinSurface", "Surface %s already added to manager", surf->GetName());
2052 return;
2053 }
2054 fSkinSurfaces->Add(surf);
2055}
2056
2057////////////////////////////////////////////////////////////////////////////////
2058/// Get border surface with a given name;
2059
2061{
2063}
2064
2065////////////////////////////////////////////////////////////////////////////////
2066/// Add border surface;
2068{
2069 if (GetBorderSurface(surf->GetName())) {
2070 Error("AddBorderSurface", "Surface %s already added to manager", surf->GetName());
2071 return;
2072 }
2073 fBorderSurfaces->Add(surf);
2074}
2075
2076////////////////////////////////////////////////////////////////////////////////
2077/// Fill node copy numbers of current branch into an array.
2078
2079void TGeoManager::GetBranchNumbers(Int_t *copyNumbers, Int_t *volumeNumbers) const
2080{
2081 GetCurrentNavigator()->GetBranchNumbers(copyNumbers, volumeNumbers);
2082}
2083
2084////////////////////////////////////////////////////////////////////////////////
2085/// Fill node copy numbers of current branch into an array.
2086
2088{
2090}
2091
2092////////////////////////////////////////////////////////////////////////////////
2093/// Retrieve cartesian and radial bomb factors.
2094
2095void TGeoManager::GetBombFactors(Double_t &bombx, Double_t &bomby, Double_t &bombz, Double_t &bombr) const
2096{
2097 if (fPainter) {
2098 fPainter->GetBombFactors(bombx, bomby, bombz, bombr);
2099 return;
2100 }
2101 bombx = bomby = bombz = bombr = 1.3;
2102}
2103
2104////////////////////////////////////////////////////////////////////////////////
2105/// Return maximum number of daughters of a volume used in the geometry.
2106
2108{
2109 return fgMaxDaughters;
2110}
2111
2112////////////////////////////////////////////////////////////////////////////////
2113/// Return maximum number of levels used in the geometry.
2114
2116{
2117 return fgMaxLevel;
2118}
2119
2120////////////////////////////////////////////////////////////////////////////////
2121/// Return maximum number of vertices for an xtru shape used.
2122
2124{
2125 return fgMaxXtruVert;
2126}
2127
2128////////////////////////////////////////////////////////////////////////////////
2129/// Returns number of threads that were set to use geometry.
2130
2132{
2133 return fgNumThreads;
2134}
2135
2136////////////////////////////////////////////////////////////////////////////////
2137/// Return stored current matrix (global matrix of the next touched node).
2138
2140{
2141 if (!GetCurrentNavigator())
2142 return nullptr;
2143 return GetCurrentNavigator()->GetHMatrix();
2144}
2145
2146////////////////////////////////////////////////////////////////////////////////
2147/// Returns current depth to which geometry is drawn.
2148
2150{
2151 return fVisLevel;
2152}
2153
2154////////////////////////////////////////////////////////////////////////////////
2155/// Returns current depth to which geometry is drawn.
2156
2158{
2159 return fVisOption;
2160}
2161
2162////////////////////////////////////////////////////////////////////////////////
2163/// Find level of virtuality of current overlapping node (number of levels
2164/// up having the same tracking media.
2165
2167{
2169}
2170
2171////////////////////////////////////////////////////////////////////////////////
2172/// Search the track hierarchy to find the track with the
2173/// given id
2174///
2175/// if 'primsFirst' is true, then:
2176/// first tries TGeoManager::GetTrackOfId, then does a
2177/// recursive search if that fails. this would be faster
2178/// if the track is somehow known to be a primary
2179
2181{
2182 TVirtualGeoTrack *trk = nullptr;
2183 trk = GetTrackOfId(id);
2184 if (trk)
2185 return trk;
2186 // need recursive search
2187 TIter next(fTracks);
2188 TVirtualGeoTrack *prim;
2189 while ((prim = (TVirtualGeoTrack *)next())) {
2190 trk = prim->FindTrackWithId(id);
2191 if (trk)
2192 return trk;
2193 }
2194 return nullptr;
2195}
2196
2197////////////////////////////////////////////////////////////////////////////////
2198/// Get track with a given ID.
2199
2201{
2202 TVirtualGeoTrack *track;
2203 for (Int_t i = 0; i < fNtracks; i++) {
2204 if ((track = (TVirtualGeoTrack *)fTracks->UncheckedAt(i))) {
2205 if (track->GetId() == id)
2206 return track;
2207 }
2208 }
2209 return nullptr;
2210}
2211
2212////////////////////////////////////////////////////////////////////////////////
2213/// Get parent track with a given ID.
2214
2216{
2218 while ((track = track->GetMother())) {
2219 if (track->GetId() == id)
2220 return track;
2221 }
2222 return nullptr;
2223}
2224
2225////////////////////////////////////////////////////////////////////////////////
2226/// Get index for track id, -1 if not found.
2227
2229{
2230 TVirtualGeoTrack *track;
2231 for (Int_t i = 0; i < fNtracks; i++) {
2232 if ((track = (TVirtualGeoTrack *)fTracks->UncheckedAt(i))) {
2233 if (track->GetId() == id)
2234 return i;
2235 }
2236 }
2237 return -1;
2238}
2239
2240////////////////////////////////////////////////////////////////////////////////
2241/// Go upwards the tree until a non-overlapping node
2242
2244{
2246}
2247
2248////////////////////////////////////////////////////////////////////////////////
2249/// Go upwards the tree until a non-overlapping node
2250
2252{
2254}
2255
2256////////////////////////////////////////////////////////////////////////////////
2257/// Set default volume colors according to A of material
2258
2260{
2261 const Int_t nmax = 110;
2262 Int_t col[nmax];
2263 for (Int_t i = 0; i < nmax; i++)
2264 col[i] = kGray;
2265
2266 // here we should create a new TColor with the same rgb as in the default
2267 // ROOT colors used below
2268 col[3] = kYellow - 10;
2269 col[4] = col[5] = kGreen - 10;
2270 col[6] = col[7] = kBlue - 7;
2271 col[8] = col[9] = kMagenta - 3;
2272 col[10] = col[11] = kRed - 10;
2273 col[12] = kGray + 1;
2274 col[13] = kBlue - 10;
2275 col[14] = kOrange + 7;
2276 col[16] = kYellow + 1;
2277 col[20] = kYellow - 10;
2278 col[24] = col[25] = col[26] = kBlue - 8;
2279 col[29] = kOrange + 9;
2280 col[79] = kOrange - 2;
2281
2282 TGeoVolume *vol;
2283 TIter next(fVolumes);
2284 while ((vol = (TGeoVolume *)next())) {
2285 TGeoMedium *med = vol->GetMedium();
2286 if (!med)
2287 continue;
2288 TGeoMaterial *mat = med->GetMaterial();
2289 Int_t matZ = (Int_t)mat->GetZ();
2290 vol->SetLineColor(col[matZ]);
2291 if (mat->GetDensity() < 0.1)
2292 vol->SetTransparency(60);
2293 }
2294}
2295
2296////////////////////////////////////////////////////////////////////////////////
2297/// Compute safe distance from the current point. This represent the distance
2298/// from POINT to the closest boundary.
2299
2301{
2302 return GetCurrentNavigator()->Safety(inside);
2303}
2304
2305////////////////////////////////////////////////////////////////////////////////
2306/// Set volume attributes in G3 style.
2307
2308void TGeoManager::SetVolumeAttribute(const char *name, const char *att, Int_t val)
2309{
2310 TGeoVolume *volume;
2311 Bool_t all = kFALSE;
2312 if (strstr(name, "*"))
2313 all = kTRUE;
2314 Int_t ivo = 0;
2315 TIter next(fVolumes);
2316 TString chatt = att;
2317 chatt.ToLower();
2318 while ((volume = (TGeoVolume *)next())) {
2319 if (strcmp(volume->GetName(), name) && !all)
2320 continue;
2321 ivo++;
2322 if (chatt.Contains("colo"))
2323 volume->SetLineColor(val);
2324 if (chatt.Contains("lsty"))
2325 volume->SetLineStyle(val);
2326 if (chatt.Contains("lwid"))
2327 volume->SetLineWidth(val);
2328 if (chatt.Contains("fill"))
2329 volume->SetFillColor(val);
2330 if (chatt.Contains("seen"))
2331 volume->SetVisibility(val);
2332 }
2333 TIter next1(fGVolumes);
2334 while ((volume = (TGeoVolume *)next1())) {
2335 if (strcmp(volume->GetName(), name) && !all)
2336 continue;
2337 ivo++;
2338 if (chatt.Contains("colo"))
2339 volume->SetLineColor(val);
2340 if (chatt.Contains("lsty"))
2341 volume->SetLineStyle(val);
2342 if (chatt.Contains("lwid"))
2343 volume->SetLineWidth(val);
2344 if (chatt.Contains("fill"))
2345 volume->SetFillColor(val);
2346 if (chatt.Contains("seen"))
2347 volume->SetVisibility(val);
2348 }
2349 if (!ivo) {
2350 Warning("SetVolumeAttribute", "volume: %s does not exist", name);
2351 }
2352}
2353
2354////////////////////////////////////////////////////////////////////////////////
2355/// Set factors that will "bomb" all translations in cartesian and cylindrical coordinates.
2356
2358{
2359 if (fPainter)
2360 fPainter->SetBombFactors(bombx, bomby, bombz, bombr);
2361}
2362
2363////////////////////////////////////////////////////////////////////////////////
2364/// Set a user-defined shape as clipping for ray tracing.
2365
2367{
2369 if (shape) {
2370 if (fClippingShape && (fClippingShape != shape))
2372 fClippingShape = shape;
2373 }
2374 painter->SetClippingShape(shape);
2375}
2376
2377////////////////////////////////////////////////////////////////////////////////
2378/// set the maximum number of visible nodes.
2379
2381{
2382 fMaxVisNodes = maxnodes;
2383 if (maxnodes > 0 && fgVerboseLevel > 0)
2384 Info("SetMaxVisNodes", "Automatic visible depth for %d visible nodes", maxnodes);
2385 if (!fPainter)
2386 return;
2388 Int_t level = fPainter->GetVisLevel();
2389 if (level != fVisLevel)
2390 fVisLevel = level;
2391}
2392
2393////////////////////////////////////////////////////////////////////////////////
2394/// make top volume visible on screen
2395
2397{
2399 fPainter->SetTopVisible(vis);
2400}
2401
2402////////////////////////////////////////////////////////////////////////////////
2403/// Assign a given node to be checked for overlaps. Any other overlaps will be ignored.
2404
2406{
2408}
2409
2410////////////////////////////////////////////////////////////////////////////////
2411/// Set the number of points to be generated on the shape outline when checking
2412/// for overlaps.
2413
2415{
2416 GetGeomPainter()->SetNmeshPoints(npoints);
2417}
2418
2419////////////////////////////////////////////////////////////////////////////////
2420/// set drawing mode :
2421/// - option=0 (default) all nodes drawn down to vislevel
2422/// - option=1 leaves and nodes at vislevel drawn
2423/// - option=2 path is drawn
2424/// - option=4 visibility changed
2425
2427{
2428 if ((option >= 0) && (option < 3))
2430 if (fPainter)
2432}
2433
2434////////////////////////////////////////////////////////////////////////////////
2435/// Set visualization option (leaves only OR all volumes)
2436
2438{
2439 if (flag)
2440 SetVisOption(1);
2441 else
2442 SetVisOption(0);
2443}
2444
2445////////////////////////////////////////////////////////////////////////////////
2446/// Set density threshold. Volumes with densities lower than this become
2447/// transparent.
2448
2450{
2451 fVisDensity = density;
2452 if (fPainter)
2454}
2455
2456////////////////////////////////////////////////////////////////////////////////
2457/// set default level down to which visualization is performed
2458
2460{
2461 if (level > 0) {
2462 fVisLevel = level;
2463 fMaxVisNodes = 0;
2464 if (fgVerboseLevel > 0)
2465 Info("SetVisLevel", "Automatic visible depth disabled");
2466 if (fPainter)
2468 } else {
2470 }
2471}
2472
2473////////////////////////////////////////////////////////////////////////////////
2474/// Sort overlaps by decreasing overlap distance. Extrusions comes first.
2475
2477{
2478 fOverlaps->Sort();
2479}
2480
2481////////////////////////////////////////////////////////////////////////////////
2482/// Optimize voxelization type for all volumes. Save best choice in a macro.
2483
2485{
2486 if (!fTopNode) {
2487 Error("OptimizeVoxels", "Geometry must be closed first");
2488 return;
2489 }
2490 std::ofstream out;
2491 TString fname = filename;
2492 if (fname.IsNull())
2493 fname = "tgeovox.C";
2494 out.open(fname, std::ios::out);
2495 if (!out.good()) {
2496 Error("OptimizeVoxels", "cannot open file");
2497 return;
2498 }
2499 // write header
2500 TDatime t;
2501 TString sname(fname);
2502 sname.ReplaceAll(".C", "");
2503 out << sname.Data() << "()" << std::endl;
2504 out << "{" << std::endl;
2505 out << "//=== Macro generated by ROOT version " << gROOT->GetVersion() << " : " << t.AsString() << std::endl;
2506 out << "//=== Voxel optimization for " << GetTitle() << " geometry" << std::endl;
2507 out << "//===== <run this macro JUST BEFORE closing the geometry>" << std::endl;
2508 out << " TGeoVolume *vol = 0;" << std::endl;
2509 out << " // parse all voxelized volumes" << std::endl;
2510 TGeoVolume *vol = nullptr;
2511 Bool_t cyltype;
2512 TIter next(fVolumes);
2513 while ((vol = (TGeoVolume *)next())) {
2514 if (!vol->GetVoxels())
2515 continue;
2516 out << " vol = gGeoManager->GetVolume(\"" << vol->GetName() << "\");" << std::endl;
2517 cyltype = vol->OptimizeVoxels();
2518 if (cyltype) {
2519 out << " vol->SetCylVoxels();" << std::endl;
2520 } else {
2521 out << " vol->SetCylVoxels(kFALSE);" << std::endl;
2522 }
2523 }
2524 out << "}" << std::endl;
2525 out.close();
2526}
2527////////////////////////////////////////////////////////////////////////////////
2528/// Parse a string boolean expression and do a syntax check. Find top
2529/// level boolean operator and returns its type. Fill the two
2530/// substrings to which this operator applies. The returned integer is :
2531/// - -1 : parse error
2532/// - 0 : no boolean operator
2533/// - 1 : union - represented as '+' in expression
2534/// - 2 : difference (subtraction) - represented as '-' in expression
2535/// - 3 : intersection - represented as '*' in expression.
2536/// Parentheses should be used to avoid ambiguities. For instance :
2537/// - A+B-C will be interpreted as (A+B)-C which is not the same as A+(B-C)
2538/// eliminate not needed parentheses
2539
2540Int_t TGeoManager::Parse(const char *expr, TString &expr1, TString &expr2, TString &expr3)
2541{
2542 TString startstr(expr);
2543 Int_t len = startstr.Length();
2544 Int_t i;
2545 TString e0 = "";
2546 expr3 = "";
2547 // eliminate blanks
2548 for (i = 0; i < len; i++) {
2549 if (startstr(i) == ' ')
2550 continue;
2551 e0 += startstr(i, 1);
2552 }
2553 Int_t level = 0;
2554 Int_t levmin = 999;
2555 Int_t boolop = 0;
2556 Int_t indop = 0;
2557 Int_t iloop = 1;
2558 Int_t lastop = 0;
2559 Int_t lastdp = 0;
2560 Int_t lastpp = 0;
2561 Bool_t foundmat = kFALSE;
2562 // check/eliminate parentheses
2563 while (iloop == 1) {
2564 iloop = 0;
2565 lastop = 0;
2566 lastdp = 0;
2567 lastpp = 0;
2568 len = e0.Length();
2569 for (i = 0; i < len; i++) {
2570 if (e0(i) == '(') {
2571 if (!level)
2572 iloop++;
2573 level++;
2574 continue;
2575 }
2576 if (e0(i) == ')') {
2577 level--;
2578 if (level == 0)
2579 lastpp = i;
2580 continue;
2581 }
2582 if ((e0(i) == '+') || (e0(i) == '-') || (e0(i) == '*')) {
2583 lastop = i;
2584 if (level < levmin) {
2585 levmin = level;
2586 indop = i;
2587 }
2588 continue;
2589 }
2590 if ((e0(i) == ':') && (level == 0)) {
2591 lastdp = i;
2592 continue;
2593 }
2594 }
2595 if (level != 0) {
2596 if (gGeoManager)
2597 gGeoManager->Error("Parse", "parentheses does not match");
2598 return -1;
2599 }
2600 if (iloop == 1 && (e0(0) == '(') && (e0(len - 1) == ')')) {
2601 // eliminate extra parentheses
2602 e0 = e0(1, len - 2);
2603 continue;
2604 }
2605 if (foundmat)
2606 break;
2607 if (((lastop == 0) && (lastdp > 0)) || ((lastpp > 0) && (lastdp > lastpp) && (indop < lastpp))) {
2608 expr3 = e0(lastdp + 1, len - lastdp);
2609 e0 = e0(0, lastdp);
2610 foundmat = kTRUE;
2611 iloop = 1;
2612 continue;
2613 } else
2614 break;
2615 }
2616 // loop expression and search parentheses/operators
2617 levmin = 999;
2618 for (i = 0; i < len; i++) {
2619 if (e0(i) == '(') {
2620 level++;
2621 continue;
2622 }
2623 if (e0(i) == ')') {
2624 level--;
2625 continue;
2626 }
2627 // Take LAST operator at lowest level (revision 28/07/08)
2628 if (level <= levmin) {
2629 if (e0(i) == '+') {
2630 boolop = 1; // union
2631 levmin = level;
2632 indop = i;
2633 }
2634 if (e0(i) == '-') {
2635 boolop = 2; // difference
2636 levmin = level;
2637 indop = i;
2638 }
2639 if (e0(i) == '*') {
2640 boolop = 3; // intersection
2641 levmin = level;
2642 indop = i;
2643 }
2644 }
2645 }
2646 if (indop == 0) {
2647 expr1 = e0;
2648 return indop;
2649 }
2650 expr1 = e0(0, indop);
2651 expr2 = e0(indop + 1, len - indop);
2652 return boolop;
2653}
2654
2655////////////////////////////////////////////////////////////////////////////////
2656/// Save current attributes in a macro
2657
2659{
2660 if (!fTopNode) {
2661 Error("SaveAttributes", "geometry must be closed first\n");
2662 return;
2663 }
2664 std::ofstream out;
2665 TString fname(filename);
2666 if (fname.IsNull())
2667 fname = "tgeoatt.C";
2668 out.open(fname, std::ios::out);
2669 if (!out.good()) {
2670 Error("SaveAttributes", "cannot open file");
2671 return;
2672 }
2673 // write header
2674 TDatime t;
2675 TString sname(fname);
2676 sname.ReplaceAll(".C", "");
2677 out << sname.Data() << "()" << std::endl;
2678 out << "{" << std::endl;
2679 out << "//=== Macro generated by ROOT version " << gROOT->GetVersion() << " : " << t.AsString() << std::endl;
2680 out << "//=== Attributes for " << GetTitle() << " geometry" << std::endl;
2681 out << "//===== <run this macro AFTER loading the geometry in memory>" << std::endl;
2682 // save current top volume
2683 out << " TGeoVolume *top = gGeoManager->GetVolume(\"" << fTopVolume->GetName() << "\");" << std::endl;
2684 out << " TGeoVolume *vol = 0;" << std::endl;
2685 out << " TGeoNode *node = 0;" << std::endl;
2686 out << " // clear all volume attributes and get painter" << std::endl;
2687 out << " gGeoManager->ClearAttributes();" << std::endl;
2688 out << " gGeoManager->GetGeomPainter();" << std::endl;
2689 out << " // set visualization modes and bomb factors" << std::endl;
2690 out << " gGeoManager->SetVisOption(" << GetVisOption() << ");" << std::endl;
2691 out << " gGeoManager->SetVisLevel(" << GetVisLevel() << ");" << std::endl;
2692 out << " gGeoManager->SetExplodedView(" << GetBombMode() << ");" << std::endl;
2693 Double_t bombx, bomby, bombz, bombr;
2694 GetBombFactors(bombx, bomby, bombz, bombr);
2695 out << " gGeoManager->SetBombFactors(" << bombx << "," << bomby << "," << bombz << "," << bombr << ");"
2696 << std::endl;
2697 out << " // iterate volumes container and set new attributes" << std::endl;
2698 // out << " TIter next(gGeoManager->GetListOfVolumes());"<<std::endl;
2699 TGeoVolume *vol = nullptr;
2701
2702 TIter next(fVolumes);
2703 while ((vol = (TGeoVolume *)next())) {
2704 vol->SetVisStreamed(kFALSE);
2705 }
2706 out << " // draw top volume with new settings" << std::endl;
2707 out << " top->Draw();" << std::endl;
2708 out << " gPad->x3d();" << std::endl;
2709 out << "}" << std::endl;
2710 out.close();
2711}
2712
2713////////////////////////////////////////////////////////////////////////////////
2714/// Returns the deepest node containing fPoint, which must be set a priori.
2715
2717{
2718 return GetCurrentNavigator()->SearchNode(downwards, skipnode);
2719}
2720
2721////////////////////////////////////////////////////////////////////////////////
2722/// Cross next boundary and locate within current node
2723/// The current point must be on the boundary of fCurrentNode.
2724
2726{
2727 return GetCurrentNavigator()->CrossBoundaryAndLocate(downwards, skipnode);
2728}
2729
2730////////////////////////////////////////////////////////////////////////////////
2731/// Compute distance to next boundary within STEPMAX. If no boundary is found,
2732/// propagate current point along current direction with fStep=STEPMAX. Otherwise
2733/// propagate with fStep=SNEXT (distance to boundary) and locate/return the next
2734/// node.
2735
2737{
2738 return GetCurrentNavigator()->FindNextBoundaryAndStep(stepmax, compsafe);
2739}
2740
2741////////////////////////////////////////////////////////////////////////////////
2742/// Find distance to next boundary and store it in fStep. Returns node to which this
2743/// boundary belongs. If PATH is specified, compute only distance to the node to which
2744/// PATH points. If STEPMAX is specified, compute distance only in case fSafety is smaller
2745/// than this value. STEPMAX represent the step to be made imposed by other reasons than
2746/// geometry (usually physics processes). Therefore in this case this method provides the
2747/// answer to the question : "Is STEPMAX a safe step ?" returning a NULL node and filling
2748/// fStep with a big number.
2749/// In case frombdr=kTRUE, the isotropic safety is set to zero.
2750///
2751/// Note : safety distance for the current point is computed ONLY in case STEPMAX is
2752/// specified, otherwise users have to call explicitly TGeoManager::Safety() if
2753/// they want this computed for the current point.
2754
2755TGeoNode *TGeoManager::FindNextBoundary(Double_t stepmax, const char *path, Bool_t frombdr)
2756{
2757 // convert current point and direction to local reference
2758 return GetCurrentNavigator()->FindNextBoundary(stepmax, path, frombdr);
2759}
2760
2761////////////////////////////////////////////////////////////////////////////////
2762/// Computes as fStep the distance to next daughter of the current volume.
2763/// The point and direction must be converted in the coordinate system of the current volume.
2764/// The proposed step limit is fStep.
2765
2767{
2768 return GetCurrentNavigator()->FindNextDaughterBoundary(point, dir, idaughter, compmatrix);
2769}
2770
2771////////////////////////////////////////////////////////////////////////////////
2772/// Reset current state flags.
2773
2775{
2777}
2778
2779////////////////////////////////////////////////////////////////////////////////
2780/// Returns deepest node containing current point.
2781
2783{
2784 return GetCurrentNavigator()->FindNode(safe_start);
2785}
2786
2787////////////////////////////////////////////////////////////////////////////////
2788/// Returns deepest node containing current point.
2789
2791{
2792 return GetCurrentNavigator()->FindNode(x, y, z);
2793}
2794
2795////////////////////////////////////////////////////////////////////////////////
2796/// Computes fast normal to next crossed boundary, assuming that the current point
2797/// is close enough to the boundary. Works only after calling FindNextBoundary.
2798
2800{
2802}
2803
2804////////////////////////////////////////////////////////////////////////////////
2805/// Computes normal vector to the next surface that will be or was already
2806/// crossed when propagating on a straight line from a given point/direction.
2807/// Returns the normal vector cosines in the MASTER coordinate system. The dot
2808/// product of the normal and the current direction is positive defined.
2809
2811{
2812 return GetCurrentNavigator()->FindNormal(forward);
2813}
2814
2815////////////////////////////////////////////////////////////////////////////////
2816/// Checks if point (x,y,z) is still in the current node.
2817
2819{
2820 return GetCurrentNavigator()->IsSameLocation(x, y, z, change);
2821}
2822
2823////////////////////////////////////////////////////////////////////////////////
2824/// Check if a new point with given coordinates is the same as the last located one.
2825
2827{
2828 return GetCurrentNavigator()->IsSamePoint(x, y, z);
2829}
2830
2831////////////////////////////////////////////////////////////////////////////////
2832/// True if current node is in phi range
2833
2835{
2836 if (!fPhiCut)
2837 return kTRUE;
2838 const Double_t *origin;
2840 return kFALSE;
2841 origin = ((TGeoBBox *)GetCurrentNavigator()->GetCurrentVolume()->GetShape())->GetOrigin();
2842 Double_t point[3];
2843 LocalToMaster(origin, &point[0]);
2844 Double_t phi = TMath::ATan2(point[1], point[0]) * TMath::RadToDeg();
2845 if (phi < 0)
2846 phi += 360.;
2847 if ((phi >= fPhimin) && (phi <= fPhimax))
2848 return kFALSE;
2849 return kTRUE;
2850}
2851
2852////////////////////////////////////////////////////////////////////////////////
2853/// Initialize current point and current direction vector (normalized)
2854/// in MARS. Return corresponding node.
2855
2857{
2858 return GetCurrentNavigator()->InitTrack(point, dir);
2859}
2860
2861////////////////////////////////////////////////////////////////////////////////
2862/// Initialize current point and current direction vector (normalized)
2863/// in MARS. Return corresponding node.
2864
2866{
2867 return GetCurrentNavigator()->InitTrack(x, y, z, nx, ny, nz);
2868}
2869
2870////////////////////////////////////////////////////////////////////////////////
2871/// Inspects path and all flags for the current state.
2872
2874{
2876}
2877
2878////////////////////////////////////////////////////////////////////////////////
2879/// Get path to the current node in the form /node0/node1/...
2880
2881const char *TGeoManager::GetPath() const
2882{
2883 return GetCurrentNavigator()->GetPath();
2884}
2885
2886////////////////////////////////////////////////////////////////////////////////
2887/// Get total size of geometry in bytes.
2888
2890{
2891 Int_t count = 0;
2892 TIter next(fVolumes);
2893 TGeoVolume *vol;
2894 while ((vol = (TGeoVolume *)next()))
2895 count += vol->GetByteCount();
2896 TIter next1(fMatrices);
2897 TGeoMatrix *matrix;
2898 while ((matrix = (TGeoMatrix *)next1()))
2899 count += matrix->GetByteCount();
2900 TIter next2(fMaterials);
2901 TGeoMaterial *mat;
2902 while ((mat = (TGeoMaterial *)next2()))
2903 count += mat->GetByteCount();
2904 TIter next3(fMedia);
2905 TGeoMedium *med;
2906 while ((med = (TGeoMedium *)next3()))
2907 count += med->GetByteCount();
2908 if (fgVerboseLevel > 0)
2909 Info("GetByteCount", "Total size of logical tree : %i bytes", count);
2910 return count;
2911}
2912
2913////////////////////////////////////////////////////////////////////////////////
2914/// Make a default painter if none present. Returns pointer to it.
2915
2917{
2918 if (!fPainter) {
2919 const char *kind = "root";
2920 if (gROOT->IsWebDisplay() && !gROOT->IsWebDisplayBatch())
2921 kind = "web";
2922 if (auto h = gROOT->GetPluginManager()->FindHandler("TVirtualGeoPainter", kind)) {
2923 if (h->LoadPlugin() == -1) {
2924 Error("GetGeomPainter", "could not load plugin for %s geo_painter", kind);
2925 return nullptr;
2926 }
2927 fPainter = (TVirtualGeoPainter *)h->ExecPlugin(1, this);
2928 if (!fPainter) {
2929 Error("GetGeomPainter", "could not create %s geo_painter", kind);
2930 return nullptr;
2931 }
2932 } else {
2933 Error("GetGeomPainter", "not found plugin %s for geo_painter", kind);
2934 }
2935 }
2936 return fPainter;
2937}
2938
2939////////////////////////////////////////////////////////////////////////////////
2940/// Search for a named volume. All trailing blanks stripped.
2941
2943{
2944 TString sname = name;
2945 sname = sname.Strip();
2946 TGeoVolume *vol = (TGeoVolume *)fVolumes->FindObject(sname.Data());
2947 return vol;
2948}
2949
2950////////////////////////////////////////////////////////////////////////////////
2951/// Fast search for a named volume. All trailing blanks stripped.
2952
2954{
2955 if (!fHashVolumes) {
2956 Int_t nvol = fVolumes->GetEntriesFast();
2957 Int_t ngvol = fGVolumes->GetEntriesFast();
2958 fHashVolumes = new THashList(nvol + 1);
2959 fHashGVolumes = new THashList(ngvol + 1);
2960 Int_t i;
2961 for (i = 0; i < ngvol; i++)
2963 for (i = 0; i < nvol; i++)
2965 }
2966 TString sname = name;
2967 sname = sname.Strip();
2968 THashList *list = fHashVolumes;
2969 if (multi)
2970 list = fHashGVolumes;
2971 TGeoVolume *vol = (TGeoVolume *)list->FindObject(sname.Data());
2972 return vol;
2973}
2974
2975////////////////////////////////////////////////////////////////////////////////
2976/// Retrieve unique id for a volume name. Return -1 if name not found.
2977
2978Int_t TGeoManager::GetUID(const char *volname) const
2979{
2980 TGeoManager *geom = (TGeoManager *)this;
2981 TGeoVolume *vol = geom->FindVolumeFast(volname, kFALSE);
2982 if (!vol)
2983 vol = geom->FindVolumeFast(volname, kTRUE);
2984 if (!vol)
2985 return -1;
2986 return vol->GetNumber();
2987}
2988
2989////////////////////////////////////////////////////////////////////////////////
2990/// Find if a given material duplicates an existing one.
2991
2993{
2994 Int_t index = fMaterials->IndexOf(mat);
2995 if (index <= 0)
2996 return nullptr;
2997 TGeoMaterial *other;
2998 for (Int_t i = 0; i < index; i++) {
2999 other = (TGeoMaterial *)fMaterials->At(i);
3000 if (other == mat)
3001 continue;
3002 if (other->IsEq(mat))
3003 return other;
3004 }
3005 return nullptr;
3006}
3007
3008////////////////////////////////////////////////////////////////////////////////
3009/// Search for a named material. All trailing blanks stripped.
3010
3011TGeoMaterial *TGeoManager::GetMaterial(const char *matname) const
3012{
3013 TString sname = matname;
3014 sname = sname.Strip();
3016 return mat;
3017}
3018
3019////////////////////////////////////////////////////////////////////////////////
3020/// Search for a named tracking medium. All trailing blanks stripped.
3021
3022TGeoMedium *TGeoManager::GetMedium(const char *medium) const
3023{
3024 TString sname = medium;
3025 sname = sname.Strip();
3026 TGeoMedium *med = (TGeoMedium *)fMedia->FindObject(sname.Data());
3027 return med;
3028}
3029
3030////////////////////////////////////////////////////////////////////////////////
3031/// Search for a tracking medium with a given ID.
3032
3034{
3035 TIter next(fMedia);
3036 TGeoMedium *med;
3037 while ((med = (TGeoMedium *)next())) {
3038 if (med->GetId() == numed)
3039 return med;
3040 }
3041 return nullptr;
3042}
3043
3044////////////////////////////////////////////////////////////////////////////////
3045/// Return material at position id.
3046
3048{
3049 if (id < 0 || id >= fMaterials->GetSize())
3050 return nullptr;
3051 TGeoMaterial *mat = (TGeoMaterial *)fMaterials->At(id);
3052 return mat;
3053}
3054
3055////////////////////////////////////////////////////////////////////////////////
3056/// Return index of named material.
3057
3058Int_t TGeoManager::GetMaterialIndex(const char *matname) const
3059{
3060 TIter next(fMaterials);
3061 TGeoMaterial *mat;
3062 Int_t id = 0;
3063 TString sname = matname;
3064 sname = sname.Strip();
3065 while ((mat = (TGeoMaterial *)next())) {
3066 if (!strcmp(mat->GetName(), sname.Data()))
3067 return id;
3068 id++;
3069 }
3070 return -1; // fail
3071}
3072
3073////////////////////////////////////////////////////////////////////////////////
3074/// Randomly shoot nrays and plot intersections with surfaces for current
3075/// top node.
3076
3077void TGeoManager::RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol,
3078 Bool_t check_norm)
3079{
3080 GetGeomPainter()->RandomRays(nrays, startx, starty, startz, target_vol, check_norm);
3081}
3082
3083////////////////////////////////////////////////////////////////////////////////
3084/// Remove material at given index.
3085
3087{
3088 TObject *obj = fMaterials->At(index);
3089 if (obj)
3090 fMaterials->Remove(obj);
3091}
3092
3093////////////////////////////////////////////////////////////////////////////////
3094/// Sets all pointers TGeoVolume::fField to NULL. User data becomes decoupled
3095/// from geometry. Deletion has to be managed by users.
3096
3098{
3099 TIter next(fVolumes);
3100 TGeoVolume *vol;
3101 while ((vol = (TGeoVolume *)next()))
3102 vol->SetField(nullptr);
3103}
3104
3105////////////////////////////////////////////////////////////////////////////////
3106/// Change raytracing mode.
3107
3109{
3111 if (fPainter && fPainter->IsRaytracing())
3112 ModifiedPad();
3113}
3114
3115////////////////////////////////////////////////////////////////////////////////
3116/// Restore the master volume of the geometry.
3117
3119{
3121 return;
3122 if (fMasterVolume)
3124}
3125
3126////////////////////////////////////////////////////////////////////////////////
3127/// Voxelize all non-divided volumes.
3128
3130{
3131 TGeoVolume *vol;
3132 // TGeoVoxelFinder *vox = 0;
3133 if (!fStreamVoxels && fgVerboseLevel > 0)
3134 Info("Voxelize", "Voxelizing...");
3135 // Int_t nentries = fVolumes->GetSize();
3136 TIter next(fVolumes);
3137 while ((vol = (TGeoVolume *)next())) {
3138 if (!fIsGeomReading)
3139 vol->SortNodes();
3140 if (!fStreamVoxels) {
3141 vol->Voxelize(option);
3142 }
3143 if (!fIsGeomReading)
3144 vol->FindOverlaps();
3145 }
3146}
3147
3148////////////////////////////////////////////////////////////////////////////////
3149/// Send "Modified" signal to painter.
3150
3152{
3153 if (!fPainter)
3154 return;
3156}
3157
3158////////////////////////////////////////////////////////////////////////////////
3159/// Make an TGeoArb8 volume.
3160
3162{
3163 return TGeoBuilder::Instance(this)->MakeArb8(name, medium, dz, vertices);
3164}
3165
3166////////////////////////////////////////////////////////////////////////////////
3167/// Make in one step a volume pointing to a box shape with given medium.
3168
3170{
3171 return TGeoBuilder::Instance(this)->MakeBox(name, medium, dx, dy, dz);
3172}
3173
3174////////////////////////////////////////////////////////////////////////////////
3175/// Make in one step a volume pointing to a parallelepiped shape with given medium.
3176
3178 Double_t alpha, Double_t theta, Double_t phi)
3179{
3180 return TGeoBuilder::Instance(this)->MakePara(name, medium, dx, dy, dz, alpha, theta, phi);
3181}
3182
3183////////////////////////////////////////////////////////////////////////////////
3184/// Make in one step a volume pointing to a sphere shape with given medium
3185
3187 Double_t themax, Double_t phimin, Double_t phimax)
3188{
3189 return TGeoBuilder::Instance(this)->MakeSphere(name, medium, rmin, rmax, themin, themax, phimin, phimax);
3190}
3191
3192////////////////////////////////////////////////////////////////////////////////
3193/// Make in one step a volume pointing to a torus shape with given medium.
3194
3196 Double_t phi1, Double_t dphi)
3197{
3198 return TGeoBuilder::Instance(this)->MakeTorus(name, medium, r, rmin, rmax, phi1, dphi);
3199}
3200
3201////////////////////////////////////////////////////////////////////////////////
3202/// Make in one step a volume pointing to a tube shape with given medium.
3203
3205{
3206 return TGeoBuilder::Instance(this)->MakeTube(name, medium, rmin, rmax, dz);
3207}
3208
3209////////////////////////////////////////////////////////////////////////////////
3210/// Make in one step a volume pointing to a tube segment shape with given medium.
3211/// The segment will be from phiStart to phiEnd, the angles are expressed in degree
3212
3214 Double_t phiStart, Double_t phiEnd)
3215{
3216 return TGeoBuilder::Instance(this)->MakeTubs(name, medium, rmin, rmax, dz, phiStart, phiEnd);
3217}
3218
3219////////////////////////////////////////////////////////////////////////////////
3220/// Make in one step a volume pointing to a tube shape with given medium
3221
3223{
3224 return TGeoBuilder::Instance(this)->MakeEltu(name, medium, a, b, dz);
3225}
3226
3227////////////////////////////////////////////////////////////////////////////////
3228/// Make in one step a volume pointing to a tube shape with given medium
3229
3231 Double_t stout, Double_t dz)
3232{
3233 return TGeoBuilder::Instance(this)->MakeHype(name, medium, rin, stin, rout, stout, dz);
3234}
3235
3236////////////////////////////////////////////////////////////////////////////////
3237/// Make in one step a volume pointing to a tube shape with given medium
3238
3240{
3241 return TGeoBuilder::Instance(this)->MakeParaboloid(name, medium, rlo, rhi, dz);
3242}
3243
3244////////////////////////////////////////////////////////////////////////////////
3245/// Make in one step a volume pointing to a tube segment shape with given medium
3246
3248 Double_t phi1, Double_t phi2, Double_t lx, Double_t ly, Double_t lz, Double_t tx,
3249 Double_t ty, Double_t tz)
3250{
3251 return TGeoBuilder::Instance(this)->MakeCtub(name, medium, rmin, rmax, dz, phi1, phi2, lx, ly, lz, tx, ty, tz);
3252}
3253
3254////////////////////////////////////////////////////////////////////////////////
3255/// Make in one step a volume pointing to a cone shape with given medium.
3256
3258 Double_t rmin2, Double_t rmax2)
3259{
3260 return TGeoBuilder::Instance(this)->MakeCone(name, medium, dz, rmin1, rmax1, rmin2, rmax2);
3261}
3262
3263////////////////////////////////////////////////////////////////////////////////
3264/// Make in one step a volume pointing to a cone segment shape with given medium
3265
3267 Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
3268{
3269 return TGeoBuilder::Instance(this)->MakeCons(name, medium, dz, rmin1, rmax1, rmin2, rmax2, phi1, phi2);
3270}
3271
3272////////////////////////////////////////////////////////////////////////////////
3273/// Make in one step a volume pointing to a polycone shape with given medium.
3274
3276{
3277 return TGeoBuilder::Instance(this)->MakePcon(name, medium, phi, dphi, nz);
3278}
3279
3280////////////////////////////////////////////////////////////////////////////////
3281/// Make in one step a volume pointing to a polygone shape with given medium.
3282
3283TGeoVolume *
3284TGeoManager::MakePgon(const char *name, TGeoMedium *medium, Double_t phi, Double_t dphi, Int_t nedges, Int_t nz)
3285{
3286 return TGeoBuilder::Instance(this)->MakePgon(name, medium, phi, dphi, nedges, nz);
3287}
3288
3289////////////////////////////////////////////////////////////////////////////////
3290/// Make in one step a volume pointing to a TGeoTrd1 shape with given medium.
3291
3292TGeoVolume *
3294{
3295 return TGeoBuilder::Instance(this)->MakeTrd1(name, medium, dx1, dx2, dy, dz);
3296}
3297
3298////////////////////////////////////////////////////////////////////////////////
3299/// Make in one step a volume pointing to a TGeoTrd2 shape with given medium.
3300
3302 Double_t dy2, Double_t dz)
3303{
3304 return TGeoBuilder::Instance(this)->MakeTrd2(name, medium, dx1, dx2, dy1, dy2, dz);
3305}
3306
3307////////////////////////////////////////////////////////////////////////////////
3308/// Make in one step a volume pointing to a trapezoid shape with given medium.
3309
3311 Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
3312 Double_t tl2, Double_t alpha2)
3313{
3314 return TGeoBuilder::Instance(this)->MakeTrap(name, medium, dz, theta, phi, h1, bl1, tl1, alpha1, h2, bl2, tl2,
3315 alpha2);
3316}
3317
3318////////////////////////////////////////////////////////////////////////////////
3319/// Make in one step a volume pointing to a twisted trapezoid shape with given medium.
3320
3322 Double_t twist, Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2,
3323 Double_t bl2, Double_t tl2, Double_t alpha2)
3324{
3325 return TGeoBuilder::Instance(this)->MakeGtra(name, medium, dz, theta, phi, twist, h1, bl1, tl1, alpha1, h2, bl2, tl2,
3326 alpha2);
3327}
3328
3329////////////////////////////////////////////////////////////////////////////////
3330/// Make a TGeoXtru-shaped volume with nz planes
3331
3333{
3334 return TGeoBuilder::Instance(this)->MakeXtru(name, medium, nz);
3335}
3336
3337////////////////////////////////////////////////////////////////////////////////
3338/// Creates an alignable object with unique name corresponding to a path
3339/// and adds it to the list of alignables. An optional unique ID can be
3340/// provided, in which case PN entries can be searched fast by uid.
3341
3342TGeoPNEntry *TGeoManager::SetAlignableEntry(const char *unique_name, const char *path, Int_t uid)
3343{
3344 if (!CheckPath(path))
3345 return nullptr;
3346 if (!fHashPNE)
3347 fHashPNE = new THashList(256, 3);
3348 if (!fArrayPNE)
3349 fArrayPNE = new TObjArray(256);
3350 TGeoPNEntry *entry = GetAlignableEntry(unique_name);
3351 if (entry) {
3352 Error("SetAlignableEntry", "An alignable object with name %s already existing. NOT ADDED !", unique_name);
3353 return nullptr;
3354 }
3355 entry = new TGeoPNEntry(unique_name, path);
3356 Int_t ientry = fHashPNE->GetSize();
3357 fHashPNE->Add(entry);
3358 fArrayPNE->AddAtAndExpand(entry, ientry);
3359 if (uid >= 0) {
3360 Bool_t added = InsertPNEId(uid, ientry);
3361 if (!added)
3362 Error("SetAlignableEntry", "A PN entry: has already uid=%i", uid);
3363 }
3364 return entry;
3365}
3366
3367////////////////////////////////////////////////////////////////////////////////
3368/// Retrieves an existing alignable object.
3369
3371{
3372 if (!fHashPNE)
3373 return nullptr;
3374 return (TGeoPNEntry *)fHashPNE->FindObject(name);
3375}
3376
3377////////////////////////////////////////////////////////////////////////////////
3378/// Retrieves an existing alignable object at a given index.
3379
3381{
3382 if (!fArrayPNE && !InitArrayPNE())
3383 return nullptr;
3384 return (TGeoPNEntry *)fArrayPNE->At(index);
3385}
3386
3387////////////////////////////////////////////////////////////////////////////////
3388/// Retrieves an existing alignable object having a preset UID.
3389
3391{
3392 if (!fNPNEId || (!fArrayPNE && !InitArrayPNE()))
3393 return nullptr;
3395 if (index < 0 || fKeyPNEId[index] != uid)
3396 return nullptr;
3398}
3399
3400////////////////////////////////////////////////////////////////////////////////
3401/// Retrieves number of PN entries with or without UID.
3402
3404{
3405 if (!fHashPNE)
3406 return 0;
3407 if (with_uid)
3408 return fNPNEId;
3409 return fHashPNE->GetSize();
3410}
3411
3412////////////////////////////////////////////////////////////////////////////////
3413/// Insert a PN entry in the sorted array of indexes.
3414
3416{
3417 if (!fSizePNEId) {
3418 // Create the arrays.
3419 fSizePNEId = 128;
3420 fKeyPNEId = new Int_t[fSizePNEId];
3421 memset(fKeyPNEId, 0, fSizePNEId * sizeof(Int_t));
3423 memset(fValuePNEId, 0, fSizePNEId * sizeof(Int_t));
3424 fKeyPNEId[fNPNEId] = uid;
3425 fValuePNEId[fNPNEId++] = ientry;
3426 return kTRUE;
3427 }
3428 // Search id in the existing array and return false if it already exists.
3430 if (index > 0 && fKeyPNEId[index] == uid)
3431 return kFALSE;
3432 // Resize the arrays and insert the value
3433 Bool_t resize = (fNPNEId == fSizePNEId) ? kTRUE : kFALSE;
3434 if (resize) {
3435 // Double the size of the array
3436 fSizePNEId *= 2;
3437 // Create new arrays of keys and values
3438 Int_t *keys = new Int_t[fSizePNEId];
3439 memset(keys, 0, fSizePNEId * sizeof(Int_t));
3440 Int_t *values = new Int_t[fSizePNEId];
3441 memset(values, 0, fSizePNEId * sizeof(Int_t));
3442 // Copy all keys<uid in the new keys array (0 to index)
3443 memcpy(keys, fKeyPNEId, (index + 1) * sizeof(Int_t));
3444 memcpy(values, fValuePNEId, (index + 1) * sizeof(Int_t));
3445 // Insert current key at index+1
3446 keys[index + 1] = uid;
3447 values[index + 1] = ientry;
3448 // Copy all remaining keys from the old to new array
3449 memcpy(&keys[index + 2], &fKeyPNEId[index + 1], (fNPNEId - index - 1) * sizeof(Int_t));
3450 memcpy(&values[index + 2], &fValuePNEId[index + 1], (fNPNEId - index - 1) * sizeof(Int_t));
3451 delete[] fKeyPNEId;
3452 fKeyPNEId = keys;
3453 delete[] fValuePNEId;
3454 fValuePNEId = values;
3455 fNPNEId++;
3456 return kTRUE;
3457 }
3458 // Insert the value in the existing arrays
3459 Int_t i;
3460 for (i = fNPNEId - 1; i > index; i--) {
3461 fKeyPNEId[i + 1] = fKeyPNEId[i];
3462 fValuePNEId[i + 1] = fValuePNEId[i];
3463 }
3464 fKeyPNEId[index + 1] = uid;
3465 fValuePNEId[index + 1] = ientry;
3466 fNPNEId++;
3467 return kTRUE;
3468}
3469
3470////////////////////////////////////////////////////////////////////////////////
3471/// Make a physical node from the path pointed by an alignable object with a given name.
3472
3474{
3476 if (!entry) {
3477 Error("MakeAlignablePN", "No alignable object named %s found !", name);
3478 return nullptr;
3479 }
3480 return MakeAlignablePN(entry);
3481}
3482
3483////////////////////////////////////////////////////////////////////////////////
3484/// Make a physical node from the path pointed by a given alignable object.
3485
3487{
3488 if (!entry) {
3489 Error("MakeAlignablePN", "No alignable object specified !");
3490 return nullptr;
3491 }
3492 const char *path = entry->GetTitle();
3493 if (!cd(path)) {
3494 Error("MakeAlignablePN", "Alignable object %s poins to invalid path: %s", entry->GetName(), path);
3495 return nullptr;
3496 }
3497 TGeoPhysicalNode *node = MakePhysicalNode(path);
3498 entry->SetPhysicalNode(node);
3499 return node;
3500}
3501
3502////////////////////////////////////////////////////////////////////////////////
3503/// Makes a physical node corresponding to a path. If PATH is not specified,
3504/// makes physical node matching current modeller state.
3505
3507{
3508 TGeoPhysicalNode *node;
3509 if (path) {
3510 if (!CheckPath(path)) {
3511 Error("MakePhysicalNode", "path: %s not valid", path);
3512 return nullptr;
3513 }
3514 node = new TGeoPhysicalNode(path);
3515 } else {
3516 node = new TGeoPhysicalNode(GetPath());
3517 }
3518 fPhysicalNodes->Add(node);
3519 return node;
3520}
3521
3522////////////////////////////////////////////////////////////////////////////////
3523/// Refresh physical nodes to reflect the actual geometry paths after alignment
3524/// was applied. Optionally locks physical nodes (default).
3525
3527{
3529 TGeoPhysicalNode *pn;
3530 while ((pn = (TGeoPhysicalNode *)next()))
3531 pn->Refresh();
3534 if (lock)
3535 LockGeometry();
3536}
3537
3538////////////////////////////////////////////////////////////////////////////////
3539/// Clear the current list of physical nodes, so that we can start over with a new list.
3540/// If MUSTDELETE is true, delete previous nodes.
3541
3543{
3544 if (mustdelete)
3546 else
3548}
3549
3550////////////////////////////////////////////////////////////////////////////////
3551/// Make an assembly of volumes.
3552
3554{
3556}
3557
3558////////////////////////////////////////////////////////////////////////////////
3559/// Make a TGeoVolumeMulti handling a list of volumes.
3560
3562{
3563 return TGeoBuilder::Instance(this)->MakeVolumeMulti(name, medium);
3564}
3565
3566////////////////////////////////////////////////////////////////////////////////
3567/// Set type of exploding view (see TGeoPainter::SetExplodedView())
3568
3570{
3571 if ((ibomb >= 0) && (ibomb < 4))
3572 fExplodedView = ibomb;
3573 if (fPainter)
3574 fPainter->SetExplodedView(ibomb);
3575}
3576
3577////////////////////////////////////////////////////////////////////////////////
3578/// Set cut phi range
3579
3581{
3582 if ((phimin == 0) && (phimax == 360)) {
3583 fPhiCut = kFALSE;
3584 return;
3585 }
3586 fPhiCut = kTRUE;
3587 fPhimin = phimin;
3588 fPhimax = phimax;
3589}
3590
3591////////////////////////////////////////////////////////////////////////////////
3592/// Set number of segments for approximating circles in drawing.
3593
3595{
3596 if (fNsegments == nseg)
3597 return;
3598 if (nseg > 2)
3599 fNsegments = nseg;
3600 if (fPainter)
3601 fPainter->SetNsegments(nseg);
3602}
3603
3604////////////////////////////////////////////////////////////////////////////////
3605/// Get number of segments approximating circles
3606
3608{
3609 return fNsegments;
3610}
3611
3612////////////////////////////////////////////////////////////////////////////////
3613/// Now just a shortcut for GetElementTable.
3614
3616{
3619}
3620
3621////////////////////////////////////////////////////////////////////////////////
3622/// Returns material table. Creates it if not existing.
3623
3625{
3626 if (!fElementTable)
3628 return fElementTable;
3629}
3630
3631////////////////////////////////////////////////////////////////////////////////
3632/// Make a rectilinear step of length fStep from current point (fPoint) on current
3633/// direction (fDirection). If the step is imposed by geometry, is_geom flag
3634/// must be true (default). The cross flag specifies if the boundary should be
3635/// crossed in case of a geometry step (default true). Returns new node after step.
3636/// Set also on boundary condition.
3637
3639{
3640 return GetCurrentNavigator()->Step(is_geom, cross);
3641}
3642
3643////////////////////////////////////////////////////////////////////////////////
3644/// shoot npoints randomly in a box of 1E-5 around current point.
3645/// return minimum distance to points outside
3646
3647TGeoNode *TGeoManager::SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path)
3648{
3649 return GetGeomPainter()->SamplePoints(npoints, dist, epsil, g3path);
3650}
3651
3652////////////////////////////////////////////////////////////////////////////////
3653/// Set the top volume and corresponding node as starting point of the geometry.
3654
3656{
3657 if (fTopVolume == vol)
3658 return;
3659
3660 TSeqCollection *brlist = gROOT->GetListOfBrowsers();
3661 TIter next(brlist);
3662 TBrowser *browser = nullptr;
3663
3664 if (fTopVolume)
3665 fTopVolume->SetTitle("");
3666 fTopVolume = vol;
3667 vol->SetTitle("Top volume");
3668 if (fTopNode) {
3669 TGeoNode *topn = fTopNode;
3670 fTopNode = nullptr;
3671 while ((browser = (TBrowser *)next()))
3672 browser->RecursiveRemove(topn);
3673 delete topn;
3674 } else {
3675 fMasterVolume = vol;
3678 if (fgVerboseLevel > 0)
3679 Info("SetTopVolume", "Top volume is %s. Master volume is %s", fTopVolume->GetName(), fMasterVolume->GetName());
3680 }
3681 // fMasterVolume->FindMatrixOfDaughterVolume(vol);
3682 // fCurrentMatrix->Print();
3684 fTopNode->SetName(TString::Format("%s_1", vol->GetName()));
3685 fTopNode->SetNumber(1);
3686 fTopNode->SetTitle("Top logical node");
3687 fNodes->AddAt(fTopNode, 0);
3688 if (!GetCurrentNavigator()) {
3690 return;
3691 }
3692 Int_t nnavigators = 0;
3694 if (!arr)
3695 return;
3696 nnavigators = arr->GetEntriesFast();
3697 for (Int_t i = 0; i < nnavigators; i++) {
3698 TGeoNavigator *nav = (TGeoNavigator *)arr->At(i);
3699 nav->ResetAll();
3700 if (fClosed)
3701 nav->GetCache()->BuildInfoBranch();
3702 }
3703}
3704
3705////////////////////////////////////////////////////////////////////////////////
3706/// Define different tracking media.
3707
3709{
3710 /*
3711 Int_t nmat = fMaterials->GetSize();
3712 if (!nmat) {printf(" No materials !\n"); return;}
3713 Int_t *media = new Int_t[nmat];
3714 memset(media, 0, nmat*sizeof(Int_t));
3715 Int_t imedia = 1;
3716 TGeoMaterial *mat, *matref;
3717 mat = (TGeoMaterial*)fMaterials->At(0);
3718 if (mat->GetMedia()) {
3719 for (Int_t i=0; i<nmat; i++) {
3720 mat = (TGeoMaterial*)fMaterials->At(i);
3721 mat->Print();
3722 }
3723 return;
3724 }
3725 mat->SetMedia(imedia);
3726 media[0] = imedia++;
3727 mat->Print();
3728 for (Int_t i=0; i<nmat; i++) {
3729 mat = (TGeoMaterial*)fMaterials->At(i);
3730 for (Int_t j=0; j<i; j++) {
3731 matref = (TGeoMaterial*)fMaterials->At(j);
3732 if (mat->IsEq(matref)) {
3733 mat->SetMedia(media[j]);
3734 break;
3735 }
3736 if (j==(i-1)) {
3737 // different material
3738 mat->SetMedia(imedia);
3739 media[i] = imedia++;
3740 mat->Print();
3741 }
3742 }
3743 }
3744 */
3745}
3746
3747////////////////////////////////////////////////////////////////////////////////
3748/// Check pushes and pulls needed to cross the next boundary with respect to the
3749/// position given by FindNextBoundary. If radius is not mentioned the full bounding
3750/// box will be sampled.
3751
3753{
3754 GetGeomPainter()->CheckBoundaryErrors(ntracks, radius);
3755}
3756
3757////////////////////////////////////////////////////////////////////////////////
3758/// Check the boundary errors reference file created by CheckBoundaryErrors method.
3759/// The shape for which the crossing failed is drawn with the starting point in red
3760/// and the extrapolated point to boundary (+/- failing push/pull) in yellow.
3761
3763{
3765}
3766
3767////////////////////////////////////////////////////////////////////////////////
3768/// Classify a given point. See TGeoChecker::CheckPoint().
3769
3771{
3772 GetGeomPainter()->CheckPoint(x, y, z, option, safety);
3773}
3774
3775////////////////////////////////////////////////////////////////////////////////
3776/// Test for shape navigation methods. Summary for test numbers:
3777/// - 1: DistFromInside/Outside. Sample points inside the shape. Generate
3778/// directions randomly in cos(theta). Compute DistFromInside and move the
3779/// point with bigger distance. Compute DistFromOutside back from new point.
3780/// Plot d-(d1+d2)
3781///
3782
3784{
3785 GetGeomPainter()->CheckShape(shape, testNo, nsamples, option);
3786}
3787
3788////////////////////////////////////////////////////////////////////////////////
3789/// Geometry checking.
3790/// - if option contains 'o': Optional overlap checkings (by sampling and by mesh).
3791/// - if option contains 'b': Optional boundary crossing check + timing per volume.
3792///
3793/// STAGE 1: extensive overlap checking by sampling per volume. Stdout need to be
3794/// checked by user to get report, then TGeoVolume::CheckOverlaps(0.01, "s") can
3795/// be called for the suspicious volumes.
3796///
3797/// STAGE 2: normal overlap checking using the shapes mesh - fills the list of
3798/// overlaps.
3799///
3800/// STAGE 3: shooting NRAYS rays from VERTEX and counting the total number of
3801/// crossings per volume (rays propagated from boundary to boundary until
3802/// geometry exit). Timing computed and results stored in a histo.
3803///
3804/// STAGE 4: shooting 1 mil. random rays inside EACH volume and calling
3805/// FindNextBoundary() + Safety() for each call. The timing is normalized by the
3806/// number of crossings computed at stage 2 and presented as percentage.
3807/// One can get a picture on which are the most "burned" volumes during
3808/// transportation from geometry point of view. Another plot of the timing per
3809/// volume vs. number of daughters is produced.
3810
3812{
3813 TString opt(option);
3814 opt.ToLower();
3815 if (!opt.Length()) {
3816 Error("CheckGeometryFull", "The option string must contain a letter. See method documentation.");
3817 return;
3818 }
3819 Bool_t checkoverlaps = opt.Contains("o");
3820 Bool_t checkcrossings = opt.Contains("b");
3821 Double_t vertex[3];
3822 vertex[0] = vx;
3823 vertex[1] = vy;
3824 vertex[2] = vz;
3825 GetGeomPainter()->CheckGeometryFull(checkoverlaps, checkcrossings, ntracks, vertex);
3826}
3827
3828////////////////////////////////////////////////////////////////////////////////
3829/// Perform last checks on the geometry
3830
3832{
3833 if (fgVerboseLevel > 0)
3834 Info("CheckGeometry", "Fixing runtime shapes...");
3835 TIter next(fShapes);
3836 TIter nextv(fVolumes);
3837 TGeoShape *shape;
3838 TGeoVolume *vol;
3839 Bool_t has_runtime = kFALSE;
3840 while ((shape = (TGeoShape *)next())) {
3841 if (shape->IsRunTimeShape()) {
3842 has_runtime = kTRUE;
3843 }
3844 if (fIsGeomReading)
3845 shape->AfterStreamer();
3848 shape->ComputeBBox();
3849 }
3850 if (has_runtime)
3852 else if (fgVerboseLevel > 0)
3853 Info("CheckGeometry", "...Nothing to fix");
3854 // Compute bounding box for assemblies
3856 while ((vol = (TGeoVolume *)nextv())) {
3857 if (vol->IsAssembly())
3858 vol->GetShape()->ComputeBBox();
3859 else if (vol->GetMedium() == dummy) {
3860 Warning("CheckGeometry", "Volume \"%s\" has no medium: assigned dummy medium and material", vol->GetName());
3861 vol->SetMedium(dummy);
3862 }
3863 }
3864}
3865
3866////////////////////////////////////////////////////////////////////////////////
3867/// Check all geometry for illegal overlaps within a limit OVLP.
3868
3870{
3871 if (!fTopNode) {
3872 Error("CheckOverlaps", "Top node not set");
3873 return;
3874 }
3876}
3877
3878////////////////////////////////////////////////////////////////////////////////
3879/// Prints the current list of overlaps.
3880
3882{
3883 if (!fOverlaps)
3884 return;
3885 Int_t novlp = fOverlaps->GetEntriesFast();
3886 if (!novlp)
3887 return;
3888 TGeoManager *geom = (TGeoManager *)this;
3889 geom->GetGeomPainter()->PrintOverlaps();
3890}
3891
3892////////////////////////////////////////////////////////////////////////////////
3893/// Estimate weight of volume VOL with a precision SIGMA(W)/W better than PRECISION.
3894/// Option can be "v" - verbose (default)
3895
3897{
3899 TString opt(option);
3900 opt.ToLower();
3901 Double_t weight;
3902 TGeoVolume *volume = fTopVolume;
3903 if (opt.Contains("v")) {
3904 if (opt.Contains("a")) {
3905 if (fgVerboseLevel > 0)
3906 Info("Weight", "Computing analytically weight of %s", volume->GetName());
3907 weight = volume->WeightA();
3908 if (fgVerboseLevel > 0)
3909 Info("Weight", "Computed weight: %f [kg]\n", weight);
3910 return weight;
3911 }
3912 if (fgVerboseLevel > 0) {
3913 Info("Weight", "Estimating weight of %s with %g %% precision", fTopVolume->GetName(), 100. * precision);
3914 printf(" event weight err\n");
3915 printf("========================================\n");
3916 }
3917 }
3918 weight = fPainter->Weight(precision, option);
3919 return weight;
3920}
3921
3922////////////////////////////////////////////////////////////////////////////////
3923/// computes the total size in bytes of the branch starting with node.
3924/// The option can specify if all the branch has to be parsed or only the node
3925
3926ULong_t TGeoManager::SizeOf(const TGeoNode * /*node*/, Option_t * /*option*/)
3927{
3928 return 0;
3929}
3930
3931////////////////////////////////////////////////////////////////////////////////
3932/// Stream an object of class TGeoManager.
3933
3935{
3936 if (R__b.IsReading()) {
3937 R__b.ReadClassBuffer(TGeoManager::Class(), this);
3939 CloseGeometry();
3942 } else {
3944 }
3945}
3946
3947////////////////////////////////////////////////////////////////////////////////
3948/// Execute mouse actions on this manager.
3949
3951{
3952 if (!fPainter)
3953 return;
3954 fPainter->ExecuteManagerEvent(this, event, px, py);
3955}
3956
3957////////////////////////////////////////////////////////////////////////////////
3958/// Export this geometry to a file
3959///
3960/// - Case 1: root file or root/xml file
3961/// if filename end with ".root". The key will be named name
3962/// By default the geometry is saved without the voxelisation info.
3963/// Use option 'v" to save the voxelisation info.
3964/// if filename end with ".xml" a root/xml file is produced.
3965///
3966/// - Case 2: C++ script
3967/// if filename end with ".C"
3968///
3969/// - Case 3: gdml file
3970/// if filename end with ".gdml"
3971/// NOTE that to use this option, the PYTHONPATH must be defined like
3972/// export PYTHONPATH=$ROOTSYS/lib:$ROOTSYS/geom/gdml
3973///
3974
3976{
3977 TString sfile(filename);
3978 if (sfile.Contains(".C")) {
3979 // Save geometry as a C++ script
3980 if (fgVerboseLevel > 0)
3981 Info("Export", "Exporting %s %s as C++ code", GetName(), GetTitle());
3983 return 1;
3984 }
3985 if (sfile.Contains(".gdml")) {
3986 // Save geometry as a gdml file
3987 if (fgVerboseLevel > 0)
3988 Info("Export", "Exporting %s %s as gdml code", GetName(), GetTitle());
3989 // C++ version
3990 TString cmd;
3991 cmd = TString::Format("TGDMLWrite::StartGDMLWriting(gGeoManager,\"%s\",\"%s\")", filename, option);
3992 gROOT->ProcessLineFast(cmd);
3993 return 1;
3994 }
3995 if (sfile.Contains(".root") || sfile.Contains(".xml")) {
3996 // Save geometry as a root file
3997 TFile *f = TFile::Open(filename, "recreate");
3998 if (!f || f->IsZombie()) {
3999 Error("Export", "Cannot open file");
4000 return 0;
4001 }
4002 TString keyname = name;
4003 if (keyname.IsNull())
4004 keyname = GetName();
4005 TString opt = option;
4006 opt.ToLower();
4007 if (opt.Contains("v")) {
4009 if (fgVerboseLevel > 0)
4010 Info("Export", "Exporting %s %s as root file. Optimizations streamed.", GetName(), GetTitle());
4011 } else {
4013 if (fgVerboseLevel > 0)
4014 Info("Export", "Exporting %s %s as root file. Optimizations not streamed.", GetName(), GetTitle());
4015 }
4016
4017 const char *precision_dbl = TBufferText::GetDoubleFormat();
4018 const char *precision_flt = TBufferText::GetFloatFormat();
4019 TString new_format_dbl = TString::Format("%%.%dg", TGeoManager::GetExportPrecision());
4020 if (sfile.Contains(".xml")) {
4021 TBufferText::SetDoubleFormat(new_format_dbl.Data());
4022 TBufferText::SetFloatFormat(new_format_dbl.Data());
4023 }
4024 Int_t nbytes = Write(keyname);
4025 if (sfile.Contains(".xml")) {
4026 TBufferText::SetFloatFormat(precision_dbl);
4027 TBufferText::SetDoubleFormat(precision_flt);
4028 }
4029
4031 delete f;
4032 return nbytes;
4033 }
4034 return 0;
4035}
4036
4037////////////////////////////////////////////////////////////////////////////////
4038/// Lock current geometry so that no other geometry can be imported.
4039
4041{
4042 fgLock = kTRUE;
4043}
4044
4045////////////////////////////////////////////////////////////////////////////////
4046/// Unlock current geometry.
4047
4049{
4050 fgLock = kFALSE;
4051}
4052
4053////////////////////////////////////////////////////////////////////////////////
4054/// Check lock state.
4055
4057{
4058 return fgLock;
4059}
4060
4061////////////////////////////////////////////////////////////////////////////////
4062/// Set verbosity level (static function).
4063/// - 0 - suppress messages related to geom-painter visibility level
4064/// - 1 - default value
4065
4067{
4068 return fgVerboseLevel;
4069}
4070
4071////////////////////////////////////////////////////////////////////////////////
4072/// Return current verbosity level (static function).
4073
4075{
4076 fgVerboseLevel = vl;
4077}
4078
4079////////////////////////////////////////////////////////////////////////////////
4080/// static function
4081/// Import a geometry from a gdml or ROOT file
4082///
4083/// - Case 1: gdml
4084/// if filename ends with ".gdml" the foreign geometry described with gdml
4085/// is imported executing some python scripts in $ROOTSYS/gdml.
4086/// NOTE that to use this option, the PYTHONPATH must be defined like
4087/// export PYTHONPATH=$ROOTSYS/lib:$ROOTSYS/gdml
4088///
4089/// - Case 2: root file (.root) or root/xml file (.xml)
4090/// Import in memory from filename the geometry with key=name.
4091/// if name="" (default), the first TGeoManager object in the file is returned.
4092///
4093/// Note that this function deletes the current gGeoManager (if one)
4094/// before importing the new object.
4095
4096TGeoManager *TGeoManager::Import(const char *filename, const char *name, Option_t * /*option*/)
4097{
4098 if (fgLock) {
4099 ::Warning("TGeoManager::Import", "TGeoMananager in lock mode. NOT IMPORTING new geometry");
4100 return nullptr;
4101 }
4102 if (!filename)
4103 return nullptr;
4104 if (fgVerboseLevel > 0)
4105 ::Info("TGeoManager::Import", "Reading geometry from file: %s", filename);
4106
4107 if (gGeoManager)
4108 delete gGeoManager;
4109 gGeoManager = nullptr;
4110
4111 if (strstr(filename, ".gdml")) {
4112 // import from a gdml file
4113 new TGeoManager("GDMLImport", "Geometry imported from GDML");
4114 TString cmd = TString::Format("TGDMLParse::StartGDML(\"%s\")", filename);
4115 TGeoVolume *world = (TGeoVolume *)gROOT->ProcessLineFast(cmd);
4116
4117 if (world == nullptr) {
4118 delete gGeoManager;
4119 gGeoManager = nullptr;
4120 ::Error("TGeoManager::Import", "Cannot read file %s", filename);
4121 } else {
4122 gGeoManager->SetTopVolume(world);
4125 }
4126 } else {
4127 // import from a root file
4129 // in case a web file is specified, use the cacheread option to cache
4130 // this file in the cache directory
4131 TFile *f = nullptr;
4132 if (strstr(filename, "http"))
4133 f = TFile::Open(filename, "CACHEREAD");
4134 else
4136 if (!f || f->IsZombie()) {
4137 ::Error("TGeoManager::Import", "Cannot open file");
4138 return nullptr;
4139 }
4140 if (name && strlen(name) > 0) {
4141 gGeoManager = (TGeoManager *)f->Get(name);
4142 } else {
4143 TIter next(f->GetListOfKeys());
4144 TKey *key;
4145 while ((key = (TKey *)next())) {
4146 if (strcmp(key->GetClassName(), "TGeoManager") != 0)
4147 continue;
4148 gGeoManager = (TGeoManager *)key->ReadObj();
4149 break;
4150 }
4151 }
4152 delete f;
4153 }
4154 if (!gGeoManager)
4155 return nullptr;
4156 if (!gROOT->GetListOfGeometries()->FindObject(gGeoManager))
4157 gROOT->GetListOfGeometries()->Add(gGeoManager);
4158 if (!gROOT->GetListOfBrowsables()->FindObject(gGeoManager))
4159 gROOT->GetListOfBrowsables()->Add(gGeoManager);
4161 return gGeoManager;
4162}
4163
4164////////////////////////////////////////////////////////////////////////////////
4165/// Update element flags when geometry is loaded from a file.
4166
4168{
4169 if (!fElementTable)
4170 return;
4171 TIter next(fMaterials);
4172 TGeoMaterial *mat;
4173 TGeoMixture *mix;
4174 TGeoElement *elem, *elem_table;
4175 Int_t i, nelem;
4176 while ((mat = (TGeoMaterial *)next())) {
4177 if (mat->IsMixture()) {
4178 mix = (TGeoMixture *)mat;
4179 nelem = mix->GetNelements();
4180 for (i = 0; i < nelem; i++) {
4181 elem = mix->GetElement(i);
4182 if (!elem)
4183 continue;
4184 elem_table = fElementTable->GetElement(elem->Z());
4185 if (!elem_table)
4186 continue;
4187 if (elem != elem_table) {
4188 elem_table->SetDefined(elem->IsDefined());
4189 elem_table->SetUsed(elem->IsUsed());
4190 } else {
4191 elem_table->SetDefined();
4192 }
4193 }
4194 } else {
4195 elem = mat->GetElement();
4196 if (!elem)
4197 continue;
4198 elem_table = fElementTable->GetElement(elem->Z());
4199 if (!elem_table)
4200 continue;
4201 if (elem != elem_table) {
4202 elem_table->SetDefined(elem->IsDefined());
4203 elem_table->SetUsed(elem->IsUsed());
4204 } else {
4205 elem_table->SetUsed();
4206 }
4207 }
4208 }
4209}
4210
4211////////////////////////////////////////////////////////////////////////////////
4212/// Initialize PNE array for fast access via index and unique-id.
4213
4215{
4216 if (fHashPNE) {
4218 TIter next(fHashPNE);
4219 TObject *obj;
4220 while ((obj = next())) {
4221 fArrayPNE->Add(obj);
4222 }
4223 return kTRUE;
4224 }
4225 return kFALSE;
4226}
4227
4228////////////////////////////////////////////////////////////////////////////////
4229/// Get time cut for drawing tracks.
4230
4232{
4233 tmin = fTmin;
4234 tmax = fTmax;
4235 return fTimeCut;
4236}
4237
4238////////////////////////////////////////////////////////////////////////////////
4239/// Set time cut interval for drawing tracks. If called with no arguments, time
4240/// cut will be disabled.
4241
4243{
4244 fTmin = tmin;
4245 fTmax = tmax;
4246 if (tmin == 0 && tmax == 999)
4247 fTimeCut = kFALSE;
4248 else
4249 fTimeCut = kTRUE;
4250 if (fTracks && !IsAnimatingTracks())
4251 ModifiedPad();
4252}
4253
4254////////////////////////////////////////////////////////////////////////////////
4255/// Convert coordinates from master volume frame to top.
4256
4257void TGeoManager::MasterToTop(const Double_t *master, Double_t *top) const
4258{
4259 GetCurrentNavigator()->MasterToLocal(master, top);
4260}
4261
4262////////////////////////////////////////////////////////////////////////////////
4263/// Convert coordinates from top volume frame to master.
4264
4265void TGeoManager::TopToMaster(const Double_t *top, Double_t *master) const
4266{
4267 GetCurrentNavigator()->LocalToMaster(top, master);
4268}
4269
4270////////////////////////////////////////////////////////////////////////////////
4271/// Create a parallel world for prioritised navigation. This can be populated
4272/// with physical nodes and can be navigated independently using its API.
4273/// In case the flag SetUseParallelWorldNav is set, any navigation query in the
4274/// main geometry is checked against the parallel geometry, which gets priority
4275/// in case of overlaps with the main geometry volumes.
4276
4278{
4280 return fParallelWorld;
4281}
4282
4283////////////////////////////////////////////////////////////////////////////////
4284/// Activate/deactivate usage of parallel world navigation. Can only be done if
4285/// there is a parallel world. Activating navigation will automatically close
4286/// the parallel geometry.
4287
4289{
4290 if (!fParallelWorld) {
4291 Error("SetUseParallelWorldNav", "No parallel world geometry defined. Use CreateParallelWorld.");
4292 return;
4293 }
4294 if (!flag) {
4295 fUsePWNav = flag;
4296 return;
4297 }
4298 if (!fClosed) {
4299 Error("SetUseParallelWorldNav", "The geometry must be closed first");
4300 return;
4301 }
4302 // Closing the parallel world geometry is mandatory
4304 fUsePWNav = kTRUE;
4305}
4306
4308{
4309 Bool_t val = gGeometryLocked;
4310 gGeometryLocked = new_value;
4311 return val;
4312}
4313
4315{
4316 return fgDefaultUnits;
4317}
4318
4320{
4321 if (fgDefaultUnits == new_value) {
4322 gGeometryLocked = true;
4323 return;
4324 } else if (gGeometryLocked) {
4325 ::Fatal("TGeoManager", "The system of units may only be changed once, \n"
4326 "BEFORE any elements and materials are created! \n"
4327 "Alternatively unlock the default units at own risk.");
4328 } else if (new_value == kG4Units) {
4329 ::Info("TGeoManager", "Changing system of units to Geant4 units (mm, ns, MeV).");
4330 } else if (new_value == kRootUnits) {
4331 ::Info("TGeoManager", "Changing system of units to ROOT units (cm, s, GeV).");
4332 }
4333 fgDefaultUnits = new_value;
4334}
4335
4337{
4338 fgExportPrecision = prec;
4339}
4340
4342{
4343 return fgExportPrecision;
4344}
#define SafeDelete(p)
Definition RConfig.hxx:550
#define b(i)
Definition RSha256.hxx:100
#define f(i)
Definition RSha256.hxx:104
#define a(i)
Definition RSha256.hxx:99
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
unsigned char UChar_t
Definition RtypesCore.h:38
unsigned long ULong_t
Definition RtypesCore.h:55
unsigned int UInt_t
Definition RtypesCore.h:46
float Float_t
Definition RtypesCore.h:57
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 BIT(n)
Definition Rtypes.h:85
#define ClassImp(name)
Definition Rtypes.h:377
@ kGray
Definition Rtypes.h:65
@ kRed
Definition Rtypes.h:66
@ kOrange
Definition Rtypes.h:67
@ kGreen
Definition Rtypes.h:66
@ kMagenta
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
@ kYellow
Definition Rtypes.h:66
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
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 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 filename
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 const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
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 UChar_t len
Option_t Option_t TPoint TPoint const char mode
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 property
char name[80]
Definition TGX11.cxx:110
TGeoManager * gGeoManager
static Bool_t gGeometryLocked
R__EXTERN TGeoManager * gGeoManager
R__EXTERN TGeoIdentity * gGeoIdentity
Definition TGeoMatrix.h:537
int nentries
#define gROOT
Definition TROOT.h:406
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
#define gPad
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
Using a TBrowser one can browse all ROOT objects.
Definition TBrowser.h:37
void RecursiveRemove(TObject *obj) override
Recursively remove obj from browser.
Definition TBrowser.cxx:408
static const char * GetFloatFormat()
return current printf format for float members, default "%e"
static void SetFloatFormat(const char *fmt="%e")
set printf format for float/double members, default "%e" to change format only for doubles,...
static const char * GetDoubleFormat()
return current printf format for double members, default "%.14e"
static void SetDoubleFormat(const char *fmt="%.14e")
set printf format for double members, default "%.14e" use it after SetFloatFormat,...
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=nullptr)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
@ kRealNew
Definition TClass.h:108
@ kDummyNew
Definition TClass.h:108
static ENewType IsCallingNew()
Static method returning the defConstructor flag passed to TClass::New().
Definition TClass.cxx:5902
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition TDatime.h:37
const char * AsString() const
Return the date & time as a string (ctime() format).
Definition TDatime.cxx:102
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
A ROOT file is composed of a header, followed by consecutive data records (TKey instances) with a wel...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4070
This class is used in the process of reading and writing the GDML "matrix" tag.
Definition TGDMLMatrix.h:33
Bool_t IsVisTouched() const
Definition TGeoAtt.h:91
void SetVisStreamed(Bool_t vis=kTRUE)
Mark attributes as "streamed to file".
Definition TGeoAtt.cxx:128
void SetVisTouched(Bool_t vis=kTRUE)
Mark visualization attributes as "modified".
Definition TGeoAtt.cxx:138
void SetVisBranch()
Set branch type visibility.
Definition TGeoAtt.cxx:66
Box class.
Definition TGeoBBox.h:17
void Node(const char *name, Int_t nr, const char *mother, Double_t x, Double_t y, Double_t z, Int_t irot, Bool_t isOnly, Float_t *upar, Int_t npar=0)
Create a node called <name_nr> pointing to the volume called <name> as daughter of the volume called ...
TGeoVolume * MakePgon(const char *name, TGeoMedium *medium, Double_t phi, Double_t dphi, Int_t nedges, Int_t nz)
Make in one step a volume pointing to a polygone shape with given medium.
TGeoVolume * MakeGtra(const char *name, TGeoMedium *medium, Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2, Double_t tl2, Double_t alpha2)
Make in one step a volume pointing to a twisted trapezoid shape with given medium.
TGeoVolume * MakeXtru(const char *name, TGeoMedium *medium, Int_t nz)
Make a TGeoXtru-shaped volume with nz planes.
TGeoVolume * MakePcon(const char *name, TGeoMedium *medium, Double_t phi, Double_t dphi, Int_t nz)
Make in one step a volume pointing to a polycone shape with given medium.
TGeoVolume * MakeTrap(const char *name, TGeoMedium *medium, Double_t dz, Double_t theta, Double_t phi, Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2, Double_t tl2, Double_t alpha2)
Make in one step a volume pointing to a trapezoid shape with given medium.
Int_t AddShape(TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
TGeoVolume * MakeTorus(const char *name, TGeoMedium *medium, Double_t r, Double_t rmin, Double_t rmax, Double_t phi1=0, Double_t dphi=360)
Make in one step a volume pointing to a torus shape with given medium.
static TGeoBuilder * Instance(TGeoManager *geom)
Return pointer to singleton.
TGeoVolume * MakeEltu(const char *name, TGeoMedium *medium, Double_t a, Double_t b, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
TGeoVolume * MakeHype(const char *name, TGeoMedium *medium, Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
TGeoVolume * MakeArb8(const char *name, TGeoMedium *medium, Double_t dz, Double_t *vertices=nullptr)
Make an TGeoArb8 volume.
TGeoMaterial * Material(const char *name, Double_t a, Double_t z, Double_t dens, Int_t uid, Double_t radlen=0, Double_t intlen=0)
Create material with given A, Z and density, having an unique id.
TGeoVolume * MakeTube(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
TGeoVolume * MakePara(const char *name, TGeoMedium *medium, Double_t dx, Double_t dy, Double_t dz, Double_t alpha, Double_t theta, Double_t phi)
Make in one step a volume pointing to a parallelepiped shape with given medium.
Int_t AddTransformation(TGeoMatrix *matrix)
Add a matrix to the list. Returns index of the matrix in list.
TGeoVolume * MakeSphere(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t themin=0, Double_t themax=180, Double_t phimin=0, Double_t phimax=360)
Make in one step a volume pointing to a sphere shape with given medium.
TGeoMedium * Medium(const char *name, Int_t numed, Int_t nmat, Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd, Double_t stemax, Double_t deemax, Double_t epsil, Double_t stmin)
Create tracking medium.
TGeoVolume * MakeTubs(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
Make in one step a volume pointing to a tube segment shape with given medium.
TGeoVolume * MakeBox(const char *name, TGeoMedium *medium, Double_t dx, Double_t dy, Double_t dz)
Make in one step a volume pointing to a box shape with given medium.
void Matrix(Int_t index, Double_t theta1, Double_t phi1, Double_t theta2, Double_t phi2, Double_t theta3, Double_t phi3)
Create rotation matrix named 'mat<index>'.
TGeoVolume * Division(const char *name, const char *mother, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
Create a new volume by dividing an existing one (GEANT3 like)
TGeoVolume * MakeParaboloid(const char *name, TGeoMedium *medium, Double_t rlo, Double_t rhi, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
TGeoVolume * MakeTrd1(const char *name, TGeoMedium *medium, Double_t dx1, Double_t dx2, Double_t dy, Double_t dz)
Make in one step a volume pointing to a TGeoTrd1 shape with given medium.
TGeoVolumeAssembly * MakeVolumeAssembly(const char *name)
Make an assembly of volumes.
TGeoVolume * MakeCons(const char *name, TGeoMedium *medium, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
Make in one step a volume pointing to a cone segment shape with given medium.
Int_t AddMaterial(TGeoMaterial *material)
Add a material to the list. Returns index of the material in list.
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
TGeoVolume * MakeTrd2(const char *name, TGeoMedium *medium, Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2, Double_t dz)
Make in one step a volume pointing to a TGeoTrd2 shape with given medium.
TGeoVolume * Volume(const char *name, const char *shape, Int_t nmed, Float_t *upar, Int_t npar=0)
Create a volume in GEANT3 style.
TGeoMaterial * Mixture(const char *name, Float_t *a, Float_t *z, Double_t dens, Int_t nelem, Float_t *wmat, Int_t uid)
Create mixture OR COMPOUND IMAT as composed by THE BASIC nelem materials defined by arrays A,...
TGeoVolume * MakeCone(const char *name, TGeoMedium *medium, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Make in one step a volume pointing to a cone shape with given medium.
void RegisterMatrix(TGeoMatrix *matrix)
Register a matrix to the list of matrices.
TGeoVolume * MakeCtub(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
Make in one step a volume pointing to a tube segment shape with given medium.
Class describing rotation + translation.
Definition TGeoMatrix.h:317
Composite shapes are Boolean combinations of two or more shape components.
Table of elements.
TGeoElement * GetElement(Int_t z)
Base class for chemical elements.
Definition TGeoElement.h:36
Bool_t IsDefined() const
Definition TGeoElement.h:81
void SetDefined(Bool_t flag=kTRUE)
Definition TGeoElement.h:85
Int_t Z() const
Definition TGeoElement.h:68
void SetUsed(Bool_t flag=kTRUE)
Definition TGeoElement.h:86
Bool_t IsUsed() const
Definition TGeoElement.h:83
Matrix class used for computing global transformations Should NOT be used for node definition.
Definition TGeoMatrix.h:458
An identity transformation.
Definition TGeoMatrix.h:406
A geometry iterator.
Definition TGeoNode.h:248
Int_t GetLevel() const
Definition TGeoNode.h:294
The manager class for any TGeo geometry.
Definition TGeoManager.h:44
static void UnlockGeometry()
Unlock current geometry.
Double_t fPhimax
lowest range for phi cut
Definition TGeoManager.h:63
TGeoVolume * MakeCone(const char *name, TGeoMedium *medium, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2)
Make in one step a volume pointing to a cone shape with given medium.
void AnimateTracks(Double_t tmin=0, Double_t tmax=5E-8, Int_t nframes=200, Option_t *option="/*")
Draw animation of tracks.
void AddSkinSurface(TGeoSkinSurface *surf)
Add skin surface;.
TGeoVolume * MakeXtru(const char *name, TGeoMedium *medium, Int_t nz)
Make a TGeoXtru-shaped volume with nz planes.
Double_t * FindNormalFast()
Computes fast normal to next crossed boundary, assuming that the current point is close enough to the...
TGeoVolume * MakePcon(const char *name, TGeoMedium *medium, Double_t phi, Double_t dphi, Int_t nz)
Make in one step a volume pointing to a polycone shape with given medium.
Int_t fRaytraceMode
Flag for multi-threading.
Double_t fVisDensity
particles to be drawn
Definition TGeoManager.h:69
TGeoNavigator * AddNavigator()
Add a navigator in the list of navigators.
TVirtualGeoTrack * GetTrackOfId(Int_t id) const
Get track with a given ID.
TGeoMaterial * FindDuplicateMaterial(const TGeoMaterial *mat) const
Find if a given material duplicates an existing one.
TGeoVolume * Division(const char *name, const char *mother, Int_t iaxis, Int_t ndiv, Double_t start, Double_t step, Int_t numed=0, Option_t *option="")
Create a new volume by dividing an existing one (GEANT3 like)
TGeoVolume * Volume(const char *name, const char *shape, Int_t nmed, Float_t *upar, Int_t npar=0)
Create a volume in GEANT3 style.
Int_t ReplaceVolume(TGeoVolume *vorig, TGeoVolume *vnew)
Replaces all occurrences of VORIG with VNEW in the geometry tree.
void DoRestoreState()
Restore a backed-up state without affecting the cache stack.
Int_t GetCurrentNodeId() const
Get the unique ID of the current node.
TGeoPNEntry * GetAlignableEntry(const char *name) const
Retrieves an existing alignable object.
TGeoVolume * fMasterVolume
top physical node
TVirtualGeoTrack * FindTrackWithId(Int_t id) const
Search the track hierarchy to find the track with the given id.
TObjArray * fArrayPNE
void TestOverlaps(const char *path="")
Geometry overlap checker based on sampling.
static EDefaultUnits GetDefaultUnits()
void RemoveMaterial(Int_t index)
Remove material at given index.
void Matrix(Int_t index, Double_t theta1, Double_t phi1, Double_t theta2, Double_t phi2, Double_t theta3, Double_t phi3)
Create rotation matrix named 'mat<index>'.
TGeoElementTable * GetElementTable()
Returns material table. Creates it if not existing.
Int_t fNtracks
Definition TGeoManager.h:74
THashList * fHashPNE
hash list of group volumes providing fast search
static Int_t fgVerboseLevel
Lock preventing a second geometry to be loaded.
Definition TGeoManager.h:51
void Init()
Initialize manager class.
Bool_t InitArrayPNE() const
Initialize PNE array for fast access via index and unique-id.
TObjArray * fPhysicalNodes
Definition TGeoManager.h:96
virtual ULong_t SizeOf(const TGeoNode *node, Option_t *option)
computes the total size in bytes of the branch starting with node.
TObjArray * fUniqueVolumes
static UInt_t fgExportPrecision
Maximum number of Xtru vertices.
Definition TGeoManager.h:55
TObjArray * fRegions
void Node(const char *name, Int_t nr, const char *mother, Double_t x, Double_t y, Double_t z, Int_t irot, Bool_t isOnly, Float_t *upar, Int_t npar=0)
Create a node called <name_nr> pointing to the volume called <name> as daughter of the volume called ...
TObjArray * fGShapes
Definition TGeoManager.h:97
TGeoVolume * fPaintVolume
TGeoSkinSurface * GetSkinSurface(const char *name) const
Get skin surface with a given name;.
void UpdateElements()
Update element flags when geometry is loaded from a file.
TGeoManager()
Default constructor.
static TClass * Class()
ConstPropMap_t fProperties
TGeoVolume * MakeTube(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
void CdUp()
Go one level up in geometry.
void DoBackupState()
Backup the current state without affecting the cache stack.
TList * fMaterials
void CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.)
Check pushes and pulls needed to cross the next boundary with respect to the position given by FindNe...
TObjArray * fVolumes
Definition TGeoManager.h:95
Int_t * fValuePNEId
TGeoPNEntry * GetAlignableEntryByUID(Int_t uid) const
Retrieves an existing alignable object having a preset UID.
void AddGDMLMatrix(TGDMLMatrix *mat)
Add GDML matrix;.
Bool_t fTimeCut
Definition TGeoManager.h:85
static void SetExportPrecision(UInt_t prec)
void AddBorderSurface(TGeoBorderSurface *surf)
Add border surface;.
void SetClippingShape(TGeoShape *clip)
Set a user-defined shape as clipping for ray tracing.
TGeoVolume * fCurrentVolume
current navigator
void ClearOverlaps()
Clear the list of overlaps.
TGeoVolume * MakeCons(const char *name, TGeoMedium *medium, Double_t dz, Double_t rmin1, Double_t rmax1, Double_t rmin2, Double_t rmax2, Double_t phi1, Double_t phi2)
Make in one step a volume pointing to a cone segment shape with given medium.
THashList * fHashGVolumes
hash list of volumes providing fast search
Int_t fVisOption
Definition TGeoManager.h:71
static std::mutex fgMutex
Definition TGeoManager.h:49
Bool_t IsInPhiRange() const
True if current node is in phi range.
virtual Bool_t cd(const char *path="")
Browse the tree of nodes starting from fTopNode according to pathname.
TGeoNode * SearchNode(Bool_t downwards=kFALSE, const TGeoNode *skipnode=nullptr)
Returns the deepest node containing fPoint, which must be set a priori.
TGeoMaterial * Material(const char *name, Double_t a, Double_t z, Double_t dens, Int_t uid, Double_t radlen=0, Double_t intlen=0)
Create material with given A, Z and density, having an unique id.
void LocalToMaster(const Double_t *local, Double_t *master) const
Double_t fPhimin
Definition TGeoManager.h:62
static Bool_t fgLockNavigators
Number of registered threads.
void SaveAttributes(const char *filename="tgeoatt.C")
Save current attributes in a macro.
void RestoreMasterVolume()
Restore the master volume of the geometry.
Bool_t fDrawExtra
Definition TGeoManager.h:86
TGeoVolume * MakeArb8(const char *name, TGeoMedium *medium, Double_t dz, Double_t *vertices=nullptr)
Make an TGeoArb8 volume.
virtual Int_t Export(const char *filename, const char *name="", Option_t *option="vg")
Export this geometry to a file.
TGeoNode * FindNextDaughterBoundary(Double_t *point, Double_t *dir, Int_t &idaughter, Bool_t compmatrix=kFALSE)
Computes as fStep the distance to next daughter of the current volume.
Int_t GetUID(const char *volname) const
Retrieve unique id for a volume name. Return -1 if name not found.
TGeoShape * fClippingShape
TGeoNavigator * GetCurrentNavigator() const
Returns current navigator for the calling thread.
THashList * fHashVolumes
TObjArray * fMatrices
current painter
Definition TGeoManager.h:93
static Int_t GetNumThreads()
Returns number of threads that were set to use geometry.
TGeoVolumeMulti * MakeVolumeMulti(const char *name, TGeoMedium *medium)
Make a TGeoVolumeMulti handling a list of volumes.
void ClearNavigators()
Clear all navigators.
Int_t AddTransformation(const TGeoMatrix *matrix)
Add a matrix to the list. Returns index of the matrix in list.
TObjArray * fOpticalSurfaces
TVirtualGeoTrack * GetParentTrackOfId(Int_t id) const
Get parent track with a given ID.
void CdNode(Int_t nodeid)
Change current path to point to the node having this id.
UChar_t * fBits
static Int_t GetMaxLevels()
Return maximum number of levels used in the geometry.
Double_t fTmin
highest range for phi cut
Definition TGeoManager.h:64
static Bool_t IsLocked()
Check lock state.
TGeoVolume * fTopVolume
current volume
TGeoVolume * fUserPaintVolume
volume currently painted
TVirtualGeoPainter * GetGeomPainter()
Make a default painter if none present. Returns pointer to it.
void GetBranchOnlys(Int_t *isonly) const
Fill node copy numbers of current branch into an array.
TGeoNode * GetCurrentNode() const
Int_t AddTrack(Int_t id, Int_t pdgcode, TObject *particle=nullptr)
Add a track to the list of tracks.
void SetVisOption(Int_t option=0)
set drawing mode :
void SetPdgName(Int_t pdg, const char *name)
Set a name for a particle having a given pdg.
TObjArray * fBorderSurfaces
Int_t GetNAlignable(Bool_t with_uid=kFALSE) const
Retrieves number of PN entries with or without UID.
void RefreshPhysicalNodes(Bool_t lock=kTRUE)
Refresh physical nodes to reflect the actual geometry paths after alignment was applied.
static Bool_t fgLock
mutex for navigator booking in MT mode
Definition TGeoManager.h:50
TGeoVolume * MakePara(const char *name, TGeoMedium *medium, Double_t dx, Double_t dy, Double_t dz, Double_t alpha, Double_t theta, Double_t phi)
Make in one step a volume pointing to a parallelepiped shape with given medium.
void TopToMaster(const Double_t *top, Double_t *master) const
Convert coordinates from top volume frame to master.
TObjArray * fShapes
Definition TGeoManager.h:94
void AddOpticalSurface(TGeoOpticalSurface *optsurf)
Add optical surface;.
static void SetDefaultUnits(EDefaultUnits new_value)
Bool_t fLoopVolumes
flag that geometry is closed
Definition TGeoManager.h:80
Int_t AddMaterial(const TGeoMaterial *material)
Add a material to the list. Returns index of the material in list.
void ClearAttributes()
Reset all attributes to default ones.
static Int_t fgMaxDaughters
Maximum level in geometry.
Definition TGeoManager.h:53
Bool_t fUsePWNav
Raytrace mode: 0=normal, 1=pass through, 2=transparent.
void SetRTmode(Int_t mode)
Change raytracing mode.
Bool_t CheckPath(const char *path) const
Check if a geometry path is valid without changing the state of the current navigator.
void InspectState() const
Inspects path and all flags for the current state.
void ConvertReflections()
Convert all reflections in geometry to normal rotations + reflected shapes.
void SetVisLevel(Int_t level=3)
set default level down to which visualization is performed
TGeoNode * FindNextBoundary(Double_t stepmax=TGeoShape::Big(), const char *path="", Bool_t frombdr=kFALSE)
Find distance to next boundary and store it in fStep.
static TGeoManager * Import(const char *filename, const char *name="", Option_t *option="")
static function Import a geometry from a gdml or ROOT file
TGeoPhysicalNode * MakePhysicalNode(const char *path=nullptr)
Makes a physical node corresponding to a path.
void CountLevels()
Count maximum number of nodes per volume, maximum depth and maximum number of xtru vertices.
Int_t fMaxThreads
Bool_t fIsGeomReading
Definition TGeoManager.h:82
TGeoVolume * MakeTorus(const char *name, TGeoMedium *medium, Double_t r, Double_t rmin, Double_t rmax, Double_t phi1=0, Double_t dphi=360)
Make in one step a volume pointing to a torus shape with given medium.
TGeoHMatrix * GetHMatrix()
Return stored current matrix (global matrix of the next touched node).
TGeoParallelWorld * fParallelWorld
void RegisterMatrix(const TGeoMatrix *matrix)
Register a matrix to the list of matrices.
TVirtualGeoTrack * GetTrack(Int_t index)
static Int_t GetMaxDaughters()
Return maximum number of daughters of a volume used in the geometry.
static void ClearThreadsMap()
Clear the current map of threads.
Int_t AddVolume(TGeoVolume *volume)
Add a volume to the list. Returns index of the volume in list.
TVirtualGeoPainter * fPainter
flag that nodes are the selected objects in pad rather than volumes
Definition TGeoManager.h:91
void SetVolumeAttribute(const char *name, const char *att, Int_t val)
Set volume attributes in G3 style.
const char * GetPdgName(Int_t pdg) const
Get name for given pdg code;.
void CheckGeometryFull(Int_t ntracks=1000000, Double_t vx=0., Double_t vy=0., Double_t vz=0., Option_t *option="ob")
Geometry checking.
Bool_t fIsNodeSelectable
switch ON/OFF volume activity (default OFF - all volumes active))
Definition TGeoManager.h:90
TGeoNode * Step(Bool_t is_geom=kTRUE, Bool_t cross=kTRUE)
Make a rectilinear step of length fStep from current point (fPoint) on current direction (fDirection)...
Bool_t GotoSafeLevel()
Go upwards the tree until a non-overlapping node.
Bool_t fActivity
flag for GL reflections
Definition TGeoManager.h:89
void GetBranchNames(Int_t *names) const
Fill volume names of current branch into an array.
static ThreadsMap_t * fgThreadId
Map between thread id's and navigator arrays.
void CloseGeometry(Option_t *option="d")
Closing geometry implies checking the geometry validity, fixing shapes with negative parameters (run-...
TVirtualGeoTrack * MakeTrack(Int_t id, Int_t pdgcode, TObject *particle)
Makes a primary track but do not attach it to the list of tracks.
Int_t GetTrackIndex(Int_t id) const
Get index for track id, -1 if not found.
Int_t fNNodes
upper time limit for tracks drawing
Definition TGeoManager.h:66
Int_t fNLevel
table of elements
void OptimizeVoxels(const char *filename="tgeovox.C")
Optimize voxelization type for all volumes. Save best choice in a macro.
TGeoVolume * GetVolume(const char *name) const
Search for a named volume. All trailing blanks stripped.
void SetAnimateTracks(Bool_t flag=kTRUE)
Bool_t fIsGeomCleaning
flag set when reading geometry
Definition TGeoManager.h:83
Bool_t IsSameLocation() const
void DefaultColors()
Set default volume colors according to A of material.
TGeoNode * FindNextBoundaryAndStep(Double_t stepmax=TGeoShape::Big(), Bool_t compsafe=kFALSE)
Compute distance to next boundary within STEPMAX.
TGeoVolume * MakeTrd2(const char *name, TGeoMedium *medium, Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2, Double_t dz)
Make in one step a volume pointing to a TGeoTrd2 shape with given medium.
Double_t * FindNormal(Bool_t forward=kTRUE)
Computes normal vector to the next surface that will be or was already crossed when propagating on a ...
virtual Int_t GetByteCount(Option_t *option=nullptr)
Get total size of geometry in bytes.
TGeoVolume * MakeGtra(const char *name, TGeoMedium *medium, Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2, Double_t tl2, Double_t alpha2)
Make in one step a volume pointing to a twisted trapezoid shape with given medium.
TGeoElementTable * fElementTable
clipping shape for raytracing
static void SetNavigatorsLock(Bool_t flag)
Set the lock for navigators.
static Int_t fgMaxXtruVert
Maximum number of daughters.
Definition TGeoManager.h:54
TGeoNode * FindNode(Bool_t safe_start=kTRUE)
Returns deepest node containing current point.
Int_t GetVisOption() const
Returns current depth to which geometry is drawn.
static void LockGeometry()
Lock current geometry so that no other geometry can be imported.
TGeoVolume * MakeBox(const char *name, TGeoMedium *medium, Double_t dx, Double_t dy, Double_t dz)
Make in one step a volume pointing to a box shape with given medium.
void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)
Test for shape navigation methods.
static Int_t fgMaxLevel
Verbosity level for Info messages (no IO).
Definition TGeoManager.h:52
Int_t fNpdg
current track
Definition TGeoManager.h:77
void PrintOverlaps() const
Prints the current list of overlaps.
TGeoVolume * MakeTrd1(const char *name, TGeoMedium *medium, Double_t dx1, Double_t dx2, Double_t dy, Double_t dz)
Make in one step a volume pointing to a TGeoTrd1 shape with given medium.
TGeoVolume * MakeSphere(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t themin=0, Double_t themax=180, Double_t phimin=0, Double_t phimax=360)
Make in one step a volume pointing to a sphere shape with given medium.
void ResetUserData()
Sets all pointers TGeoVolume::fField to NULL.
TGeoVolume * FindVolumeFast(const char *name, Bool_t multi=kFALSE)
Fast search for a named volume. All trailing blanks stripped.
TList * fMedia
Bool_t GetTminTmax(Double_t &tmin, Double_t &tmax) const
Get time cut for drawing tracks.
TGeoNode * InitTrack(const Double_t *point, const Double_t *dir)
Initialize current point and current direction vector (normalized) in MARS.
ThreadsMap_t::const_iterator ThreadsMapIt_t
Bool_t fMatrixTransform
flag that the list of physical nodes has to be drawn
Definition TGeoManager.h:87
void SetVisibility(TObject *obj, Bool_t vis)
Set visibility for a volume.
void SetTopVolume(TGeoVolume *vol)
Set the top volume and corresponding node as starting point of the geometry.
Bool_t fMatrixReflection
flag for using GL matrix
Definition TGeoManager.h:88
TGeoPNEntry * SetAlignableEntry(const char *unique_name, const char *path, Int_t uid=-1)
Creates an alignable object with unique name corresponding to a path and adds it to the list of align...
void ClearShape(const TGeoShape *shape)
Remove a shape from the list of shapes.
void ModifiedPad() const
Send "Modified" signal to painter.
void BombTranslation(const Double_t *tr, Double_t *bombtr)
Get the new 'bombed' translation vector according current exploded view mode.
TGeoNavigator * fCurrentNavigator
Lock existing navigators.
static Bool_t LockDefaultUnits(Bool_t new_value)
Int_t fMaxVisNodes
Definition TGeoManager.h:75
TGeoMedium * GetMedium(const char *medium) const
Search for a named tracking medium. All trailing blanks stripped.
Bool_t InsertPNEId(Int_t uid, Int_t ientry)
Insert a PN entry in the sorted array of indexes.
Int_t fVisLevel
Definition TGeoManager.h:72
void ViewLeaves(Bool_t flag=kTRUE)
Set visualization option (leaves only OR all volumes)
TGeoVolume * MakeCtub(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2, Double_t lx, Double_t ly, Double_t lz, Double_t tx, Double_t ty, Double_t tz)
Make in one step a volume pointing to a tube segment shape with given medium.
void SetTminTmax(Double_t tmin=0, Double_t tmax=999)
Set time cut interval for drawing tracks.
NavigatorsMap_t fNavigators
void GetBranchNumbers(Int_t *copyNumbers, Int_t *volumeNumbers) const
Fill node copy numbers of current branch into an array.
Bool_t fPhiCut
flag to notify that the manager is being destructed
Definition TGeoManager.h:84
TGeoNode * CrossBoundaryAndLocate(Bool_t downwards, TGeoNode *skipnode)
Cross next boundary and locate within current node The current point must be on the boundary of fCurr...
void DrawTracks(Option_t *option="")
Draw tracks over the geometry, according to option.
void Streamer(TBuffer &) override
Stream an object of class TGeoManager.
void BuildDefaultMaterials()
Now just a shortcut for GetElementTable.
void SetMaxThreads(Int_t nthreads)
Set maximum number of threads for navigation.
TGeoMedium * Medium(const char *name, Int_t numed, Int_t nmat, Int_t isvol, Int_t ifield, Double_t fieldm, Double_t tmaxfd, Double_t stemax, Double_t deemax, Double_t epsil, Double_t stmin)
Create tracking medium.
void SetExplodedView(Int_t iopt=0)
Set type of exploding view (see TGeoPainter::SetExplodedView())
Double_t Weight(Double_t precision=0.01, Option_t *option="va")
Estimate weight of volume VOL with a precision SIGMA(W)/W better than PRECISION.
void ClearPhysicalNodes(Bool_t mustdelete=kFALSE)
Clear the current list of physical nodes, so that we can start over with a new list.
static Int_t Parse(const char *expr, TString &expr1, TString &expr2, TString &expr3)
Parse a string boolean expression and do a syntax check.
void GetBombFactors(Double_t &bombx, Double_t &bomby, Double_t &bombz, Double_t &bombr) const
Retrieve cartesian and radial bomb factors.
Double_t GetProperty(const char *name, Bool_t *error=nullptr) const
Get a user-defined property.
TObjArray * fTracks
list of runtime volumes
Definition TGeoManager.h:99
Bool_t IsAnimatingTracks() const
const char * GetPath() const
Get path to the current node in the form /node0/node1/...
static Int_t fgNumThreads
Thread id's map.
TObjArray * fGDMLMatrices
TGeoPhysicalNode * MakeAlignablePN(const char *name)
Make a physical node from the path pointed by an alignable object with a given name.
void SetCheckedNode(TGeoNode *node)
Assign a given node to be checked for overlaps. Any other overlaps will be ignored.
Int_t AddOverlap(const TNamed *ovlp)
Add an illegal overlap/extrusion to the list.
void CreateThreadData() const
Create thread private data for all geometry objects.
Int_t fNsegments
Definition TGeoManager.h:73
TObjArray * fOverlaps
TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil=1E-5, const char *g3path="")
shoot npoints randomly in a box of 1E-5 around current point.
Bool_t IsMultiThread() const
TGDMLMatrix * GetGDMLMatrix(const char *name) const
Get GDML matrix with a given name;.
Double_t fTmax
lower time limit for tracks drawing
Definition TGeoManager.h:65
Int_t TransformVolumeToAssembly(const char *vname)
Transform all volumes named VNAME to assemblies. The volumes must be virtual.
Bool_t fMultiThread
Max number of threads.
TGeoVolume * MakePgon(const char *name, TGeoMedium *medium, Double_t phi, Double_t dphi, Int_t nedges, Int_t nz)
Make in one step a volume pointing to a polygone shape with given medium.
TGeoVolume * MakeTrap(const char *name, TGeoMedium *medium, Double_t dz, Double_t theta, Double_t phi, Double_t h1, Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2, Double_t tl2, Double_t alpha2)
Make in one step a volume pointing to a trapezoid shape with given medium.
void DrawCurrentPoint(Int_t color=2)
Draw current point in the same view.
static void SetVerboseLevel(Int_t vl)
Return current verbosity level (static function).
TGeoOpticalSurface * GetOpticalSurface(const char *name) const
Get optical surface with a given name;.
void SetNsegments(Int_t nseg)
Set number of segments for approximating circles in drawing.
static UInt_t GetExportPrecision()
Bool_t IsSamePoint(Double_t x, Double_t y, Double_t z) const
Check if a new point with given coordinates is the same as the last located one.
void SetNmeshPoints(Int_t npoints=1000)
Set the number of points to be generated on the shape outline when checking for overlaps.
void CheckBoundaryReference(Int_t icheck=-1)
Check the boundary errors reference file created by CheckBoundaryErrors method.
static Int_t GetVerboseLevel()
Set verbosity level (static function).
Int_t GetVisLevel() const
Returns current depth to which geometry is drawn.
static EDefaultUnits fgDefaultUnits
Precision to be used in ASCII exports.
Definition TGeoManager.h:56
virtual void Edit(Option_t *option="")
Append a pad for this geometry.
Bool_t AddProperty(const char *property, Double_t value)
Add a user-defined property. Returns true if added, false if existing.
~TGeoManager() override
Destructor.
TObjArray * fNodes
Int_t CountNodes(const TGeoVolume *vol=nullptr, Int_t nlevels=10000, Int_t option=0)
Count the total number of nodes starting from a volume, nlevels down.
TGeoMaterial * GetMaterial(const char *matname) const
Search for a named material. All trailing blanks stripped.
void CheckGeometry(Option_t *option="")
Perform last checks on the geometry.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute mouse actions on this manager.
TGeoVolumeAssembly * MakeVolumeAssembly(const char *name)
Make an assembly of volumes.
Int_t GetBombMode() const
Int_t AddRegion(TGeoRegion *region)
Add a new region of volumes.
void SelectTrackingMedia()
Define different tracking media.
void CdNext()
Do a cd to the node found next by FindNextBoundary.
void CdTop()
Make top level node the current node.
Double_t Safety(Bool_t inside=kFALSE)
Compute safe distance from the current point.
Int_t * fKeyPNEId
void DefaultAngles()
Set default angles for a given view.
TGeoMaterial * Mixture(const char *name, Float_t *a, Float_t *z, Double_t dens, Int_t nelem, Float_t *wmat, Int_t uid)
Create mixture OR COMPOUND IMAT as composed by THE BASIC nelem materials defined by arrays A,...
std::map< std::thread::id, Int_t > ThreadsMap_t
void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="", Double_t safety=0.)
Classify a given point. See TGeoChecker::CheckPoint().
void SetUseParallelWorldNav(Bool_t flag)
Activate/deactivate usage of parallel world navigation.
void Browse(TBrowser *b) override
Describe how to browse this object.
void Test(Int_t npoints=1000000, Option_t *option="")
Check time of finding "Where am I" for n points.
Int_t GetSafeLevel() const
Go upwards the tree until a non-overlapping node.
TObjArray * fGVolumes
list of runtime shapes
Definition TGeoManager.h:98
void SetBombFactors(Double_t bombx=1.3, Double_t bomby=1.3, Double_t bombz=1.3, Double_t bombr=1.3)
Set factors that will "bomb" all translations in cartesian and cylindrical coordinates.
TGeoNode * fTopNode
top level volume in geometry
void UnbombTranslation(const Double_t *tr, Double_t *bombtr)
Get the new 'unbombed' translation vector according current exploded view mode.
void ResetState()
Reset current state flags.
void RandomPoints(const TGeoVolume *vol, Int_t npoints=10000, Option_t *option="")
Draw random points in the bounding box of a volume.
TGeoParallelWorld * CreateParallelWorld(const char *name)
Create a parallel world for prioritised navigation.
Int_t GetMaterialIndex(const char *matname) const
Return index of named material.
void CdDown(Int_t index)
Make a daughter of current node current.
TGeoNavigatorArray * GetListOfNavigators() const
Get list of navigators for the calling thread.
static Int_t GetMaxXtruVert()
Return maximum number of vertices for an xtru shape used.
void RandomRays(Int_t nrays=1000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=nullptr, Bool_t check_norm=kFALSE)
Randomly shoot nrays and plot intersections with surfaces for current top node.
void SetAllIndex()
Assigns uid's for all materials,media and matrices.
TObjArray * fSkinSurfaces
void SetVisDensity(Double_t dens=0.01)
Set density threshold.
Int_t fExplodedView
Definition TGeoManager.h:70
Bool_t fClosed
Definition TGeoManager.h:79
Int_t GetNsegments() const
Get number of segments approximating circles.
void SetPhiRange(Double_t phimin=0., Double_t phimax=360.)
Set cut phi range.
TGeoHMatrix * fGLMatrix
TVirtualGeoTrack * fCurrentTrack
Definition TGeoManager.h:76
TObjArray * fPdgNames
void DrawPath(const char *path, Option_t *option="")
Draw current path.
TObjArray * GetListOfPhysicalNodes()
static Int_t ThreadId()
Translates the current thread id to an ordinal number.
Bool_t SetCurrentNavigator(Int_t index)
Switch to another existing navigator for the calling thread.
void SetTopVisible(Bool_t vis=kTRUE)
make top volume visible on screen
TGeoVolume * MakeHype(const char *name, TGeoMedium *medium, Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
TGeoVolume * MakeParaboloid(const char *name, TGeoMedium *medium, Double_t rlo, Double_t rhi, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
Int_t AddShape(const TGeoShape *shape)
Add a shape to the list. Returns index of the shape in list.
void SetMaxVisNodes(Int_t maxnodes=10000)
set the maximum number of visible nodes.
void CleanGarbage()
Clean temporary volumes and shapes from garbage collection.
void Voxelize(Option_t *option=nullptr)
Voxelize all non-divided volumes.
Int_t GetVirtualLevel()
Find level of virtuality of current overlapping node (number of levels up having the same tracking me...
TGeoBorderSurface * GetBorderSurface(const char *name) const
Get border surface with a given name;.
void ClearThreadData() const
Int_t fSizePNEId
array of physical node entries
void CheckOverlaps(Double_t ovlp=0.1, Option_t *option="")
Check all geometry for illegal overlaps within a limit OVLP.
TGeoVolume * MakeTubs(const char *name, TGeoMedium *medium, Double_t rmin, Double_t rmax, Double_t dz, Double_t phi1, Double_t phi2)
Make in one step a volume pointing to a tube segment shape with given medium.
Int_t fPdgId[1024]
Definition TGeoManager.h:78
void SortOverlaps()
Sort overlaps by decreasing overlap distance. Extrusions comes first.
TGeoVolume * MakeEltu(const char *name, TGeoMedium *medium, Double_t a, Double_t b, Double_t dz)
Make in one step a volume pointing to a tube shape with given medium.
void RemoveNavigator(const TGeoNavigator *nav)
Clear a single navigator.
void MasterToTop(const Double_t *master, Double_t *top) const
Convert coordinates from master volume frame to top.
Bool_t fStreamVoxels
flag volume lists loop
Definition TGeoManager.h:81
Base class describing materials.
virtual Bool_t IsMixture() const
virtual Int_t GetByteCount() const
TGeoElement * GetElement() const
Get a pointer to the element this material is made of.
virtual Bool_t IsEq(const TGeoMaterial *other) const
return true if the other material has the same physical properties
virtual Double_t GetDensity() const
virtual Double_t GetZ() const
Geometrical transformation package.
Definition TGeoMatrix.h:38
@ kGeoSavePrimitive
Definition TGeoMatrix.h:48
virtual void ReflectZ(Bool_t leftside, Bool_t rotonly=kFALSE)
Multiply by a reflection respect to XY.
Bool_t IsReflection() const
Definition TGeoMatrix.h:66
virtual void RegisterYourself()
Register the matrix in the current manager, which will become the owner.
virtual Int_t GetByteCount() const
Get total size in bytes of this.
Media are used to store properties related to tracking and which are useful only when using geometry ...
Definition TGeoMedium.h:23
Int_t GetId() const
Definition TGeoMedium.h:45
@ kMedSavePrimitive
Definition TGeoMedium.h:25
TGeoMaterial * GetMaterial() const
Definition TGeoMedium.h:49
virtual Int_t GetByteCount() const
Definition TGeoMedium.h:44
Mixtures of elements.
Int_t GetNelements() const override
TGeoElement * GetElement(Int_t i=0) const override
Retrieve the pointer to the element corresponding to component I.
TGeoNavigator * AddNavigator()
Add a new navigator to the array.
TGeoNavigator * GetCurrentNavigator() const
TGeoNavigator * SetCurrentNavigator(Int_t inav)
Class providing navigation API for TGeo geometries.
void CdUp()
Go one level up in geometry.
void DoBackupState()
Backup the current state without affecting the cache stack.
void DoRestoreState()
Restore a backed-up state without affecting the cache stack.
TGeoNode * CrossBoundaryAndLocate(Bool_t downwards, TGeoNode *skipnode)
Cross next boundary and locate within current node The current point must be on the boundary of fCurr...
TGeoHMatrix * GetHMatrix()
Return stored current matrix (global matrix of the next touched node).
void LocalToMaster(const Double_t *local, Double_t *master) const
void CdNext()
Do a cd to the node found next by FindNextBoundary.
Double_t Safety(Bool_t inside=kFALSE)
Compute safe distance from the current point.
Bool_t GotoSafeLevel()
Go upwards the tree until a non-overlapping node.
Bool_t cd(const char *path="")
Browse the tree of nodes starting from top node according to pathname.
Bool_t IsSameLocation(Double_t x, Double_t y, Double_t z, Bool_t change=kFALSE)
Checks if point (x,y,z) is still in the current node.
void MasterToLocal(const Double_t *master, Double_t *local) const
Int_t GetVirtualLevel()
Find level of virtuality of current overlapping node (number of levels up having the same tracking me...
TGeoNode * InitTrack(const Double_t *point, const Double_t *dir)
Initialize current point and current direction vector (normalized) in MARS.
void InspectState() const
Inspects path and all flags for the current state.
TGeoNode * Step(Bool_t is_geom=kTRUE, Bool_t cross=kTRUE)
Make a rectiliniar step of length fStep from current point (fPoint) on current direction (fDirection)...
TGeoVolume * GetCurrentVolume() const
void ResetState()
Reset current state flags.
TGeoNode * FindNextDaughterBoundary(Double_t *point, Double_t *dir, Int_t &idaughter, Bool_t compmatrix=kFALSE)
Computes as fStep the distance to next daughter of the current volume.
void GetBranchNumbers(Int_t *copyNumbers, Int_t *volumeNumbers) const
Fill node copy numbers of current branch into an array.
Bool_t CheckPath(const char *path) const
Check if a geometry path is valid without changing the state of the navigator.
TGeoNode * FindNextBoundary(Double_t stepmax=TGeoShape::Big(), const char *path="", Bool_t frombdr=kFALSE)
Find distance to next boundary and store it in fStep.
TGeoNode * FindNode(Bool_t safe_start=kTRUE)
Returns deepest node containing current point.
TGeoNode * FindNextBoundaryAndStep(Double_t stepmax=TGeoShape::Big(), Bool_t compsafe=kFALSE)
Compute distance to next boundary within STEPMAX.
void CdTop()
Make top level node the current node.
Int_t GetCurrentNodeId() const
Double_t * FindNormalFast()
Computes fast normal to next crossed boundary, assuming that the current point is close enough to the...
void GetBranchOnlys(Int_t *isonly) const
Fill node copy numbers of current branch into an array.
TGeoNode * SearchNode(Bool_t downwards=kFALSE, const TGeoNode *skipnode=nullptr)
Returns the deepest node containing fPoint, which must be set a priori.
void CdNode(Int_t nodeid)
Change current path to point to the node having this id.
TGeoNodeCache * GetCache() const
void ResetAll()
Reset the navigator.
Bool_t IsSamePoint(Double_t x, Double_t y, Double_t z) const
Check if a new point with given coordinates is the same as the last located one.
void CdDown(Int_t index)
Make a daughter of current node current.
const char * GetPath() const
Get path to the current node in the form /node0/node1/...
Int_t GetSafeLevel() const
Go upwards the tree until a non-overlapping node.
Double_t * FindNormal(Bool_t forward=kTRUE)
Computes normal vector to the next surface that will be or was already crossed when propagating on a ...
void GetBranchNames(Int_t *names) const
Fill volume names of current branch into an array.
void BuildIdArray()
Builds node id array.
void BuildInfoBranch()
Builds info branch. Navigation is possible only after this step.
A node containing local transformation.
Definition TGeoNode.h:154
void SetMatrix(const TGeoMatrix *matrix)
Matrix setter.
Definition TGeoNode.cxx:827
A node represent a volume positioned inside another.They store links to both volumes and to the TGeoM...
Definition TGeoNode.h:39
Bool_t IsOverlapping() const
Definition TGeoNode.h:107
TGeoVolume * GetVolume() const
Definition TGeoNode.h:99
void SaveAttributes(std::ostream &out)
save attributes for this node
Definition TGeoNode.cxx:439
void SetVolume(TGeoVolume *volume)
Definition TGeoNode.h:117
void CheckShapes()
check for wrong parameters in shapes
Definition TGeoNode.cxx:329
void SetOverlapping(Bool_t flag=kTRUE)
Definition TGeoNode.h:120
Int_t GetNdaughters() const
Definition TGeoNode.h:91
virtual TGeoMatrix * GetMatrix() const =0
void SetVisibility(Bool_t vis=kTRUE) override
Set visibility of the node (obsolete).
Definition TGeoNode.cxx:716
void SetMotherVolume(TGeoVolume *mother)
Definition TGeoNode.h:125
static TClass * Class()
TGeoVolume * GetMotherVolume() const
Definition TGeoNode.h:90
void SetNumber(Int_t number)
Definition TGeoNode.h:118
void CheckOverlaps(Double_t ovlp=0.1, Option_t *option="")
Check overlaps bigger than OVLP hierarchically, starting with this node.
Definition TGeoNode.cxx:192
This is a wrapper class to G4OpticalSurface.
The knowledge of the path to the objects that need to be misaligned is essential since there is no ot...
void SetPhysicalNode(TGeoPhysicalNode *node)
Setter for the corresponding physical node.
Base class for a flat parallel geometry.
Bool_t CloseGeometry()
The main geometry must be closed.
void RefreshPhysicalNodes()
Refresh the node pointers and re-voxelize.
Bool_t IsClosed() const
Physical nodes are the actual 'touchable' objects in the geometry, representing a path of positioned ...
void Refresh()
Refresh this physical node.
Regions are groups of volumes having a common set of user tracking cuts.
Definition TGeoRegion.h:36
Base abstract class for all shapes.
Definition TGeoShape.h:25
virtual Bool_t IsComposite() const
Definition TGeoShape.h:130
Bool_t IsRunTimeShape() const
Definition TGeoShape.h:142
virtual void ComputeBBox()=0
virtual void AfterStreamer()
Definition TGeoShape.h:93
@ kGeoClosedShape
Definition TGeoShape.h:59
TClass * IsA() const override
Definition TGeoShape.h:171
Bool_t TestShapeBit(UInt_t f) const
Definition TGeoShape.h:167
Volume assemblies.
Definition TGeoVolume.h:316
static TGeoVolumeAssembly * MakeAssemblyFromVolume(TGeoVolume *vol)
Make a clone of volume VOL but which is an assembly.
Volume families.
Definition TGeoVolume.h:266
TGeoVolume, TGeoVolumeMulti, TGeoVolumeAssembly are the volume classes.
Definition TGeoVolume.h:43
Double_t WeightA() const
Analytical computation of the weight.
virtual void ClearThreadData() const
void SetVisibility(Bool_t vis=kTRUE) override
set visibility of this volume
void SaveAs(const char *filename, Option_t *option="") const override
Save geometry having this as top volume as a C++ macro.
void SetNumber(Int_t number)
Definition TGeoVolume.h:245
void SetLineWidth(Width_t lwidth) override
Set the line width.
TGeoMedium * GetMedium() const
Definition TGeoVolume.h:175
Int_t GetRefCount() const
Definition TGeoVolume.h:131
void SortNodes()
sort nodes by decreasing volume of the bounding box.
void Voxelize(Option_t *option)
build the voxels for this volume
Bool_t IsRunTime() const
Definition TGeoVolume.h:109
virtual void CreateThreadData(Int_t nthreads)
virtual Int_t GetByteCount() const
get the total size in bytes for this volume
Bool_t OptimizeVoxels()
Perform an extensive sampling to find which type of voxelization is most efficient.
virtual Bool_t IsVolumeMulti() const
Definition TGeoVolume.h:110
Int_t CountNodes(Int_t nlevels=1000, Int_t option=0)
Count total number of subnodes starting from this volume, nlevels down.
void UnmarkSaved()
Reset SavePrimitive bits.
void SetFinder(TGeoPatternFinder *finder)
Definition TGeoVolume.h:244
Int_t GetNdaughters() const
Definition TGeoVolume.h:362
void Grab()
Definition TGeoVolume.h:136
static TClass * Class()
void SetTransparency(Char_t transparency=0)
Definition TGeoVolume.h:376
void Release()
Definition TGeoVolume.h:137
void FindOverlaps() const
loop all nodes marked as overlaps and find overlapping brothers
TGeoNode * GetNode(const char *name) const
get the pointer to a daughter node
virtual void SetMedium(TGeoMedium *medium)
Definition TGeoVolume.h:242
TGeoVoxelFinder * GetVoxels() const
Getter for optimization structure.
static TGeoMedium * DummyMedium()
void SetLineColor(Color_t lcolor) override
Set the line color.
Int_t GetNumber() const
Definition TGeoVolume.h:184
TGeoShape * GetShape() const
Definition TGeoVolume.h:190
void SetField(TObject *field)
Definition TGeoVolume.h:231
static void CreateDummyMedium()
Create a dummy medium.
void SetLineStyle(Style_t lstyle) override
Set the line style.
virtual Bool_t IsAssembly() const
Returns true if the volume is an assembly or a scaled assembly.
TGeoVolume * MakeReflectedVolume(const char *newname="") const
Make a copy of this volume which is reflected with respect to XY plane.
virtual Bool_t IsVisible() const
Definition TGeoVolume.h:155
Finder class handling voxels.
void SetNeedRebuild(Bool_t flag=kTRUE)
A TGeoXtru shape is represented by the extrusion of an arbitrary polygon with fixed outline between s...
Definition TGeoXtru.h:22
Int_t GetNvert() const
Definition TGeoXtru.h:96
static TClass * Class()
THashList implements a hybrid collection class consisting of a hash table and a list to store TObject...
Definition THashList.h:34
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
void Clear(Option_t *option="") override
Remove all objects from the list.
TObject * FindObject(const char *name) const override
Find object using its name.
void AddLast(TObject *obj) override
Add object at the end of the list.
Definition THashList.cxx:95
Book space in a file, create I/O buffers, to fill them, (un)compress them.
Definition TKey.h:28
virtual const char * GetClassName() const
Definition TKey.h:75
virtual TObject * ReadObj()
To read a TObject* from the file.
Definition TKey.cxx:758
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:578
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:822
void Delete(Option_t *option="") override
Remove all objects from the list AND delete all heap based objects.
Definition TList.cxx:470
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
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
TNamed()
Definition TNamed.h:36
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
An array of TObjects.
Definition TObjArray.h:31
Int_t GetEntriesFast() const
Definition TObjArray.h:58
Int_t IndexOf(const TObject *obj) const override
void AddAt(TObject *obj, Int_t idx) override
Add object at position ids.
virtual void Sort(Int_t upto=kMaxInt)
If objects in array are sortable (i.e.
void Clear(Option_t *option="") override
Remove all objects from the array.
virtual void AddAtAndExpand(TObject *obj, Int_t idx)
Add object at position idx.
Int_t GetEntries() const override
Return the number of objects in array (i.e.
void Delete(Option_t *option="") override
Remove all objects from the array AND delete all heap based objects.
TObject * At(Int_t idx) const override
Definition TObjArray.h:164
TObject * UncheckedAt(Int_t i) const
Definition TObjArray.h:84
TObject * Remove(TObject *obj) override
Remove object from array.
TObject * FindObject(const char *name) const override
Find an object in this collection using its name.
void Add(TObject *obj) override
Definition TObjArray.h:68
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:201
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:962
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:184
virtual Int_t Write(const char *name=nullptr, Int_t option=0, Int_t bufsize=0)
Write this object to the current directory.
Definition TObject.cxx:880
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:976
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1004
virtual void SetUniqueID(UInt_t uid)
Set the unique object id.
Definition TObject.cxx:791
virtual TClass * IsA() const
Definition TObject.h:245
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:274
void ResetBit(UInt_t f)
Definition TObject.h:200
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:950
Bool_t Connect(const char *signal, const char *receiver_class, void *receiver, const char *slot)
Non-static method is used to connect from the signal of this object to the receiver slot.
Definition TQObject.cxx:869
Sequenceable collection abstract base class.
virtual Int_t IndexOf(const TObject *obj) const
Return index of object in collection.
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:419
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
TSubString Strip(EStripType s=kTrailing, char c=' ') const
Return a substring of self stripped at beginning and/or end.
Definition TString.cxx:1163
const char * Data() const
Definition TString.h:378
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:706
Bool_t IsNull() const
Definition TString.h:416
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:2378
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:634
Abstract class for geometry painters.
virtual TGeoNode * SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil, const char *g3path)=0
virtual void SetTopVisible(Bool_t vis=kTRUE)=0
virtual Int_t GetVisLevel() const =0
virtual Double_t Weight(Double_t precision, Option_t *option="v")=0
virtual void CheckGeometryFull(Bool_t checkoverlaps=kTRUE, Bool_t checkcrossings=kTRUE, Int_t nrays=10000, const Double_t *vertex=nullptr)=0
virtual void DrawPath(const char *path, Option_t *option="")=0
virtual void ModifiedPad(Bool_t update=kFALSE) const =0
virtual void SetCheckedNode(TGeoNode *node)=0
virtual void GetViewAngles(Double_t &, Double_t &, Double_t &)
virtual TVirtualGeoTrack * AddTrack(Int_t id, Int_t pdgcode, TObject *particle)=0
virtual void SetExplodedView(Int_t iopt=0)=0
virtual Bool_t IsRaytracing() const =0
virtual void DrawCurrentPoint(Int_t color)=0
virtual void SetClippingShape(TGeoShape *shape)=0
virtual void TestOverlaps(const char *path)=0
virtual void PrintOverlaps() const =0
virtual void GrabFocus(Int_t nfr=0, Double_t dlong=0, Double_t dlat=0, Double_t dpsi=0)=0
virtual void Test(Int_t npoints, Option_t *option)=0
virtual void SetVisOption(Int_t option=0)=0
virtual Double_t * GetViewBox()=0
virtual void EstimateCameraMove(Double_t, Double_t, Double_t *, Double_t *)
virtual void CheckBoundaryReference(Int_t icheck=-1)=0
virtual void SetNsegments(Int_t nseg=20)=0
virtual void RandomRays(Int_t nrays, Double_t startx, Double_t starty, Double_t startz, const char *target_vol, Bool_t check_norm)=0
virtual void RandomPoints(const TGeoVolume *vol, Int_t npoints, Option_t *option="")=0
virtual Int_t CountVisibleNodes()=0
virtual void DefaultAngles()=0
virtual void SetNmeshPoints(Int_t npoints)=0
virtual void UnbombTranslation(const Double_t *tr, Double_t *bombtr)=0
virtual void EditGeometry(Option_t *option="")=0
virtual void GetBombFactors(Double_t &bombx, Double_t &bomby, Double_t &bombz, Double_t &bombr) const =0
virtual void CheckPoint(Double_t x=0, Double_t y=0, Double_t z=0, Option_t *option="", Double_t safety=0.)=0
virtual void BombTranslation(const Double_t *tr, Double_t *bombtr)=0
virtual void CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option)=0
virtual void ExecuteManagerEvent(TGeoManager *geom, Int_t event, Int_t px, Int_t py)=0
virtual void SetBombFactors(Double_t bombx=1.3, Double_t bomby=1.3, Double_t bombz=1.3, Double_t bombr=1.3)=0
virtual void CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.)=0
Base class for user-defined tracks attached to a geometry.
Int_t GetId() const
TVirtualGeoTrack * GetMother() const
virtual TVirtualGeoTrack * FindTrackWithId(Int_t id) const
Recursively search through this track for a daughter particle (at any depth) with the specified id.
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition fillpatterns.C:1
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
void EnableThreadSafety()
Enable support for multi-threading within the ROOT code in particular, enables the global mutex to ma...
Definition TROOT.cxx:499
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:646
Long64_t BinarySearch(Long64_t n, const T *array, T value)
Binary search in an array of n values to locate value.
Definition TMathBase.h:347
constexpr Double_t RadToDeg()
Conversion from radian to degree: .
Definition TMath.h:72