]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Add special treatment for TPC and TRD, they use different trees than other detectors
authorjchudoba <jchudoba@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 15 Nov 2001 09:00:11 +0000 (09:00 +0000)
committerjchudoba <jchudoba@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 15 Nov 2001 09:00:11 +0000 (09:00 +0000)
STEER/AliRunDigitizer.cxx
STEER/AliRunDigitizer.h

index 9f2edac74f25d1bd722a18d503a8c5cd12ed5484..44b9658e10a46a0f9ecbbbf0924031513f974b39 100644 (file)
@@ -15,6 +15,9 @@
 
 /*
 $Log$
+Revision 1.8  2001/10/21 18:38:43  hristov
+Several pointers were set to zero in the default constructors to avoid memory management problems
+
 Revision 1.7  2001/10/04 15:56:07  jchudoba
 TTask inheritance
 
@@ -238,6 +241,9 @@ Bool_t AliRunDigitizer::ConnectInputTrees()
       sprintf(treeName,"TreeS_75x40_100x60_%d",serialNr);
       tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
       fArrayTreeTPCS[i] = tree;
+      sprintf(treeName,"TreeS%d_TRD",serialNr);
+      tree = static_cast<TTree*>(iStream->CurrentFile()->Get(treeName));
+      fArrayTreeTRDS[i] = tree;
     } else if (delta[i] != 0) {
       Error("ConnectInputTrees","Only delta 0 or 1 is implemented");
       return kFALSE;
@@ -306,6 +312,23 @@ void AliRunDigitizer::InitEvent()
     fTreeD = new TTree(treeName,"Digits");
     fTreeD->Write(0,TObject::kOverwrite);
   }
+
+// special tree for TPC
+  sprintf(treeName,"TreeD_75x40_100x60_%d",fEvent);
+  fTreeDTPC = static_cast<TTree*>(fOutput->Get(treeName));
+  if (!fTreeDTPC) {
+    fTreeDTPC = new TTree(treeName,"TPC_Digits");
+    fTreeDTPC->Write(0,TObject::kOverwrite);
+  }
+
+// special tree for TRD
+  sprintf(treeName,"TreeD%d_TRD",fEvent);
+  fTreeDTRD = static_cast<TTree*>(fOutput->Get(treeName));
+  if (!fTreeDTRD) {
+    fTreeDTRD = new TTree(treeName,"TRD_Digits");
+    fTreeDTRD->Write(0,TObject::kOverwrite);
+  }
+
 }
 
 ////////////////////////////////////////////////////////////////////////
@@ -328,6 +351,14 @@ void AliRunDigitizer::FinishEvent()
     delete fTreeD;
     fTreeD = 0;
   }
+  if (fTreeDTPC) {
+    delete fTreeDTPC;
+    fTreeDTPC = 0;
+  }
+  if (fTreeDTRD) {
+    delete fTreeDTRD;
+    fTreeDTRD = 0;
+  }
 }
 ////////////////////////////////////////////////////////////////////////
 void AliRunDigitizer::FinishGlobal()
index 63ec6269e987a843be0f39711828a74495fe76b7..e5a99ee3226919f7d858a6067c9383374dd59798 100644 (file)
@@ -53,7 +53,10 @@ public:
   TTree*    GetInputTreeS(Int_t i) const {return fArrayTreeS[i];}
   TTree*    GetInputTreeH(Int_t i) const {return fArrayTreeH[i];}
   TTree*    GetInputTreeTPCS(Int_t i) const {return fArrayTreeTPCS[i];}
+  TTree*    GetInputTreeTRDS(Int_t i) const {return fArrayTreeTRDS[i];}
   TTree*    GetTreeD() const {return fTreeD;}
+  TTree*    GetTreeDTPC() const {return fTreeDTPC;} 
+  TTree*    GetTreeDTRD() const {return fTreeDTRD;} 
   void      Digitize(Option_t* option = 0);
   void      Exec(Option_t *option) {this->Digitize();}
   void      ExecuteTask(Option_t* option = 0);
@@ -100,12 +103,15 @@ private:
   Int_t             fCopyTreesFromInput;  // from which input file the trees
                                           // should be copied, -1 for no copies
   TTree *           fTreeD;               //! output TreeD
+  TTree *           fTreeDTPC;            //! output TreeD for TPC
+  TTree *           fTreeDTRD;            //! output TreeD for TRD
   Int_t             fNinputs;             // nr of input streams - can be taken from the TClonesArray dimension
   Int_t             fNinputsGiven;        // nr of input streams given by user
   TClonesArray *    fInputStreams;        // input streams
-  TFile *          fInputFiles[MAXSTREAMSTOMERGE];   //! p. to current input files
+  TFile *           fInputFiles[MAXSTREAMSTOMERGE];   //! p. to current input files
   TTree *           fArrayTreeS[MAXSTREAMSTOMERGE];   //! array with p. to TreeS
   TTree *           fArrayTreeTPCS[MAXSTREAMSTOMERGE];   //! array with p. to TreeD_75x40_100x60_x (TPC Sdigits)
+  TTree *           fArrayTreeTRDS[MAXSTREAMSTOMERGE];   //! array with p. to TreeSx_TRD (TRD Sdigits)
   TTree *           fArrayTreeH[MAXSTREAMSTOMERGE];   //! array with p. to TreeH
   AliMergeCombi *   fCombi;               // pointer to the combination object
   TArrayI           fCombination;         //! combination of events from
@@ -119,7 +125,7 @@ private:
   void              FinishGlobal();
   Int_t             fDebug;                //! specifies debug level, 0 is min
   
-  ClassDef(AliRunDigitizer,2)
+  ClassDef(AliRunDigitizer,3)
 };
 
 #endif // ALIRUNDIGITIZER_H