- remove legacy code for old raw stream
authorcblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 19 Feb 2011 17:35:11 +0000 (17:35 +0000)
committercblume <cblume@f7af4fe6-9843-0410-8265-dc069ae4e863>
Sat, 19 Feb 2011 17:35:11 +0000 (17:35 +0000)
- prepare tracklet/track writing in case of digits conversion

TRD/AliTRDrawData.cxx
TRD/AliTRDrawData.h

index e31406a..794dc12 100644 (file)
@@ -57,7 +57,8 @@ AliTRDrawData::AliTRDrawData()
   ,fFee(NULL)
   ,fNumberOfDDLs(0)
   ,fTrackletTree(NULL)
-  ,fTrackletContainer(NULL)
+  ,fTracklets(NULL)
+  ,fTracks(NULL)
   ,fSMindexPos(0)
   ,fStackindexPos(0)
   ,fEventCounter(0)
@@ -81,7 +82,8 @@ AliTRDrawData::AliTRDrawData(const AliTRDrawData &r)
   ,fFee(NULL)
   ,fNumberOfDDLs(0)
   ,fTrackletTree(NULL)
-  ,fTrackletContainer(NULL)
+  ,fTracklets(NULL)
+  ,fTracks(NULL)
   ,fSMindexPos(0)
   ,fStackindexPos(0)
   ,fEventCounter(0)
@@ -104,9 +106,14 @@ AliTRDrawData::~AliTRDrawData()
   // Destructor
   //
 
