ROOT logo
// @(#)root/geom:$Id$
// Author: Andrei Gheata   25/10/01

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#ifndef ROOT_TGeoManager
#define ROOT_TGeoManager

#ifndef ROOT_TObjArray
#include "TObjArray.h"
#endif
#ifndef ROOT_TGeoNavigator
#include "TGeoNavigator.h"
#endif

// forward declarations
class TVirtualGeoTrack;
class TGeoNode;
class TGeoPhysicalNode;
class TGeoPNEntry;
class TGeoVolume;
class TGeoVolumeMulti;
class TGeoMatrix;
class TGeoHMatrix;
class TGeoMaterial;
class TGeoMedium;
class TGeoShape;
class TVirtualGeoPainter;
class THashList;
class TGeoParallelWorld;

////////////////////////////////////////////////////////////////////////////
//                                                                        //
// TGeoManager - The manager class for any TGeo geometry. Provides user   //
//    interface for geometry creation, navigation, state querying,        //
//    visualization, IO, geometry checking and other utilities.           //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

class TGeoManager : public TNamed
{
protected:
   static Bool_t         fgLock;            //! Lock preventing a second geometry to be loaded
   static Int_t          fgVerboseLevel;    //! Verbosity level for Info messages (no IO).
   static Int_t          fgMaxLevel;        //! Maximum level in geometry
   static Int_t          fgMaxDaughters;    //! Maximum number of daughters
   static Int_t          fgMaxXtruVert;     //! Maximum number of Xtru vertices

   TGeoManager(const TGeoManager&); 
   TGeoManager& operator=(const TGeoManager&); 

private :
   Double_t              fPhimin;           //! lowest range for phi cut
   Double_t              fPhimax;           //! highest range for phi cut
   Double_t              fTmin;             //! lower time limit for tracks drawing
   Double_t              fTmax;             //! upper time limit for tracks drawing
   Int_t                 fNNodes;           // total number of physical nodes
   TString               fPath;             //! path to current node
   TString               fParticleName;     //! particles to be drawn
   Double_t              fVisDensity;       // transparency threshold by density
   Int_t                 fExplodedView;     // exploded view mode
   Int_t                 fVisOption;        // global visualization option
   Int_t                 fVisLevel;         // maximum visualization depth
   Int_t                 fNsegments;        // number of segments to approximate circles
   Int_t                 fNtracks;          // number of tracks
   Int_t                 fMaxVisNodes;      // maximum number of visible nodes
   TVirtualGeoTrack     *fCurrentTrack;     //! current track
   Int_t                 fNpdg;             // number of different pdg's stored
   Int_t                 fPdgId[1024];      // pdg conversion table
   Bool_t                fClosed;           //! flag that geometry is closed
   Bool_t                fLoopVolumes;      //! flag volume lists loop
   Bool_t                fStreamVoxels;     // flag to allow voxelization I/O
   Bool_t                fIsGeomReading;    //! flag set when reading geometry
   Bool_t                fIsGeomCleaning;   //! flag to notify that the manager is being destructed
   Bool_t                fPhiCut;           // flag for phi cuts
   Bool_t                fTimeCut;          // time cut for tracks
   Bool_t                fDrawExtra;        //! flag that the list of physical nodes has to be drawn
   Bool_t                fMatrixTransform;  //! flag for using GL matrix
   Bool_t                fMatrixReflection; //! flag for GL reflections
   Bool_t                fActivity;         //! switch ON/OFF volume activity (default OFF - all volumes active))
   Bool_t                fIsNodeSelectable; //! flag that nodes are the selected objects in pad rather than volumes
   TVirtualGeoPainter   *fPainter;          //! current painter

   TObjArray            *fMatrices;         //-> list of local transformations
   TObjArray            *fShapes;           //-> list of shapes
   TObjArray            *fVolumes;          //-> list of volumes
   TObjArray            *fPhysicalNodes;    //-> list of physical nodes
   TObjArray            *fGShapes;          //! list of runtime shapes
   TObjArray            *fGVolumes;         //! list of runtime volumes
   TObjArray            *fTracks;           //-> list of tracks attached to geometry
   TObjArray            *fPdgNames;         //-> list of pdg names for tracks
//   TObjArray            *fNavigators;       //! list of navigators
   TList                *fMaterials;        //-> list of materials
   TList                *fMedia;            //-> list of tracking media
   TObjArray            *fNodes;            //-> current branch of nodes
   TObjArray            *fOverlaps;         //-> list of geometrical overlaps
   UChar_t              *fBits;             //! bits used for voxelization
   // Map of navigatorr arrays per thread
   typedef std::map<Long_t, TGeoNavigatorArray *>   NavigatorsMap_t;
   typedef NavigatorsMap_t::iterator                NavigatorsMapIt_t;
   typedef std::map<Long_t, Int_t>                  ThreadsMap_t;
   typedef ThreadsMap_t::const_iterator             ThreadsMapIt_t;
   
   NavigatorsMap_t       fNavigators;       //! Map between thread id's and navigator arrays
   static ThreadsMap_t  *fgThreadId;        //! Thread id's map
   static Int_t          fgNumThreads;      //! Number of registered threads
   static Bool_t         fgLockNavigators;   //! Lock existing navigators
   TGeoNavigator        *fCurrentNavigator; //! current navigator
   TGeoVolume           *fCurrentVolume;    //! current volume
   TGeoVolume           *fTopVolume;        //! top level volume in geometry
   TGeoNode             *fTopNode;          //! top physical node
   TGeoVolume           *fMasterVolume;     // master volume
   TGeoHMatrix          *fGLMatrix;         // matrix to be used for view transformations
   TObjArray            *fUniqueVolumes;    //-> list of unique volumes
   TGeoShape            *fClippingShape;    //! clipping shape for raytracing
   TGeoElementTable     *fElementTable;     //! table of elements

