]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
- TrackReference related methods and data members moved from AliDetector to AliModule
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Apr 2003 08:10:50 +0000 (08:10 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 8 Apr 2003 08:10:50 +0000 (08:10 +0000)
- AddTrackReference method added to AliModule
- AliTrackReference constructor without passing pointer to VMC.

STEER/AliDetector.cxx
STEER/AliDetector.h
STEER/AliModule.cxx
STEER/AliModule.h
STEER/AliTrackReference.cxx
STEER/AliTrackReference.h

index b6b54d54c20d61e675492612d4db71a264de9601..b0fad1a3de67042dceef353320ffb9149e686628 100644 (file)
@@ -47,7 +47,6 @@
 #include "AliHit.h"
 #include "AliPoints.h"
 #include "AliRun.h"
-#include "AliTrackReference.h"
 
 
 // Static variables for the hit iterator routines
@@ -67,10 +66,7 @@ AliDetector::AliDetector():
   fHits(0),
   fDigits(0),
   fDigitsFile(0),
-  fPoints(0),
-  fTrackReferences(0),
-  fMaxIterTrackRef(0),
-  fCurrentIterTrackRef(0)
+  fPoints(0)
 {
   //
   // Default constructor for the AliDetector class
@@ -88,10 +84,7 @@ AliDetector::AliDetector(const AliDetector &det):
   fHits(0),
   fDigits(0),
   fDigitsFile(0),
-  fPoints(0),
-  fTrackReferences(0),
-  fMaxIterTrackRef(0),
-  fCurrentIterTrackRef(0)
+  fPoints(0)
 {
   det.Copy(*this);
 }
@@ -107,10 +100,7 @@ AliDetector::AliDetector(const char* name,const char *title):
   fHits(0),
   fDigits(0),
   fDigitsFile(0),
-  fPoints(0),
-  fTrackReferences(new TClonesArray("AliTrackReference", 100)),
-  fMaxIterTrackRef(0),
-  fCurrentIterTrackRef(0)
+  fPoints(0)
 {
   //
   // Normal constructor invoked by all Detectors.
@@ -293,24 +283,6 @@ void AliDetector::FinishRun()
   //
 }
 
-//_______________________________________________________________________
-void AliDetector::RemapTrackReferencesIDs(Int_t *map)
-{
-  // 
-  // Remapping track reference
-  // Called at finish primary
-  //
-  if (!fTrackReferences) return;
-  for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
-    AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
-    if (ref) {
-      Int_t newID = map[ref->GetTrack()];
-      if (newID>=0) ref->SetTrack(newID);
-      else ref->SetTrack(-1);
-      
-    }
-  }
-}
 
 //_______________________________________________________________________
 AliHit* AliDetector::FirstHit(Int_t track)
@@ -334,27 +306,6 @@ AliHit* AliDetector::FirstHit(Int_t track)
 }
 
 
-//_______________________________________________________________________
-AliTrackReference* AliDetector::FirstTrackReference(Int_t track)
-{
-  //
-  // Initialise the hit iterator
-  // Return the address of the first hit for track
-  // If track>=0 the track is read from disk
-  // while if track<0 the first hit of the current
-  // track is returned
-  // 
-  if(track>=0) {
-    gAlice->ResetTrackReferences();
-    gAlice->TreeTR()->GetEvent(track);
-  }
-  //
-  fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
-  fCurrentIterTrackRef = 0;
-  if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
-  else            return 0;
-}
-
 //_______________________________________________________________________
 AliHit* AliDetector::NextHit()
 {
@@ -372,22 +323,6 @@ AliHit* AliDetector::NextHit()
   }
 }
 
-//_______________________________________________________________________
-AliTrackReference* AliDetector::NextTrackReference()
-{
-  //
-  // Return the next hit for the current track
-  //
-  if(fMaxIterTrackRef) {
-    if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
-      return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
-    else        
-      return 0;
-  } else {
-    printf("* AliDetector::NextTrackReference * TrackReference  Iterator called without calling FistTrackReference before\n");
-    return 0;
-  }
-}
 
 //_______________________________________________________________________
 void AliDetector::LoadPoints(Int_t)