-  if (fTrackletContainer){
-    delete fTrackletContainer;
-    fTrackletContainer = NULL;
+  if (fTracklets){
+    fTracklets->Delete();
+    delete fTracklets;
+  }
+
+  if (fTracks){
+    fTracks->Delete();
+    delete fTracks;
   }
 
   delete fMcmSim;
@@ -542,24 +549,11 @@ AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
   AliTRDdigitsManager* digitsManager = new AliTRDdigitsManager();
   digitsManager->CreateArrays();
 
-  if (!fTrackletContainer) {
-    const Int_t kTrackletChmb=256;
-    fTrackletContainer = new UInt_t *[2];
-    fTrackletContainer[0] = new UInt_t[kTrackletChmb];
-    fTrackletContainer[1] = new UInt_t[kTrackletChmb];
-    memset(fTrackletContainer[0], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
-    memset(fTrackletContainer[1], 0, kTrackletChmb*sizeof(UInt_t)); //jkl
-  }
-
-  AliTRDrawStream *pinput = new AliTRDrawStream(rawReader);
-  AliTRDrawStream &input  = *pinput;
-
-  AliInfo(Form("Stream version: %s", input.IsA()->GetName()));
+  AliTRDrawStream input(rawReader);
 
   // ----- preparing tracklet output -----
   AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets");
   if (!trklLoader) {
-    //AliInfo("Could not get the tracklets data loader, adding it now!");
     trklLoader = new AliDataLoader("TRD.Tracklets.root","tracklets", "tracklets");
     AliRunLoader::Instance()->GetLoader("TRDLoader")->AddDataLoader(trklLoader);
   }
@@ -578,14 +572,15 @@ AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
   if (!trklTreeLoader->Tree())
     trklTreeLoader->MakeTree();
 
+  input.SetTrackletArray(TrackletsArray());
+  input.SetTrackArray(TracksArray());
+
   // Loop through the digits
   Int_t det    = 0;
 
   while (det >= 0)
     {
-      det = input.NextChamber(digitsManager,fTrackletContainer);
-
-      if (*(fTrackletContainer[0]) > 0 || *(fTrackletContainer[1]) > 0) WriteTracklets(det);
+      det = input.NextChamber(digitsManager);
 
       if (det >= 0)
        {
@@ -605,18 +600,7 @@ AliTRDdigitsManager *AliTRDrawData::Raw2Digits(AliRawReader *rawReader)
   trklTreeLoader->WriteData("OVERWRITE");
   trklLoader->UnloadAll();
 
-  if (fTrackletContainer){
-    delete [] fTrackletContainer[0];
-    delete [] fTrackletContainer[1];
-    delete [] fTrackletContainer;
-    fTrackletContainer = NULL;
-  }
-
-  delete pinput;
-  pinput = NULL;
-
   return digitsManager;
-
 }
 
 //_____________________________________________________________________________
@@ -678,43 +662,22 @@ void AliTRDrawData::WriteIntermediateWords(UInt_t* buf, Int_t& nw, Int_t& of, co
     if (nw < maxSize) buf[nw++] = x; else of++;
 }
 
-//_____________________________________________________________________________
-Bool_t AliTRDrawData::WriteTracklets(Int_t det)
+TClonesArray *AliTRDrawData::TrackletsArray()
 {
-  //
-  // Write the raw data tracklets into seperate file
-  //
-
-  UInt_t **leaves = new UInt_t *[2];
-  for (Int_t i=0; i<2 ;i++){
-    leaves[i] = new UInt_t[258];
-    leaves[i][0] = det; // det
-    leaves[i][1] = i;   // side
-    memcpy(leaves[i]+2, fTrackletContainer[i], sizeof(UInt_t) * 256);
-  }
+  // Returns the array of on-line tracklets
 
-  if (!fTrackletTree){
-    AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets");
-    dl->MakeTree();
-    fTrackletTree = dl->Tree();
+  if (!fTracklets) {
+    fTracklets = new TClonesArray("AliTRDtrackletWord", 2*MAXTRACKLETSPERHC);
   }
+  return fTracklets;
+}
 
-  TBranch *trkbranch = fTrackletTree->GetBranch("trkbranch");
-  if (!trkbranch) {
-    trkbranch = fTrackletTree->Branch("trkbranch",leaves[0],"det/i:side/i:tracklets[256]/i");
-  }
+TClonesArray* AliTRDrawData::TracksArray()
+{
+  // return array of GTU tracks (create TClonesArray if necessary)
 
-  for (Int_t i=0; i<2; i++){
-    if (leaves[i][2]>0) {
-      trkbranch->SetAddress(leaves[i]);
-      fTrackletTree->Fill();
-    }
+  if (!fTracks) {
+    fTracks = new TClonesArray("AliESDTrdTrack",100);
   }
-
-  //  AliDataLoader *dl = fRunLoader->GetLoader("TRDLoader")->GetDataLoader("tracklets"); //jkl: wrong
-  //  dl->WriteData("OVERWRITE"); //jkl: wrong
-  //dl->Unload();
-  delete [] leaves;
-
-  return kTRUE;
+  return fTracks;
 }
index 564330d..11275a3 100644 (file)
@@ -14,6 +14,7 @@
 #include "TObject.h"
 
 class TTree;
+class TClonesArray;
 
 class AliRunLoader;
 
@@ -39,7 +40,11 @@ class AliTRDrawData : public TObject {
   virtual Bool_t       Digits2Raw(TTree *digits, const TTree *tracks = NULL);
 
   virtual AliTRDdigitsManager *Raw2Digits(AliRawReader *rawReader);
-  Bool_t WriteTracklets(Int_t det);
+
+  virtual TClonesArray    *TrackletsArray();
+  virtual TClonesArray    *TracksArray();
+  void                    SetTrackletsOwner(TClonesArray *trkl = 0x0) { fTracklets = trkl; } // set to own the given array
+  void                    SetTracksOwner(TClonesArray *trk = 0x0) { fTracks = trk; } // set to own the given array
 
  protected:
 
@@ -56,8 +61,10 @@ class AliTRDrawData : public TObject {
   AliTRDgeometry      *fGeo;            //! Geometry
   AliTRDfeeParam      *fFee;            //! Fee Parameters
   Int_t                fNumberOfDDLs;   //  Number of DDLs
-  TTree               *fTrackletTree;        //! Tree for tracklets
-  UInt_t              **fTrackletContainer;  //! tracklet container
+  TTree               *fTrackletTree;   //! Tree for tracklets
+
+  TClonesArray        *fTracklets;      //! Array of online tracklets
+  TClonesArray        *fTracks;         //! Array of GTU tracks
 
  private:
 
@@ -71,7 +78,7 @@ class AliTRDrawData : public TObject {
   AliTRDmcmSim      *fMcmSim;         //! MCM simulation for raw data output
   AliTRDdigitsParam *fDigitsParam;    // Digits parameter
 
-  ClassDef(AliTRDrawData,7)             //  TRD raw data class
+  ClassDef(AliTRDrawData,8)             //  TRD raw data class
 
 };
 #endif