   Int_t                *fNodeIdArray;      //! array of node id's
   Int_t                 fNLevel;           // maximum accepted level in geometry
   TGeoVolume           *fPaintVolume;      //! volume currently painted
   THashList            *fHashVolumes;      //! hash list of volumes providing fast search
   THashList            *fHashGVolumes;     //! hash list of group volumes providing fast search
   THashList            *fHashPNE;          //-> hash list of phisical node entries
   mutable TObjArray    *fArrayPNE;         //! array of phisical node entries
   Int_t                 fSizePNEId;        // size of the array of unique ID's for PN entries
   Int_t                 fNPNEId;           // number of PN entries having a unique ID
   Int_t                *fKeyPNEId;         //[fSizePNEId] array of uid values for PN entries
   Int_t                *fValuePNEId;       //[fSizePNEId] array of pointers to PN entries with ID's
   Int_t                 fMaxThreads;       //! Max number of threads
   Bool_t                fMultiThread;      //! Flag for multi-threading
   Bool_t                fUsePWNav;         // Activate usage of parallel world in navigation
   TGeoParallelWorld    *fParallelWorld;    // Parallel world
//--- private methods
   Bool_t                IsLoopingVolumes() const     {return fLoopVolumes;}
   void                  Init();
   Bool_t                InitArrayPNE() const;
   Bool_t                InsertPNEId(Int_t uid, Int_t ientry);
   void                  SetLoopVolumes(Bool_t flag=kTRUE) {fLoopVolumes=flag;}
   void                  UpdateElements();
   void                  Voxelize(Option_t *option = 0);

public:
   // constructors
   TGeoManager();
   TGeoManager(const char *name, const char *title);
   // destructor
   virtual ~TGeoManager();
   //--- adding geometrical objects
   Int_t                  AddMaterial(const TGeoMaterial *material);
   Int_t                  AddOverlap(const TNamed *ovlp);
   Int_t                  AddTransformation(const TGeoMatrix *matrix);
   Int_t                  AddShape(const TGeoShape *shape);
   Int_t                  AddTrack(Int_t id, Int_t pdgcode, TObject *particle=0);
   Int_t                  AddTrack(TVirtualGeoTrack *track);
   Int_t                  AddVolume(TGeoVolume *volume);
   TGeoNavigator         *AddNavigator();
   void                   ClearOverlaps();
   void                   RegisterMatrix(const TGeoMatrix *matrix);
   void                   SortOverlaps();
   //--- browsing and tree navigation
   void                   Browse(TBrowser *b);
   void                   SetVisibility(TObject *obj, Bool_t vis);
   virtual Bool_t         cd(const char *path=""); // *MENU*
   Bool_t                 CheckPath(const char *path) const;
   void                   CdNode(Int_t nodeid);
   void                   CdDown(Int_t index);
   void                   CdUp();
   void                   CdTop();
   void                   CdNext();
   void                   GetBranchNames(Int_t *names) const;
   void                   GetBranchNumbers(Int_t *copyNumbers, Int_t *volumeNumbers) const;
   void                   GetBranchOnlys(Int_t *isonly) const;
   Int_t                  GetNmany() const {return GetCurrentNavigator()->GetNmany();}
   const char            *GetPdgName(Int_t pdg) const;
   void                   SetPdgName(Int_t pdg, const char *name);
   Bool_t                 IsFolder() const { return kTRUE; }
   //--- visualization settings
   virtual void           Edit(Option_t *option="");  // *MENU*
   void                   BombTranslation(const Double_t *tr, Double_t *bombtr);
   void                   UnbombTranslation(const Double_t *tr, Double_t *bombtr);
   void                   ClearAttributes(); // *MENU*
   void                   DefaultAngles();   // *MENU*
   void                   DefaultColors();   // *MENU*
   TGeoShape             *GetClippingShape() const {return fClippingShape;}
   Int_t                  GetNsegments() const;
   TVirtualGeoPainter    *GetGeomPainter();
   TVirtualGeoPainter    *GetPainter() const {return fPainter;}
   Int_t                  GetBombMode() const  {return fExplodedView;}
   void                   GetBombFactors(Double_t &bombx, Double_t &bomby, Double_t &bombz, Double_t &bombr) const;
   Int_t                  GetMaxVisNodes() const {return fMaxVisNodes;}
   Bool_t                 GetTminTmax(Double_t &tmin, Double_t &tmax) const;
   Double_t               GetTmax() const {return fTmax;}
   TGeoVolume            *GetPaintVolume() const {return fPaintVolume;}
   Double_t               GetVisDensity() const  {return fVisDensity;}
   Int_t                  GetVisLevel() const;
   Int_t                  GetVisOption() const;
   Bool_t                 IsInPhiRange() const;
   Bool_t                 IsDrawingExtra() const {return fDrawExtra;}
   Bool_t                 IsNodeSelectable() const {return fIsNodeSelectable;}
   Bool_t                 IsVisLeaves() const {return (fVisOption==1)?kTRUE:kFALSE;}
   void                   ModifiedPad() const;
   void                   OptimizeVoxels(const char *filename="tgeovox.C"); // *MENU*
   void                   SetClipping(Bool_t flag=kTRUE) {SetClippingShape(((flag)?fClippingShape:0));} // *MENU*
   void                   SetClippingShape(TGeoShape *clip);
   void                   SetExplodedView(Int_t iopt=0); // *MENU*
   void                   SetPhiRange(Double_t phimin=0., Double_t phimax=360.);
   void                   SetNsegments(Int_t nseg); // *MENU*
   Bool_t                 SetCurrentNavigator(Int_t index);
   void                   SetBombFactors(Double_t bombx=1.3, Double_t bomby=1.3, Double_t bombz=1.3, Double_t bombr=1.3); // *MENU*
   void                   SetPaintVolume(TGeoVolume *vol) {fPaintVolume = vol;}
   void                   SetTopVisible(Bool_t vis=kTRUE);
   void                   SetTminTmax(Double_t tmin=0, Double_t tmax=999);
   void                   SetDrawExtraPaths(Bool_t flag=kTRUE) {fDrawExtra=flag;}
   void                   SetNodeSelectable(Bool_t flag=kTRUE) {fIsNodeSelectable=flag;}
   void                   SetVisDensity(Double_t dens=0.01); // *MENU*
   void                   SetVisLevel(Int_t level=3);   // *MENU*
   void                   SetVisOption(Int_t option=0);
   void                   ViewLeaves(Bool_t flag=kTRUE); // *TOGGLE* *GETTER=IsVisLeaves
   void                   SaveAttributes(const char *filename="tgeoatt.C"); // *MENU*
   void                   RestoreMasterVolume(); // *MENU*
   void                   SetMaxVisNodes(Int_t maxnodes=10000); // *MENU*
   //--- geometry checking
   void                   AnimateTracks(Double_t tmin=0, Double_t tmax=5E-8, Int_t nframes=200, Option_t *option="/*"); // *MENU*
   void                   CheckBoundaryErrors(Int_t ntracks=1000000, Double_t radius=-1.); // *MENU*
   void                   CheckBoundaryReference(Int_t icheck=-1);
   void                   CheckGeometryFull(Int_t ntracks=1000000, Double_t vx=0., Double_t vy=0., Double_t vz=0., Option_t *option="ob"); // *MENU*
   void                   CheckGeometry(Option_t *option="");
   void                   CheckOverlaps(Double_t ovlp=0.1, Option_t *option=""); // *MENU*
   void                   CheckPoint(Double_t x=0,Double_t y=0, Double_t z=0, Option_t *option=""); // *MENU*
   void                   CheckShape(TGeoShape *shape, Int_t testNo, Int_t nsamples, Option_t *option);
   void                   ConvertReflections();
   void                   DrawCurrentPoint(Int_t color=2); // *MENU*
   void                   DrawTracks(Option_t *option=""); // *MENU*
   void                   SetParticleName(const char *pname) {fParticleName=pname;}
   const char            *GetParticleName() const {return fParticleName.Data();}
   void                   DrawPath(const char *path, Option_t *option="");
   void                   PrintOverlaps() const; // *MENU*
   void                   RandomPoints(const TGeoVolume *vol, Int_t npoints=10000, Option_t *option="");
   void                   RandomRays(Int_t nrays=1000, Double_t startx=0, Double_t starty=0, Double_t startz=0, const char *target_vol=0, Bool_t check_norm=kFALSE);
   TGeoNode              *SamplePoints(Int_t npoints, Double_t &dist, Double_t epsil=1E-5,
                                       const char *g3path="");
   void                   SetNmeshPoints(Int_t npoints=1000);
   void                   SetCheckedNode(TGeoNode *node);                                    
   void                   Test(Int_t npoints=1000000, Option_t *option=""); // *MENU*
   void                   TestOverlaps(const char* path=""); // *MENU*
   Double_t               Weight(Double_t precision=0.01, Option_t *option="va"); // *MENU*