@@ -527,16 +462,6 @@ void AliDetector::ResetHits()
   if (fHits)   fHits->Clear();
 }
 
-//_______________________________________________________________________
-void AliDetector::ResetTrackReferences()
-{
-  //
-  // Reset number of hits and the hits array
-  //
-  fMaxIterTrackRef   = 0;
-  if (fTrackReferences)   fTrackReferences->Clear();
-}
-
 //_______________________________________________________________________
 void AliDetector::ResetPoints()
 {
@@ -573,13 +498,8 @@ void AliDetector::SetTreeAddress()
     branch = treeD->GetBranch(branchname);
     if (branch) branch->SetAddress(&fDigits);
   }
-
-  // Branch address for tr  tree
-  TTree *treeTR = gAlice->TreeTR();
-  if (treeTR && fTrackReferences) {
-    branch = treeTR->GetBranch(branchname);
-    if (branch) branch->SetAddress(&fTrackReferences);
-  }
+  
+  AliModule::SetTreeAddress();
 }
 
  
index 8100fc11128411c916625e133793659a8614b3c5..bc827dbfd46eba8a403643ad6899deda06acdf85 100644 (file)
@@ -7,10 +7,8 @@
 
 #include <AliModule.h>
 class AliHit;
-class AliTrackReference;
 class TTree;
 class TBranch;
-
 class AliDetector : public AliModule {
 
 public:
@@ -27,9 +25,6 @@ public:
   virtual int   GetNhits()   const {return fNhits;}
   TClonesArray *Digits() const {return fDigits;}
   TClonesArray *Hits()   const {return fHits;}
-  TClonesArray *TrackReferences()   const {return fTrackReferences;}
-  virtual void        RemapTrackReferencesIDs(Int_t *map); //remaping track references MI
-
   TObjArray    *Points() const {return fPoints;}
   Int_t         GetIshunt() const {return fIshunt;}
   void          SetIshunt(Int_t ishunt) {fIshunt=ishunt;}
@@ -45,8 +40,6 @@ public:
   virtual void        MakeBranchTR(Option_t *opt=" ", const char *file=0 );
   virtual void        ResetDigits();
   virtual void        ResetHits();
-  virtual void        ResetTrackReferences();
-
   virtual void        ResetPoints();
   virtual void        SetTreeAddress();
   virtual void        SetTimeGate(Float_t gate) {fTimeGate=gate;}
@@ -55,8 +48,6 @@ public:
   virtual void        DrawModule() {}
   virtual AliHit*     FirstHit(Int_t track);
   virtual AliHit*     NextHit();
-  virtual AliTrackReference * FirstTrackReference(Int_t track);
-  virtual AliTrackReference * NextTrackReference();
   virtual void        SetBufferSize(Int_t bufsize=8000) {fBufferSize = bufsize;}  
   virtual TBranch*    MakeBranchInTree(TTree *tree, const char* cname, void* address, Int_t size=32000, const char *file=0);
   virtual TBranch*    MakeBranchInTree(TTree *tree, const char* cname, const char* name, void* address, Int_t size=32000, Int_t splitlevel=99, const char *file=0);
@@ -74,9 +65,7 @@ protected:
   TClonesArray *fDigits;      //List of digits for this detector
   char         *fDigitsFile;  //!File to store branches of digits tree for detector 
   TObjArray    *fPoints;      //!Array of points for each track (all tracks in memory)
-  TClonesArray *fTrackReferences; //list of track references - for one primary track only -MI
-  Int_t         fMaxIterTrackRef; //!for track refernce iterator routines
-  Int_t         fCurrentIterTrackRef; //!for track refernce iterator routines
-  ClassDef(AliDetector,2)  //Base class for ALICE detectors
+
+  ClassDef(AliDetector,3)  //Base class for ALICE detectors
 };
 #endif
index 11884bcd777b202cf4c97e703c3a5886e4f96a28..368375de4a5e0384614bc44b4f4b55a6e8209f8e 100644 (file)
 ///////////////////////////////////////////////////////////////////////////////
 #include <TNode.h>
 #include <TObjArray.h>