   //--- GEANT3-like geometry creation
   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="");
   void                   Matrix(Int_t index, Double_t theta1, Double_t phi1,
                                       Double_t theta2, Double_t phi2,
                                       Double_t theta3, Double_t phi3);
   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);
   TGeoMaterial          *Mixture(const char *name, Float_t *a, Float_t *z, Double_t dens,
                                        Int_t nelem, Float_t *wmat, Int_t uid);
   TGeoMaterial          *Mixture(const char *name, Double_t *a, Double_t *z, Double_t dens,
                                        Int_t nelem, Double_t *wmat, Int_t uid);
   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);
   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);
   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, Double_t *upar, Int_t npar=0);
   TGeoVolume            *Volume(const char *name, const char *shape, Int_t nmed,
                                       Float_t *upar, Int_t npar=0);
   TGeoVolume            *Volume(const char *name, const char *shape, Int_t nmed,
                                       Double_t *upar, Int_t npar=0);
   void                   SetVolumeAttribute(const char *name, const char *att, Int_t val);
   //--- geometry building
   void                   BuildDefaultMaterials();
   void                   CloseGeometry(Option_t *option="d");
   Bool_t                 IsClosed() const {return fClosed;}
   TGeoVolume            *MakeArb8(const char *name, TGeoMedium *medium,
                                     Double_t dz, Double_t *vertices=0);
   TGeoVolume            *MakeBox(const char *name, TGeoMedium *medium,
                                     Double_t dx, Double_t dy, Double_t dz);
   TGeoVolume            *MakeCone(const char *name, TGeoMedium *medium,
                                      Double_t dz, Double_t rmin1, Double_t rmax1,
                                      Double_t rmin2, Double_t rmax2);
   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);
   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);
   TGeoVolume            *MakeEltu(const char *name, TGeoMedium *medium,
                                      Double_t a, Double_t b, Double_t dz);
   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);
   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);
   TGeoVolume            *MakePcon(const char *name, TGeoMedium *medium,
                                      Double_t phi, Double_t dphi, Int_t nz);
   TGeoVolume            *MakeParaboloid(const char *name, TGeoMedium *medium,
                                      Double_t rlo, Double_t rhi, Double_t dz);
   TGeoVolume            *MakeHype(const char *name, TGeoMedium *medium,
                                      Double_t rin, Double_t stin, Double_t rout, Double_t stout, Double_t dz);
   TGeoVolume            *MakePgon(const char *name, TGeoMedium *medium,
                                      Double_t phi, Double_t dphi, Int_t nedges, Int_t nz);
   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);
   TGeoVolume            *MakeTorus(const char *name, TGeoMedium *medium, Double_t r,
                                    Double_t rmin, Double_t rmax, Double_t phi1=0, Double_t dphi=360);
   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);
   TGeoVolume            *MakeTrd1(const char *name, TGeoMedium *medium,
                                      Double_t dx1, Double_t dx2, Double_t dy, Double_t dz);
   TGeoVolume            *MakeTrd2(const char *name, TGeoMedium *medium,
                                      Double_t dx1, Double_t dx2, Double_t dy1, Double_t dy2,
                                      Double_t dz);
   TGeoVolume            *MakeTube(const char *name, TGeoMedium *medium,
                                      Double_t rmin, Double_t rmax, Double_t dz);
   TGeoVolume            *MakeTubs(const char *name, TGeoMedium *medium,
                                      Double_t rmin, Double_t rmax, Double_t dz,
                                      Double_t phi1, Double_t phi2);
   TGeoVolume            *MakeXtru(const char *name, TGeoMedium *medium,
                                   Int_t nz);

   TGeoPNEntry           *SetAlignableEntry(const char *unique_name, const char *path, Int_t uid=-1);
   TGeoPNEntry           *GetAlignableEntry(const char *name) const;   
   TGeoPNEntry           *GetAlignableEntry(Int_t index) const;
   TGeoPNEntry           *GetAlignableEntryByUID(Int_t uid) const;
   Int_t                  GetNAlignable(Bool_t with_uid=kFALSE) const;
   TGeoPhysicalNode      *MakeAlignablePN(const char *name);
   TGeoPhysicalNode      *MakeAlignablePN(TGeoPNEntry *entry);
   TGeoPhysicalNode      *MakePhysicalNode(const char *path=0);
   void                   ClearPhysicalNodes(Bool_t mustdelete=kFALSE);
   void                   RefreshPhysicalNodes(Bool_t lock=kTRUE);
   TVirtualGeoTrack      *MakeTrack(Int_t id, Int_t pdgcode, TObject *particle);
   TGeoVolumeAssembly    *MakeVolumeAssembly(const char *name);
   TGeoVolumeMulti       *MakeVolumeMulti(const char *name, TGeoMedium *medium);
   void                   SetTopVolume(TGeoVolume *vol);

   //--- geometry queries
   TGeoNode              *CrossBoundaryAndLocate(Bool_t downwards, TGeoNode *skipnode);
   TGeoNode              *FindNextBoundary(Double_t stepmax=TGeoShape::Big(),const char *path="", Bool_t frombdr=kFALSE);
   TGeoNode              *FindNextDaughterBoundary(Double_t *point, Double_t *dir, Int_t &idaughter, Bool_t compmatrix=kFALSE);
   TGeoNode              *FindNextBoundaryAndStep(Double_t stepmax=TGeoShape::Big(), Bool_t compsafe=kFALSE);
   TGeoNode              *FindNode(Bool_t safe_start=kTRUE);
   TGeoNode              *FindNode(Double_t x, Double_t y, Double_t z);
   Double_t              *FindNormal(Bool_t forward=kTRUE);
   Double_t              *FindNormalFast();
   TGeoNode              *InitTrack(const Double_t *point, const Double_t *dir);
   TGeoNode              *InitTrack(Double_t x, Double_t y, Double_t z, Double_t nx, Double_t ny, Double_t nz);
   void                   ResetState();
   Double_t               Safety(Bool_t inside=kFALSE);
   TGeoNode              *SearchNode(Bool_t downwards=kFALSE, const TGeoNode *skipnode=0);
   TGeoNode              *Step(Bool_t is_geom=kTRUE, Bool_t cross=kTRUE);
   void                   DisableInactiveVolumes() {fActivity=kTRUE;}
   void                   EnableInactiveVolumes()  {fActivity=kFALSE;}
   void                   SetCurrentTrack(Int_t i) {fCurrentTrack = (TVirtualGeoTrack*)fTracks->At(i);}
   void                   SetCurrentTrack(TVirtualGeoTrack *track) {fCurrentTrack=track;}
   Int_t                  GetNtracks() const {return fNtracks;}
   TVirtualGeoTrack      *GetCurrentTrack() {return fCurrentTrack;}
   TVirtualGeoTrack      *GetLastTrack() {return (TVirtualGeoTrack*)((fNtracks>0)?fTracks->At(fNtracks-1):NULL);}
   const Double_t        *GetLastPoint() const {return GetCurrentNavigator()->GetLastPoint();}
   TVirtualGeoTrack      *GetTrack(Int_t index)         {return (index<fNtracks)?(TVirtualGeoTrack*)fTracks->At(index):0;}
   Int_t                  GetTrackIndex(Int_t id) const;
   TVirtualGeoTrack      *GetTrackOfId(Int_t id) const;
   TVirtualGeoTrack      *FindTrackWithId(Int_t id) const;
   TVirtualGeoTrack      *GetParentTrackOfId(Int_t id) const;
   Int_t                  GetVirtualLevel();
   Bool_t                 GotoSafeLevel();
   Int_t                  GetSafeLevel() const;
   Double_t               GetSafeDistance() const      {return GetCurrentNavigator()->GetSafeDistance();}
   Double_t               GetLastSafety() const        {return GetCurrentNavigator()->GetLastSafety();}
   Double_t               GetStep() const              {return GetCurrentNavigator()->GetStep();}
   void                   InspectState() const;
   Bool_t                 IsAnimatingTracks() const    {return fIsGeomReading;}
   Bool_t                 IsCheckingOverlaps() const   {return GetCurrentNavigator()->IsCheckingOverlaps();}
   Bool_t                 IsMatrixTransform() const    {return fMatrixTransform;}
   Bool_t                 IsMatrixReflection() const   {return fMatrixReflection;}
   Bool_t                 IsSameLocation(Double_t x, Double_t y, Double_t z, Bool_t change=kFALSE);
   Bool_t                 IsSameLocation() const {return GetCurrentNavigator()->IsSameLocation();}
   Bool_t                 IsSamePoint(Double_t x, Double_t y, Double_t z) const;
   Bool_t                 IsStartSafe() const {return GetCurrentNavigator()->IsStartSafe();}
   void                   SetCheckingOverlaps(Bool_t flag=kTRUE) {GetCurrentNavigator()->SetCheckingOverlaps(flag);}
   void                   SetStartSafe(Bool_t flag=kTRUE)   {GetCurrentNavigator()->SetStartSafe(flag);}
   void                   SetMatrixTransform(Bool_t on=kTRUE) {fMatrixTransform = on;}
   void                   SetMatrixReflection(Bool_t flag=kTRUE) {fMatrixReflection = flag;}
   void                   SetStep(Double_t step) {GetCurrentNavigator()->SetStep(step);}
   Bool_t                 IsCurrentOverlapping() const {return GetCurrentNavigator()->IsCurrentOverlapping();}
   Bool_t                 IsEntering() const           {return GetCurrentNavigator()->IsEntering();}
   Bool_t                 IsExiting() const            {return GetCurrentNavigator()->IsExiting();}
   Bool_t                 IsStepEntering() const       {return GetCurrentNavigator()->IsStepEntering();}
   Bool_t                 IsStepExiting() const        {return GetCurrentNavigator()->IsStepExiting();}
   Bool_t                 IsOutside() const            {return GetCurrentNavigator()->IsOutside();}
   Bool_t                 IsOnBoundary() const         {return GetCurrentNavigator()->IsOnBoundary();}
   Bool_t                 IsNullStep() const           {return GetCurrentNavigator()->IsNullStep();}
   Bool_t                 IsActivityEnabled() const    {return fActivity;}
   void                   SetOutside(Bool_t flag=kTRUE) {GetCurrentNavigator()->SetOutside(flag);}
   

   //--- cleaning
   void                   CleanGarbage();
   void                   ClearShape(const TGeoShape *shape);
   void                   ClearTracks() {fTracks->Delete(); fNtracks=0;}
   void                   ClearNavigators();
   void                   RemoveMaterial(Int_t index);
   void                   RemoveNavigator(const TGeoNavigator *nav);
   void                   ResetUserData();


   //--- utilities
   Int_t                  CountNodes(const TGeoVolume *vol=0, Int_t nlevels=10000, Int_t option=0);
   void                   CountLevels();
   virtual void           ExecuteEvent(Int_t event, Int_t px, Int_t py);
   static Int_t           Parse(const char* expr, TString &expr1, TString &expr2, TString &expr3);
   Int_t                  ReplaceVolume(TGeoVolume *vorig, TGeoVolume *vnew);
   Int_t                  TransformVolumeToAssembly(const char *vname);
   UChar_t               *GetBits() {return fBits;}
   virtual Int_t          GetByteCount(Option_t *option=0);
   void                   SetAllIndex();
   static Int_t           GetMaxDaughters();
   static Int_t           GetMaxLevels();
   static Int_t           GetMaxXtruVert();
   Int_t                  GetMaxThreads() const {return fMaxThreads-1;}
   void                   SetMaxThreads(Int_t nthreads);
   void                   SetMultiThread(Bool_t flag=kTRUE) {fMultiThread = flag;}
   Bool_t                 IsMultiThread() const {return fMultiThread;}
   static void            SetNavigatorsLock(Bool_t flag);
   static Int_t           ThreadId();
   static Int_t           GetNumThreads();
   static void            ClearThreadsMap();
   void                   ClearThreadData() const;
   void                   CreateThreadData() const;

   //--- I/O
   virtual Int_t          Export(const char *filename, const char *name="", Option_t *option="vg");
   static  void           LockGeometry();
   static  void           UnlockGeometry();
   static  Int_t          GetVerboseLevel();
   static  void           SetVerboseLevel(Int_t vl);
   static TGeoManager    *Import(const char *filename, const char *name="", Option_t *option="");
   static Bool_t          IsLocked();
   Bool_t                 IsStreamingVoxels() const {return fStreamVoxels;}
   Bool_t                 IsCleaning() const {return fIsGeomCleaning;}

   //--- list getters
   TObjArray             *GetListOfNodes()              {return fNodes;}
   TObjArray             *GetListOfPhysicalNodes()      {return fPhysicalNodes;}
   TObjArray             *GetListOfOverlaps()           {return fOverlaps;}
   TObjArray             *GetListOfMatrices() const     {return fMatrices;}
   TList                 *GetListOfMaterials() const    {return fMaterials;}
   TList                 *GetListOfMedia() const        {return fMedia;}
   TObjArray             *GetListOfVolumes() const      {return fVolumes;}
   TObjArray             *GetListOfGVolumes() const     {return fGVolumes;}
   TObjArray             *GetListOfShapes() const       {return fShapes;}
   TObjArray             *GetListOfGShapes() const      {return fGShapes;}
   TObjArray             *GetListOfUVolumes() const     {return fUniqueVolumes;}
   TObjArray             *GetListOfTracks() const       {return fTracks;}
   TGeoNavigatorArray    *GetListOfNavigators() const;
   TGeoElementTable      *GetElementTable();

   //--- modeler state getters/setters
   void                   DoBackupState();
   void                   DoRestoreState();
   TGeoNode              *GetNode(Int_t level) const  {return (TGeoNode*)fNodes->UncheckedAt(level);}
   Int_t                  GetNodeId() const           {return GetCurrentNavigator()->GetNodeId();}
   TGeoNode              *GetNextNode() const         {return GetCurrentNavigator()->GetNextNode();}
   TGeoNode              *GetMother(Int_t up=1) const {return GetCurrentNavigator()->GetMother(up);}
   TGeoHMatrix           *GetMotherMatrix(Int_t up=1) const {return GetCurrentNavigator()->GetMotherMatrix(up);}
   TGeoHMatrix           *GetHMatrix();
   TGeoHMatrix           *GetCurrentMatrix() const    {return GetCurrentNavigator()->GetCurrentMatrix();}
   TGeoHMatrix           *GetGLMatrix() const         {return fGLMatrix;}
   TGeoNavigator         *GetCurrentNavigator() const;
   TGeoNode              *GetCurrentNode() const      {return GetCurrentNavigator()->GetCurrentNode();}
   Int_t                  GetCurrentNodeId() const;
   const Double_t        *GetCurrentPoint() const     {return GetCurrentNavigator()->GetCurrentPoint();}
   const Double_t        *GetCurrentDirection() const {return GetCurrentNavigator()->GetCurrentDirection();}
   TGeoVolume            *GetCurrentVolume() const {return GetCurrentNavigator()->GetCurrentVolume();}
   const Double_t        *GetCldirChecked() const  {return GetCurrentNavigator()->GetCldirChecked();}
   const Double_t        *GetCldir() const         {return GetCurrentNavigator()->GetCldir();}
   const Double_t        *GetNormal() const        {return GetCurrentNavigator()->GetNormal();}
   Int_t                  GetLevel() const         {return GetCurrentNavigator()->GetLevel();}
   Int_t                  GetMaxLevel() const      {return fNLevel;}
   const char            *GetPath() const;
   Int_t                  GetStackLevel() const    {return GetCurrentNavigator()->GetStackLevel();}
   TGeoVolume            *GetMasterVolume() const  {return fMasterVolume;}
   TGeoVolume            *GetTopVolume() const     {return fTopVolume;}
   TGeoNode              *GetTopNode() const       {return fTopNode;}
   TGeoPhysicalNode      *GetPhysicalNode(Int_t i) const {return (TGeoPhysicalNode*)fPhysicalNodes->UncheckedAt(i);}
   void                   SetCurrentPoint(Double_t *point) {GetCurrentNavigator()->SetCurrentPoint(point);}
   void                   SetCurrentPoint(Double_t x, Double_t y, Double_t z) {GetCurrentNavigator()->SetCurrentPoint(x,y,z);}
   void                   SetLastPoint(Double_t x, Double_t y, Double_t z) {GetCurrentNavigator()->SetLastPoint(x,y,z);}
   void                   SetCurrentDirection(Double_t *dir) {GetCurrentNavigator()->SetCurrentDirection(dir);}
   void                   SetCurrentDirection(Double_t nx, Double_t ny, Double_t nz) {GetCurrentNavigator()->SetCurrentDirection(nx,ny,nz);}
   void                   SetCldirChecked(Double_t *dir) {GetCurrentNavigator()->SetCldirChecked(dir);}

   //--- point/vector reference frame conversion
   void                   LocalToMaster(const Double_t *local, Double_t *master) const {GetCurrentNavigator()->LocalToMaster(local, master);}
   void                   LocalToMasterVect(const Double_t *local, Double_t *master) const {GetCurrentNavigator()->LocalToMasterVect(local, master);}
   void                   LocalToMasterBomb(const Double_t *local, Double_t *master) const {GetCurrentNavigator()->LocalToMasterBomb(local, master);}
   void                   MasterToLocal(const Double_t *master, Double_t *local) const {GetCurrentNavigator()->MasterToLocal(master, local);}
   void                   MasterToLocalVect(const Double_t *master, Double_t *local) const {GetCurrentNavigator()->MasterToLocalVect(master, local);}
   void                   MasterToLocalBomb(const Double_t *master, Double_t *local) const {GetCurrentNavigator()->MasterToLocalBomb(master, local);}
   void                   MasterToTop(const Double_t *master, Double_t *top) const;
   void                   TopToMaster(const Double_t *top, Double_t *master) const;

   //--- general use getters/setters
   TGeoMaterial          *FindDuplicateMaterial(const TGeoMaterial *mat) const;
   TGeoVolume            *FindVolumeFast(const char*name, Bool_t multi=kFALSE);
   TGeoMaterial          *GetMaterial(const char *matname) const;
   TGeoMaterial          *GetMaterial(Int_t id) const;
   TGeoMedium            *GetMedium(const char *medium) const;
   TGeoMedium            *GetMedium(Int_t numed) const;
   Int_t                  GetMaterialIndex(const char *matname) const;
//   TGeoShape             *GetShape(const char *name) const;
   TGeoVolume            *GetVolume(const char*name) const;
   TGeoVolume            *GetVolume(Int_t uid) const {return (TGeoVolume*)fUniqueVolumes->At(uid);}
   Int_t                  GetUID(const char *volname) const;
   Int_t                  GetNNodes() {if (!fNNodes) CountNodes(); return fNNodes;}
   TGeoNodeCache         *GetCache() const         {return GetCurrentNavigator()->GetCache();}
//   void                   SetCache(const TGeoNodeCache *cache) {fCache = (TGeoNodeCache*)cache;}
   void                   SetAnimateTracks(Bool_t flag=kTRUE) {fIsGeomReading=flag;}
   virtual ULong_t        SizeOf(const TGeoNode *node, Option_t *option); // size of the geometry in memory
   void                   SelectTrackingMedia();

   //--- stack manipulation
   Int_t                  PushPath(Int_t startlevel=0) {return GetCurrentNavigator()->PushPath(startlevel);}
   Bool_t                 PopPath() {return GetCurrentNavigator()->PopPath();}
   Bool_t                 PopPath(Int_t index) {return GetCurrentNavigator()->PopPath(index);}
   Int_t                  PushPoint(Int_t startlevel=0) {return GetCurrentNavigator()->PushPoint(startlevel);}
   Bool_t                 PopPoint() {return GetCurrentNavigator()->PopPoint();}
   Bool_t                 PopPoint(Int_t index) {return GetCurrentNavigator()->PopPoint(index);}
   void                   PopDummy(Int_t ipop=9999) {return GetCurrentNavigator()->PopDummy(ipop);}

   //--- parallel world navigation
   TGeoParallelWorld    *CreateParallelWorld(const char *name);
   TGeoParallelWorld    *GetParallelWorld() const {return fParallelWorld;}
   void                  SetUseParallelWorldNav(Bool_t flag);
   Bool_t                IsParallelWorldNav() const {return fUsePWNav;}

   ClassDef(TGeoManager, 14)          // geometry manager
};