+#include <TClonesArray.h>
+#include <TTree.h>
 #include <TSystem.h>
+#include <TDirectory.h>
 
 #include "AliModule.h"
 #include "AliRun.h"
 #include "AliMagF.h"
 #include "AliConfig.h"
+#include "AliTrackReference.h"
 
 ClassImp(AliModule)
  
@@ -55,7 +59,10 @@ AliModule::AliModule():
   fHistograms(0),
   fNodes(0),
   fDebug(0),
-  fEnable(1)
+  fEnable(1),
+  fTrackReferences(0),
+  fMaxIterTrackRef(0),
+  fCurrentIterTrackRef(0)
 {
   //
   // Default constructor for the AliModule class
@@ -75,7 +82,10 @@ AliModule::AliModule(const char* name,const char *title):
   fHistograms(new TList()),
   fNodes(new TList()),
   fDebug(0),
-  fEnable(1)
+  fEnable(1),
+  fTrackReferences(new TClonesArray("AliTrackReference", 100)),
+  fMaxIterTrackRef(0),
+  fCurrentIterTrackRef(0)
 {
   //
   // Normal constructor invoked by all Modules.
@@ -121,7 +131,10 @@ AliModule::AliModule(const AliModule &mod):
   fHistograms(0),
   fNodes(0),
   fDebug(0),
-  fEnable(0)
+  fEnable(0),
+  fTrackReferences(0),
+  fMaxIterTrackRef(0),
+  fCurrentIterTrackRef(0)
 {
   //
   // Copy constructor
@@ -630,5 +643,137 @@ void AliModule::ReadEuclidMedia(const char* filnam)
  L20:
   Warning("ReadEuclidMedia","reading error or premature end of file\n");
 } 
+
+//_______________________________________________________________________
+void AliModule::RemapTrackReferencesIDs(Int_t *map)
+{
+  // 
+  // Remapping track reference
+  // Called at finish primary
+  //
+  if (!fTrackReferences) return;
+  for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
+    AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
+    if (ref) {
+      Int_t newID = map[ref->GetTrack()];
+      if (newID>=0) ref->SetTrack(newID);
+      else ref->SetTrack(-1);
+      
+    }
+  }
+}
+
+
+//_______________________________________________________________________
+AliTrackReference* AliModule::FirstTrackReference(Int_t track)
+{
+  //
+  // Initialise the hit iterator
+  // Return the address of the first hit for track
+  // If track>=0 the track is read from disk
+  // while if track<0 the first hit of the current
+  // track is returned
+  // 
+  if(track>=0) {
+    gAlice->ResetTrackReferences();
+    gAlice->TreeTR()->GetEvent(track);
+  }
+  //
+  fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
+  fCurrentIterTrackRef = 0;
+  if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
+  else            return 0;
+}
+
+//_______________________________________________________________________
+AliTrackReference* AliModule::NextTrackReference()
+{
+  //
+  // Return the next hit for the current track
+  //
+  if(fMaxIterTrackRef) {
+    if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
+      return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
+    else        
+      return 0;
+  } else {
+    printf("* AliDetector::NextTrackReference * TrackReference  Iterator called without calling FistTrackReference before\n");
+    return 0;
+  }
+}
+
+
+//_______________________________________________________________________
+void AliModule::ResetTrackReferences()
+{
+  //
+  // Reset number of hits and the hits array
+  //
+  fMaxIterTrackRef   = 0;
+  if (fTrackReferences)   fTrackReferences->Clear();
+}
  
  
+
+void AliModule::SetTreeAddress()
+{
+  //
+  // Set branch address for the Hits and Digits Trees
+  //
+    TBranch *branch;
+    char branchname[20];
+    sprintf(branchname,"%s",GetName());
+    // Branch address for track reference tree
+    TTree *treeTR = gAlice->TreeTR();
+    if (treeTR && fTrackReferences) {
+       branch = treeTR->GetBranch(branchname);
+       if (branch) branch->SetAddress(&fTrackReferences);
+    }
+}
+
+void  AliModule::AddTrackReference(Int_t label){
+  //
+  // add a trackrefernce to the list
+  if (!fTrackReferences) {
+    cerr<<"Container trackrefernce not active\n";
+    return;
+  }
+  Int_t nref = fTrackReferences->GetEntriesFast();
+  TClonesArray &lref = *fTrackReferences;
+  new(lref[nref]) AliTrackReference(label);
+}
+
+
+void AliModule::MakeBranchTR(Option_t *option, const char *file)
+{ 
+    //
+    // Makes branch in treeTR
+    //  
+    char name[10];
+    sprintf(name,"%s",GetName());
+
+    if (GetDebug()>1)
+       printf("* MakeBranch * Making Branch %s \n",name);
+    
+    TDirectory *cwd = gDirectory;
+    TBranch *branch = 0;
+    TTree* tree = gAlice->TreeTR();
+    if (tree) {
+       branch = tree->Branch(name, &fTrackReferences, 1600);
+       if (file) {
+           char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
+           sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
+           branch->SetFile(outFile);
+           TIter next( branch->GetListOfBranches());
+           while ((branch=dynamic_cast<TBranch*>(next()))) {
+               branch->SetFile(outFile);
+           } 
+           delete outFile;
+           
+           cwd->cd();
+           
+           if (GetDebug()>1)
+               printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
+       }
+    }
+}
index 16c6305f0cf3a981c66d3525ae894368af5d7a73..36de549445c3f1f928007ece72c29f96e4634dec 100644 (file)
@@ -17,6 +17,7 @@ class TClonesArray;
 class TBrowser;
 class TArrayI;
 class TFile;
+class AliTrackReference;
 
 class AliModule : public TNamed , public TAttLine, public TAttMarker,
                   public AliRndm {
@@ -90,21 +91,16 @@ public:
   virtual void        FinishEvent() {}
   virtual void        FinishRun() {}
   virtual void        FinishPrimary() {}
-  virtual void        RemapTrackHitIDs(Int_t *) {}
-  virtual void        RemapTrackReferencesIDs(Int_t *) {} //remaping track references MI
-
   //virtual void        Hits2Digits() {}
   virtual void        Init() {}
   virtual void        LoadPoints(Int_t ) {}
   virtual void        MakeBranch(Option_t *, const char* =0 ) {} 
-  virtual void        MakeBranchTR(Option_t * =" ", const char * =0 ){}
   virtual void        Paint(Option_t *) {}
   virtual void        ResetDigits() {}
   virtual void        ResetSDigits() {}
   virtual void        ResetHits() {}
-  virtual void        ResetTrackReferences() {}
   virtual void        ResetPoints() {}
-  virtual void        SetTreeAddress() {}
+  virtual void        SetTreeAddress();
   virtual void        SetTimeGate(Float_t) {}
   virtual Float_t     GetTimeGate() const {return 1.e10;}
   virtual void        StepManager() {}
@@ -116,6 +112,17 @@ public:
   virtual void        SetEuclidFile(char *material,char *geometry=0);
   virtual void ReadEuclid(const char *filnam, char *topvol);
   virtual void ReadEuclidMedia(const char *filnam);
+// Track reference related
+  TClonesArray *TrackReferences()   const {return fTrackReferences;}
+  virtual void        RemapTrackHitIDs(Int_t *) {}
+  virtual void        RemapTrackReferencesIDs(Int_t *map); //remaping track references MI
+  virtual void        ResetTrackReferences();
+  virtual void  AddTrackReference(Int_t label);
+  virtual  AliTrackReference * FirstTrackReference(Int_t track);
+  virtual  AliTrackReference * NextTrackReference();
+  virtual void MakeBranchTR(Option_t *option, const char *file);
+  
+//
   AliModule& operator=(const AliModule &mod)
     {mod.Copy(*this); return (*this);}
  
@@ -137,6 +144,9 @@ protected:
   TList        *fNodes;       //List of geometry nodes
   Int_t         fDebug;       //Debug flag
   Bool_t        fEnable;      //StepManager enabling flag
-  ClassDef(AliModule,2)  //Base class for ALICE Modules
+  TClonesArray *fTrackReferences;     //list of track references - for one primary track only -MI
+  Int_t         fMaxIterTrackRef;     //!for track refernce iterator routines
+  Int_t         fCurrentIterTrackRef; //!for track refernce iterator routines
+  ClassDef(AliModule,3)  //Base class for ALICE Modules
 };
 #endif
index ffd513bc276a71ebc15fa92c250c7f0b3ee602a9..7c3fb26a5d643a0f75a297b00438fa3cc1fe3bb6 100644 (file)
@@ -15,6 +15,7 @@
 
 
 #include "AliTrackReference.h"
+#include "TVirtualMC.h"
 #include "TParticle.h"
 #include "AliRun.h"
 #include "TLorentzVector.h"
@@ -53,7 +54,7 @@ ClassImp(AliTrackReference)
 }
 
 //_______________________________________________________________________