R__EXTERN TGeoManager *gGeoManager;

#endif

 TGeoManager.h:1
 TGeoManager.h:2
 TGeoManager.h:3
 TGeoManager.h:4
 TGeoManager.h:5
 TGeoManager.h:6
 TGeoManager.h:7
 TGeoManager.h:8
 TGeoManager.h:9
 TGeoManager.h:10
 TGeoManager.h:11
 TGeoManager.h:12
 TGeoManager.h:13
 TGeoManager.h:14
 TGeoManager.h:15
 TGeoManager.h:16
 TGeoManager.h:17
 TGeoManager.h:18
 TGeoManager.h:19
 TGeoManager.h:20
 TGeoManager.h:21
 TGeoManager.h:22
 TGeoManager.h:23
 TGeoManager.h:24
 TGeoManager.h:25
 TGeoManager.h:26
 TGeoManager.h:27
 TGeoManager.h:28
 TGeoManager.h:29
 TGeoManager.h:30
 TGeoManager.h:31
 TGeoManager.h:32
 TGeoManager.h:33
 TGeoManager.h:34
 TGeoManager.h:35
 TGeoManager.h:36
 TGeoManager.h:37
 TGeoManager.h:38
 TGeoManager.h:39
 TGeoManager.h:40
 TGeoManager.h:41
 TGeoManager.h:42
 TGeoManager.h:43
 TGeoManager.h:44
 TGeoManager.h:45
 TGeoManager.h:46
 TGeoManager.h:47
 TGeoManager.h:48
 TGeoManager.h:49
 TGeoManager.h:50
 TGeoManager.h:51
 TGeoManager.h:52
 TGeoManager.h:53
 TGeoManager.h:54
 TGeoManager.h:55
 TGeoManager.h:56
 TGeoManager.h:57
 TGeoManager.h:58
 TGeoManager.h:59
 TGeoManager.h:60
 TGeoManager.h:61
 TGeoManager.h:62
 TGeoManager.h:63
 TGeoManager.h:64
 TGeoManager.h:65
 TGeoManager.h:66
 TGeoManager.h:67
 TGeoManager.h:68
 TGeoManager.h:69
 TGeoManager.h:70
 TGeoManager.h:71
 TGeoManager.h:72
 TGeoManager.h:73
 TGeoManager.h:74
 TGeoManager.h:75
 TGeoManager.h:76
 TGeoManager.h:77
 TGeoManager.h:78
 TGeoManager.h:79
 TGeoManager.h:80
 TGeoManager.h:81
 TGeoManager.h:82
 TGeoManager.h:83
 TGeoManager.h:84
 TGeoManager.h:85
 TGeoManager.h:86
 TGeoManager.h:87
 TGeoManager.h:88
 TGeoManager.h:89
 TGeoManager.h:90
 TGeoManager.h:91
 TGeoManager.h:92
 TGeoManager.h:93
 TGeoManager.h:94
 TGeoManager.h:95
 TGeoManager.h:96
 TGeoManager.h:97
 TGeoManager.h:98
 TGeoManager.h:99
 TGeoManager.h:100
 TGeoManager.h:101
 TGeoManager.h:102
 TGeoManager.h:103
 TGeoManager.h:104
 TGeoManager.h:105
 TGeoManager.h:106
 TGeoManager.h:107
 TGeoManager.h:108
 TGeoManager.h:109
 TGeoManager.h:110
 TGeoManager.h:111
 TGeoManager.h:112
 TGeoManager.h:113
 TGeoManager.h:114
 TGeoManager.h:115
 TGeoManager.h:116
 TGeoManager.h:117
 TGeoManager.h:118
 TGeoManager.h:119
 TGeoManager.h:120
 TGeoManager.h:121
 TGeoManager.h:122
 TGeoManager.h:123
 TGeoManager.h:124
 TGeoManager.h:125
 TGeoManager.h:126
 TGeoManager.h:127
 TGeoManager.h:128
 TGeoManager.h:129
 TGeoManager.h:130
 TGeoManager.h:131
 TGeoManager.h:132
 TGeoManager.h:133
 TGeoManager.h:134
 TGeoManager.h:135
 TGeoManager.h:136
 TGeoManager.h:137
 TGeoManager.h:138
 TGeoManager.h:139
 TGeoManager.h:140
 TGeoManager.h:141
 TGeoManager.h:142
 TGeoManager.h:143
 TGeoManager.h:144
 TGeoManager.h:145
 TGeoManager.h:146
 TGeoManager.h:147
 TGeoManager.h:148
 TGeoManager.h:149
 TGeoManager.h:150
 TGeoManager.h:151
 TGeoManager.h:152
 TGeoManager.h:153
 TGeoManager.h:154
 TGeoManager.h:155
 TGeoManager.h:156
 TGeoManager.h:157
 TGeoManager.h:158
 TGeoManager.h:159
 TGeoManager.h:160
 TGeoManager.h:161
 TGeoManager.h:162
 TGeoManager.h:163
 TGeoManager.h:164
 TGeoManager.h:165
 TGeoManager.h:166
 TGeoManager.h:167
 TGeoManager.h:168
 TGeoManager.h:169
 TGeoManager.h:170
 TGeoManager.h:171
 TGeoManager.h:172
 TGeoManager.h:173
 TGeoManager.h:174
 TGeoManager.h:175
 TGeoManager.h:176
 TGeoManager.h:177
 TGeoManager.h:178
 TGeoManager.h:179
 TGeoManager.h:180
 TGeoManager.h:181
 TGeoManager.h:182
 TGeoManager.h:183
 TGeoManager.h:184
 TGeoManager.h:185
 TGeoManager.h:186
 TGeoManager.h:187
 TGeoManager.h:188
 TGeoManager.h:189
 TGeoManager.h:190
 TGeoManager.h:191
 TGeoManager.h:192
 TGeoManager.h:193
 TGeoManager.h:194
 TGeoManager.h:195
 TGeoManager.h:196
 TGeoManager.h:197
 TGeoManager.h:198
 TGeoManager.h:199
 TGeoManager.h:200
 TGeoManager.h:201
 TGeoManager.h:202
 TGeoManager.h:203
 TGeoManager.h:204
 TGeoManager.h:205
 TGeoManager.h:206
 TGeoManager.h:207
 TGeoManager.h:208
 TGeoManager.h:209
 TGeoManager.h:210
 TGeoManager.h:211
 TGeoManager.h:212
 TGeoManager.h:213
 TGeoManager.h:214
 TGeoManager.h:215
 TGeoManager.h:216
 TGeoManager.h:217
 TGeoManager.h:218
 TGeoManager.h:219
 TGeoManager.h:220
 TGeoManager.h:221
 TGeoManager.h:222
 TGeoManager.h:223
 TGeoManager.h:224
 TGeoManager.h:225
 TGeoManager.h:226
 TGeoManager.h:227
 TGeoManager.h:228
 TGeoManager.h:229
 TGeoManager.h:230
 TGeoManager.h:231
 TGeoManager.h:232
 TGeoManager.h:233
 TGeoManager.h:234
 TGeoManager.h:235
 TGeoManager.h:236
 TGeoManager.h:237
 TGeoManager.h:238
 TGeoManager.h:239
 TGeoManager.h:240
 TGeoManager.h:241
 TGeoManager.h:242
 TGeoManager.h:243
 TGeoManager.h:244
 TGeoManager.h:245
 TGeoManager.h:246
 TGeoManager.h:247
 TGeoManager.h:248
 TGeoManager.h:249
 TGeoManager.h:250
 TGeoManager.h:251
 TGeoManager.h:252
 TGeoManager.h:253
 TGeoManager.h:254
 TGeoManager.h:255
 TGeoManager.h:256
 TGeoManager.h:257
 TGeoManager.h:258
 TGeoManager.h:259
 TGeoManager.h:260
 TGeoManager.h:261
 TGeoManager.h:262
 TGeoManager.h:263
 TGeoManager.h:264
 TGeoManager.h:265
 TGeoManager.h:266
 TGeoManager.h:267
 TGeoManager.h:268
 TGeoManager.h:269
 TGeoManager.h:270
 TGeoManager.h:271
 TGeoManager.h:272
 TGeoManager.h:273
 TGeoManager.h:274
 TGeoManager.h:275
 TGeoManager.h:276
 TGeoManager.h:277
 TGeoManager.h:278
 TGeoManager.h:279
 TGeoManager.h:280
 TGeoManager.h:281
 TGeoManager.h:282
 TGeoManager.h:283
 TGeoManager.h:284
 TGeoManager.h:285
 TGeoManager.h:286
 TGeoManager.h:287
 TGeoManager.h:288
 TGeoManager.h:289
 TGeoManager.h:290
 TGeoManager.h:291
 TGeoManager.h:292
 TGeoManager.h:293
 TGeoManager.h:294
 TGeoManager.h:295
 TGeoManager.h:296
 TGeoManager.h:297
 TGeoManager.h:298
 TGeoManager.h:299
 TGeoManager.h:300
 TGeoManager.h:301
 TGeoManager.h:302
 TGeoManager.h:303
 TGeoManager.h:304
 TGeoManager.h:305
 TGeoManager.h:306
 TGeoManager.h:307
 TGeoManager.h:308
 TGeoManager.h:309
 TGeoManager.h:310
 TGeoManager.h:311
 TGeoManager.h:312
 TGeoManager.h:313
 TGeoManager.h:314
 TGeoManager.h:315
 TGeoManager.h:316
 TGeoManager.h:317
 TGeoManager.h:318
 TGeoManager.h:319
 TGeoManager.h:320
 TGeoManager.h:321
 TGeoManager.h:322
 TGeoManager.h:323
 TGeoManager.h:324
 TGeoManager.h:325
 TGeoManager.h:326
 TGeoManager.h:327
 TGeoManager.h:328
 TGeoManager.h:329
 TGeoManager.h:330
 TGeoManager.h:331
 TGeoManager.h:332
 TGeoManager.h:333
 TGeoManager.h:334
 TGeoManager.h:335
 TGeoManager.h:336
 TGeoManager.h:337
 TGeoManager.h:338
 TGeoManager.h:339
 TGeoManager.h:340
 TGeoManager.h:341
 TGeoManager.h:342
 TGeoManager.h:343
 TGeoManager.h:344
 TGeoManager.h:345
 TGeoManager.h:346
 TGeoManager.h:347
 TGeoManager.h:348
 TGeoManager.h:349
 TGeoManager.h:350
 TGeoManager.h:351
 TGeoManager.h:352
 TGeoManager.h:353
 TGeoManager.h:354
 TGeoManager.h:355
 TGeoManager.h:356
 TGeoManager.h:357
 TGeoManager.h:358
 TGeoManager.h:359
 TGeoManager.h:360
 TGeoManager.h:361
 TGeoManager.h:362
 TGeoManager.h:363
 TGeoManager.h:364
 TGeoManager.h:365
 TGeoManager.h:366
 TGeoManager.h:367
 TGeoManager.h:368
 TGeoManager.h:369
 TGeoManager.h:370
 TGeoManager.h:371
 TGeoManager.h:372
 TGeoManager.h:373
 TGeoManager.h:374
 TGeoManager.h:375
 TGeoManager.h:376
 TGeoManager.h:377
 TGeoManager.h:378
 TGeoManager.h:379
 TGeoManager.h:380
 TGeoManager.h:381
 TGeoManager.h:382
 TGeoManager.h:383
 TGeoManager.h:384
 TGeoManager.h:385
 TGeoManager.h:386
 TGeoManager.h:387
 TGeoManager.h:388
 TGeoManager.h:389
 TGeoManager.h:390
 TGeoManager.h:391
 TGeoManager.h:392
 TGeoManager.h:393
 TGeoManager.h:394
 TGeoManager.h:395
 TGeoManager.h:396
 TGeoManager.h:397
 TGeoManager.h:398
 TGeoManager.h:399
 TGeoManager.h:400
 TGeoManager.h:401
 TGeoManager.h:402
 TGeoManager.h:403
 TGeoManager.h:404
 TGeoManager.h:405
 TGeoManager.h:406
 TGeoManager.h:407
 TGeoManager.h:408
 TGeoManager.h:409
 TGeoManager.h:410
 TGeoManager.h:411
 TGeoManager.h:412
 TGeoManager.h:413
 TGeoManager.h:414
 TGeoManager.h:415
 TGeoManager.h:416
 TGeoManager.h:417
 TGeoManager.h:418
 TGeoManager.h:419
 TGeoManager.h:420
 TGeoManager.h:421
 TGeoManager.h:422
 TGeoManager.h:423
 TGeoManager.h:424
 TGeoManager.h:425
 TGeoManager.h:426
 TGeoManager.h:427
 TGeoManager.h:428
 TGeoManager.h:429
 TGeoManager.h:430
 TGeoManager.h:431
 TGeoManager.h:432
 TGeoManager.h:433
 TGeoManager.h:434
 TGeoManager.h:435
 TGeoManager.h:436
 TGeoManager.h:437
 TGeoManager.h:438
 TGeoManager.h:439
 TGeoManager.h:440
 TGeoManager.h:441
 TGeoManager.h:442
 TGeoManager.h:443
 TGeoManager.h:444
 TGeoManager.h:445
 TGeoManager.h:446
 TGeoManager.h:447
 TGeoManager.h:448
 TGeoManager.h:449
 TGeoManager.h:450
 TGeoManager.h:451
 TGeoManager.h:452
 TGeoManager.h:453
 TGeoManager.h:454
 TGeoManager.h:455
 TGeoManager.h:456
 TGeoManager.h:457
 TGeoManager.h:458
 TGeoManager.h:459
 TGeoManager.h:460
 TGeoManager.h:461
 TGeoManager.h:462
 TGeoManager.h:463
 TGeoManager.h:464
 TGeoManager.h:465
 TGeoManager.h:466
 TGeoManager.h:467
 TGeoManager.h:468
 TGeoManager.h:469
 TGeoManager.h:470
 TGeoManager.h:471
 TGeoManager.h:472
 TGeoManager.h:473
 TGeoManager.h:474
 TGeoManager.h:475
 TGeoManager.h:476
 TGeoManager.h:477
 TGeoManager.h:478
 TGeoManager.h:479
 TGeoManager.h:480
 TGeoManager.h:481
 TGeoManager.h:482
 TGeoManager.h:483
 TGeoManager.h:484
 TGeoManager.h:485
 TGeoManager.h:486
 TGeoManager.h:487
 TGeoManager.h:488
 TGeoManager.h:489
 TGeoManager.h:490
 TGeoManager.h:491
 TGeoManager.h:492
 TGeoManager.h:493
 TGeoManager.h:494
 TGeoManager.h:495
 TGeoManager.h:496
 TGeoManager.h:497
 TGeoManager.h:498
 TGeoManager.h:499
 TGeoManager.h:500
 TGeoManager.h:501
 TGeoManager.h:502
 TGeoManager.h:503
 TGeoManager.h:504
 TGeoManager.h:505
 TGeoManager.h:506
 TGeoManager.h:507
 TGeoManager.h:508
 TGeoManager.h:509
 TGeoManager.h:510
 TGeoManager.h:511
 TGeoManager.h:512
 TGeoManager.h:513
 TGeoManager.h:514
 TGeoManager.h:515
 TGeoManager.h:516
 TGeoManager.h:517
 TGeoManager.h:518
 TGeoManager.h:519
 TGeoManager.h:520
 TGeoManager.h:521
 TGeoManager.h:522
 TGeoManager.h:523
 TGeoManager.h:524
 TGeoManager.h:525
 TGeoManager.h:526
 TGeoManager.h:527
 TGeoManager.h:528
 TGeoManager.h:529
 TGeoManager.h:530
 TGeoManager.h:531
 TGeoManager.h:532
 TGeoManager.h:533
 TGeoManager.h:534
 TGeoManager.h:535
 TGeoManager.h:536
 TGeoManager.h:537
 TGeoManager.h:538
 TGeoManager.h:539
 TGeoManager.h:540
 TGeoManager.h:541
 TGeoManager.h:542
 TGeoManager.h:543
 TGeoManager.h:544
 TGeoManager.h:545
 TGeoManager.h:546
 TGeoManager.h:547
 TGeoManager.h:548
 TGeoManager.h:549
 TGeoManager.h:550
 TGeoManager.h:551
 TGeoManager.h:552
 TGeoManager.h:553
 TGeoManager.h:554
 TGeoManager.h:555
 TGeoManager.h:556
 TGeoManager.h:557
 TGeoManager.h:558
 TGeoManager.h:559