-AliTrackReference::AliTrackReference(Int_t label, TVirtualMC *vMC) {
+AliTrackReference::AliTrackReference(Int_t label) {
   //
   // Create Reference object out of label and
   // data in TVirtualMC object
@@ -64,20 +65,20 @@ AliTrackReference::AliTrackReference(Int_t label, TVirtualMC *vMC) {
   // Sylwester Radomski, (S.Radomski@gsi.de)
   // GSI, Jan 31, 2003
   //
-
+    
   TLorentzVector vec;
   
   fTrack = label;
-  fLength = vMC->TrackLength();
-  fTime = vMC->TrackTime();
+  fLength = gMC->TrackLength();
+  fTime = gMC->TrackTime();
 
-  vMC->TrackPosition(vec);
+  gMC->TrackPosition(vec);
 
   fX = vec[0];
   fY = vec[1];
   fZ = vec[2];
   
-  vMC->TrackMomentum(vec);
+  gMC->TrackMomentum(vec);
   
   fPx = vec[0];
   fPy = vec[1];
@@ -88,13 +89,13 @@ AliTrackReference::AliTrackReference(Int_t label, TVirtualMC *vMC) {
 
   for(Int_t i=0; i<16; i++) ResetBit(BIT(i));
 
-  SetBit(BIT(0), vMC->IsNewTrack());
-  SetBit(BIT(1), vMC->IsTrackAlive());
-  SetBit(BIT(2), vMC->IsTrackDisappeared());
-  SetBit(BIT(3), vMC->IsTrackEntering());
-  SetBit(BIT(4), vMC->IsTrackExiting());
-  SetBit(BIT(5), vMC->IsTrackInside());
-  SetBit(BIT(6), vMC->IsTrackOut());
-  SetBit(BIT(7), vMC->IsTrackStop()); 
+  SetBit(BIT(0), gMC->IsNewTrack());
+  SetBit(BIT(1), gMC->IsTrackAlive());
+  SetBit(BIT(2), gMC->IsTrackDisappeared());
+  SetBit(BIT(3), gMC->IsTrackEntering());
+  SetBit(BIT(4), gMC->IsTrackExiting());
+  SetBit(BIT(5), gMC->IsTrackInside());
+  SetBit(BIT(6), gMC->IsTrackOut());
+  SetBit(BIT(7), gMC->IsTrackStop()); 
 }
 //_______________________________________________________________________
index 204be7b50f7f6ca3cbf39ef16089c1209ce75c20..04ebf7b41ed60fe9d14f198d0914189f9837431f 100644 (file)
@@ -6,14 +6,14 @@
 /* $Id$ */
 
 #include "TObject.h"
-#include "TVirtualMC.h"
+#include "TMath.h"
 
 class AliTrackReference : public TObject {
 
 public:
 
   AliTrackReference();
-  AliTrackReference(Int_t label, TVirtualMC *vMC);
+  AliTrackReference(Int_t label);
   virtual ~AliTrackReference() {}
 
   virtual Int_t GetTrack() const {return fTrack;}