]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
write on-line tracklets and tracks to ESD during reconstruction
authorjklein <jklein@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 6 Jul 2011 07:26:16 +0000 (07:26 +0000)
committerjklein <jklein@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 6 Jul 2011 07:26:16 +0000 (07:26 +0000)
- implement FillESD in Reconstructor
- parse tracking and trigger headers in rawStream
- fix reading of data with number of timebins != 3n

TRD/AliTRDReconstructor.cxx
TRD/AliTRDReconstructor.h
TRD/AliTRDclusterizer.cxx
TRD/AliTRDclusterizer.h
TRD/AliTRDrawStream.cxx
TRD/AliTRDtrackletBase.h
TRD/AliTRDtrackletWord.h

index 26fe457ce3259704fcbb5f6d5481f2245bd1b131..a457fb75f81ecb9d84bbc7c83977126161ac83fe 100644 (file)
 #include <TObjArray.h>
 #include <TTreeStream.h>
 #include <TDirectory.h>
+#include <TRef.h>
 
 #include "AliRawReader.h"
+#include "AliLog.h"
 
 #include "AliTRDReconstructor.h"
 #include "AliTRDclusterizer.h"
 #include "AliTRDrawStream.h"
 #include "AliTRDdigitsManager.h"
 #include "AliTRDtrackerV1.h"
+#include "AliESDEvent.h"
+#include "AliESDTrdTrack.h"
+#include "AliESDTrdTracklet.h"
+#include "AliTRDtrackletWord.h"
 
 #define SETFLG(n,f) ((n) |= f)
 #define CLRFLG(n,f) ((n) &= ~f)
@@ -45,6 +51,7 @@ ClassImp(AliTRDReconstructor)
 
 TClonesArray *AliTRDReconstructor::fgClusters = NULL;
 TClonesArray *AliTRDReconstructor::fgTracklets = NULL;
+TClonesArray *AliTRDReconstructor::fgTracks = NULL;
 Char_t const * AliTRDReconstructor::fgSteerNames[kNsteer] = {
   "DigitsConversion       "
  ,"Write Clusters         "
@@ -118,6 +125,11 @@ AliTRDReconstructor::~AliTRDReconstructor()
     delete fgTracklets;
     fgTracklets = NULL;
   }
+  if(fgTracks) {
+    fgTracks->Delete();
+    delete fgTracks;
+    fgTracks = NULL;
+  }
   if(fSteerParam&kOwner){
     for(Int_t itask = 0; itask < AliTRDrecoParam::kTRDreconstructionTasks; itask++)
       if(fDebugStream[itask]) delete fDebugStream[itask];
@@ -171,6 +183,13 @@ void AliTRDReconstructor::ConvertDigits(AliRawReader *rawReader
   manager->WriteDigits();
   delete manager;
 
+  // take over ownership of online tracklets
+  fgTracklets = rawData.TrackletsArray();
+  rawData.SetTrackletsOwner(0x0);
+
+  // take over GTU tracks
+  fgTracks = rawData.TracksArray();
+  rawData.SetTracksOwner(0x0);
 }
 
 //_____________________________________________________________________________
@@ -200,15 +219,19 @@ void AliTRDReconstructor::Reconstruct(AliRawReader *rawReader
   
   fgNTimeBins = fClusterizer->GetNTimeBins();
   
+  // take over ownership of online tracklets
+  fgTracklets = fClusterizer->TrackletsArray();
+  fClusterizer->SetTrackletsOwner(kFALSE);
+
+  // take over GTU tracks
+  fgTracks = fClusterizer->TracksArray();
+  fClusterizer->SetTracksOwner(kFALSE);
+
   if(IsWritingClusters()) return;
 
   // take over ownership of clusters
   fgClusters = fClusterizer->RecPoints();
   fClusterizer->SetClustersOwner(kFALSE);
-
-  // take over ownership of online tracklets
-  fgTracklets = fClusterizer->TrackletsArray();
-  fClusterizer->SetTrackletsOwner(kFALSE);
 }
 
 //_____________________________________________________________________________
@@ -227,6 +250,18 @@ void AliTRDReconstructor::Reconstruct(TTree *digitsTree
   clusterer.ReadDigits(digitsTree);
   clusterer.MakeClusters();
 
+  // read tracklets and tracks if not done during reading of raw data
+  if (!fgTracklets) {
+    clusterer.ReadTracklets();
+    fgTracklets = clusterer.TrackletsArray();
+    clusterer.SetTrackletsOwner(kFALSE);
+  }
+  if (!fgTracks) {
+    clusterer.ReadTracks();
+    fgTracks = clusterer.TracksArray();
+    clusterer.SetTracksOwner(kFALSE);
+  }
+
   fgNTimeBins = clusterer.GetNTimeBins();
 
   if(IsWritingClusters()) return;
@@ -234,11 +269,6 @@ void AliTRDReconstructor::Reconstruct(TTree *digitsTree
   // take over ownership of clusters
   fgClusters = clusterer.RecPoints();
   clusterer.SetClustersOwner(kFALSE);
-
-  // take over ownership of online tracklets
-  fgTracklets = clusterer.TrackletsArray();
-  clusterer.SetTrackletsOwner(kFALSE);
-
 }
 
 //_____________________________________________________________________________
@@ -258,12 +288,130 @@ AliTracker *AliTRDReconstructor::CreateTracker() const
 //_____________________________________________________________________________
 void AliTRDReconstructor::FillESD(TTree* /*digitsTree*/
         , TTree* /*clusterTree*/
-        , AliESDEvent* /*esd*/) const
+        , AliESDEvent* esd) const
 {
   //
   // Fill ESD
   //
 
+  // ----- filling tracklets -----
+  Int_t trackletIndex[1080] = { 0 };
+  AliDebug(1, Form("Filling tracklets from %p (%i)",
+                  fgTracklets, fgTracklets ? fgTracklets->GetEntriesFast() : 0));
+  if (fgTracklets) {
+    Int_t nTracklets = fgTracklets->GetEntriesFast();
+
+    Int_t lastHC = -1;
+    for (Int_t iTracklet = 0; iTracklet < nTracklets; iTracklet++) {
+      AliTRDtrackletBase *trkl = (AliTRDtrackletBase*) ((*fgTracklets)[iTracklet]);
+      Int_t hc = trkl->GetHCId();
+      //      hc += trkl->GetY() > 0 ? 1 : 0;
+      if ((hc < 0) || (hc >= 1080)) {
+       AliError(Form("HC for tracklet: 0x%08x out of range: %i", trkl->GetTrackletWord(), trkl->GetHCId()));
+       continue;
+}
+      AliDebug(5, Form("hc: %4i : 0x%08x z: %2i", hc, trkl->GetTrackletWord(), trkl->GetZbin()));
+      if (hc != lastHC) {
+       AliDebug(2, Form("set tracklet index for HC %i to %i", hc, iTracklet));
+       trackletIndex[hc] = iTracklet + 1;
+       lastHC = hc;
+      }
+    }
+
+    for (Int_t iDet = 0; iDet < 540; iDet++) {
+      Int_t trklIndexA = trackletIndex[2*iDet + 0] - 1;
+      Int_t trklIndexB = trackletIndex[2*iDet + 1] - 1;
+      Int_t trklIndex  = esd->GetNumberOfTrdTracklets();
+      AliTRDtrackletBase *trklA = trklIndexA > -1 ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexA]) : 0x0;
+      AliTRDtrackletBase *trklB = trklIndexB > -1 ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexB]) : 0x0;
+      while (trklA != 0x0 || trklB != 0x0) {
+       AliDebug(5, Form("det %i - A: %i/%i -> %p, B: %i/%i -> %p",
+                        iDet, trklIndexA, nTracklets, trklA, trklIndexB, nTracklets, trklB));
+       if (trklA == 0x0) {
+         AliESDTrdTracklet esdTracklet(trklB->GetTrackletWord(), trklB->GetHCId());
+         esd->AddTrdTracklet(&esdTracklet);
+         trklIndexB++;
+         trklB = trklIndexB < nTracklets ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexB]) : 0x0;
+         if (trklB && trklB->GetHCId() != 2*iDet + 1)
+           trklB = 0x0;
+       }
+       else if (trklB == 0x0) {
+         AliESDTrdTracklet esdTracklet(trklA->GetTrackletWord(), trklA->GetHCId());
+         esd->AddTrdTracklet(&esdTracklet);
+         trklIndexA++;
+         trklA = trklIndexA < nTracklets ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexA]) : 0x0;
+         if (trklA && trklA->GetHCId() != 2*iDet)
+           trklA = 0x0;
+       }
+       else {
+         if ((trklA->GetZbin() < trklB->GetZbin()) ||
+             ((trklA->GetZbin() == trklB->GetZbin()) && (trklA->GetYbin() < trklB->GetYbin()))) {
+           AliESDTrdTracklet esdTracklet(trklA->GetTrackletWord(), trklA->GetHCId());
+           esd->AddTrdTracklet(&esdTracklet);
+           trklIndexA++;
+           trklA = trklIndexA < nTracklets ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexA]) : 0x0;
+           if (trklA && trklA->GetHCId() != 2*iDet)
+             trklA = 0x0;
+         }
+         else {
+           AliESDTrdTracklet esdTracklet(trklB->GetTrackletWord(), trklB->GetHCId());
+           esd->AddTrdTracklet(&esdTracklet);
+           trklIndexB++;
+           trklB = trklIndexB < nTracklets ? (AliTRDtrackletBase*) ((*fgTracklets)[trklIndexB]) : 0x0;
+           if (trklB && trklB->GetHCId() != 2*iDet + 1)
+             trklB = 0x0;
+         }
+       }
+      }
+      // updating tracklet indices as in ESD
+      if (esd->GetNumberOfTrdTracklets() != trklIndex) {
+       trackletIndex[2*iDet + 0] = trackletIndex[2*iDet + 1] = trklIndex;
+      }
+      else
+       trackletIndex[2*iDet + 0] = trackletIndex[2*iDet + 1] = -1;
+    }
+  }
+
+  // ----- filling GTU tracks -----
+  AliDebug(1, Form("Now filling ESD with GTU tracks from %p (%i)",
+                  fgTracks, fgTracks ? fgTracks->GetEntriesFast() : 0));
+  if (fgTracks) {
+    for (Int_t iTrack = 0; iTrack < fgTracks->GetEntriesFast(); iTrack++) {
+      AliESDTrdTrack *trdTrack = (AliESDTrdTrack*) ((*fgTracks)[iTrack]);
+
+      UInt_t mask  = trdTrack->GetLayerMask();
+      UInt_t stack = trdTrack->GetStack();
+
+      for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
+       if (mask & (1 << iLayer)) {
+
+         Int_t det = trdTrack->GetSector()*30 + stack*6 + iLayer;
+         Int_t idx = trdTrack->GetTrackletIndex(iLayer);
+
+         if ((det < 0) || (det > 539)) {
+           AliError(Form("Invalid detector no. from track: %i", 2*det));
+           continue;
+         }
+         if (trackletIndex[2*det] >= 0) {
+           if ((trackletIndex[2*det] + idx > -1) &&
+               (trackletIndex[2*det] + idx < esd->GetNumberOfTrdTracklets())) {
+             AliESDTrdTracklet *trkl = esd->GetTrdTracklet(trackletIndex[2*det] + idx);
+             if (trkl) {
+               AliDebug(5, Form("adding tracklet: 0x%08x", trkl->GetTrackletWord()));
+               trdTrack->AddTrackletReference(trkl, iLayer);
+             }
+           }
+         }
+       }
+      }
+      // only add the track when it's complete (including tracklet references)
+      esd->AddTrdTrack(trdTrack);
+    }
+  }
+
+  // clearing variables for next event
+  fgTracklets = 0x0;
+  fgTracks = 0x0;
 }
 
 //_____________________________________________________________________________
index ea0298dfc750c32b4e9590622cb23d375c41733f..c15ff7e03f3f5d7aff3e7d1383cde24dd11f779d 100644 (file)
@@ -50,6 +50,7 @@ public:
   virtual void        FillESD(TTree *digitsTree, TTree *clusterTree, AliESDEvent *esd) const;
   static TClonesArray* GetClusters()             {return fgClusters;}
   static TClonesArray* GetTracklets()            { return fgTracklets;}
+  static TClonesArray* GetTracks()               { return fgTracks;}
   static Int_t        GetNTimeBins()             { return fgNTimeBins;}
   Int_t               GetNdEdxSlices() const     { return (Int_t)AliTRDpidUtil::GetNdEdxSlices(GetPIDMethod());}
   AliTRDpidUtil::ETRDPIDMethod       GetPIDMethod() const       { return GetRecoParam()->IsPIDNeuralNetwork() ? AliTRDpidUtil::kNN : AliTRDpidUtil::kLQ;}
@@ -77,6 +78,7 @@ public:
 
   static void         SetClusters(TClonesArray *clusters)  { fgClusters = clusters;} 
   static void         SetTracklets(TClonesArray *tracklets) { fgTracklets = tracklets;}
+  static void         SetTracks(TClonesArray *tracks) { fgTracks = tracks;}
   void               SetOption(Option_t *opt);
 
 private:
@@ -100,10 +102,11 @@ private:
  
   static TClonesArray *fgClusters;    //  list of clusters for local reconstructor
   static TClonesArray *fgTracklets;   //  list of online tracklets for local reconstructor
+  static TClonesArray *fgTracks;      //  list of GTU tracks for local reconstructor
   static Int_t         fgNTimeBins;   //  number of time bins as given by the clusterizer
   AliTRDclusterizer   *fClusterizer;  //! instance of TRD clusterizer
 
-  ClassDef(AliTRDReconstructor, 4)    //  Class for the TRD reconstruction
+  ClassDef(AliTRDReconstructor, 5)    //  Class for the TRD reconstruction
 
 };
 
index 82637dbe7a612a4ca8885fa6b39e64cf99bd64dc..08b34a0972431bf9c3804984e6e7047efe20d97e 100644 (file)
@@ -44,6 +44,9 @@
 #include "AliTRDrawStream.h"
 #include "AliTRDfeeParam.h"
 #include "AliTRDtrackletWord.h"
+#include "AliTRDtrackletMCM.h"
+#include "AliTRDtrackGTU.h"
+#include "AliESDTrdTrack.h"
 
 #include "TTreeStream.h"
 
@@ -483,6 +486,83 @@ Bool_t AliTRDclusterizer::ReadDigits(AliRawReader *rawReader)
 
 }
 
+Bool_t AliTRDclusterizer::ReadTracklets()
+{
+  //
+  // Reads simulated tracklets from the input aliroot file
+  //
+
+  AliRunLoader *runLoader = AliRunLoader::Instance();
+  if (!runLoader) {
+    AliError("No run loader available");
+    return kFALSE;
+  }
+
+  AliLoader* loader = runLoader->GetLoader("TRDLoader");
+
+  AliDataLoader *trackletLoader = loader->GetDataLoader("tracklets");
+  if (!trackletLoader) {
+      return kFALSE;
+  }
+
+  // simulated tracklets
+  trackletLoader->Load();
+  TTree *trackletTree = trackletLoader->Tree();
+
+ if (trackletTree) {
+   TBranch *trklbranch = trackletTree->GetBranch("mcmtrklbranch");
+   TClonesArray *trklArray = TrackletsArray("AliTRDtrackletMCM");
+   if (trklbranch && trklArray) {
+     AliTRDtrackletMCM *trkl = 0x0;
+     trklbranch->SetAddress(&trkl);
+     for (Int_t iTracklet = 0; iTracklet < trklbranch->GetEntries(); iTracklet++) {
+       trklbranch->GetEntry(iTracklet);
+       new ((*trklArray)[trklArray->GetEntries()]) AliTRDtrackletMCM(*trkl);
+     }
+     return kTRUE;
+   }
+ }
+ return kFALSE;
+}
+
+Bool_t AliTRDclusterizer::ReadTracks()
+{
+  //
+  // Reads simulated GTU tracks from the input aliroot file
+  //
+
+  AliRunLoader *runLoader = AliRunLoader::Instance();
+
+  if (!runLoader) {
+    AliError("No run loader available");
+    return kFALSE;
+  }
+
+  AliLoader* loader = runLoader->GetLoader("TRDLoader");
+
+  AliDataLoader *trackLoader = loader->GetDataLoader("gtutracks");
+  if (!trackLoader) {
+      return kFALSE;
+  }
+
+  trackLoader->Load();
+
+  TTree *trackTree = trackLoader->Tree();
+  if (!trackTree) {
+    return kFALSE;
+  }
+
+  TClonesArray *trackArray = TracksArray();
+  AliTRDtrackGTU *trk = 0x0;
+  trackTree->SetBranchAddress("TRDtrackGTU", &trk);
+  for (Int_t iTrack = 0; iTrack < trackTree->GetEntries(); iTrack++) {
+    trackTree->GetEntry(iTrack);
+    new ((*trackArray)[trackArray->GetEntries()]) AliESDTrdTrack(*(trk->CreateTrdTrack()));
+  }
+
+  return kTRUE;
+}
+
 //_____________________________________________________________________________
 Bool_t AliTRDclusterizer::MakeClusters()
 {
@@ -592,8 +672,7 @@ Bool_t AliTRDclusterizer::Raw2ClustersChamber(AliRawReader *rawReader)
   //}
 
   // register tracklet array for output
-  if (fReconstructor->IsProcessingTracklets())
-    fRawStream->SetTrackletArray(TrackletsArray());
+  fRawStream->SetTrackletArray(TrackletsArray("AliTRDtrackletMCM"));
   fRawStream->SetTrackArray(TracksArray());
 
   UInt_t det = 0;
@@ -1350,19 +1429,23 @@ TClonesArray *AliTRDclusterizer::RecPoints()
 }
 
 //_____________________________________________________________________________
-TClonesArray *AliTRDclusterizer::TrackletsArray(
+TClonesArray *AliTRDclusterizer::TrackletsArray(const TString &trkltype)
 {
   //
   // Returns the array of on-line tracklets
   //
 
-  if (!fTracklets && fReconstructor->IsProcessingTracklets()) {
-    fTracklets = new TClonesArray("AliTRDtrackletWord", 200);
-    //SetClustersOwner(kTRUE);
-    //AliTRDReconstructor::SetTracklets(0x0);
+  if (trkltype.Length() != 0) {
+    if (!fTracklets) {
+      fTracklets = new TClonesArray(trkltype, 200);
+  }
+    else if (TClass::GetClass(trkltype.Data()) != fTracklets->GetClass()){
+      fTracklets->Delete();
+      delete fTracklets;
+      fTracklets = new TClonesArray(trkltype, 200);
+    }
   }
   return fTracklets;
-
 }
 
 //_____________________________________________________________________________
index 83ff093c05a35c64f02e3e304878e1bd60bab7b2..55a97d84337d3647727b3b7f29804abfd1bba943 100644 (file)
@@ -83,10 +83,13 @@ class AliTRDclusterizer : public TNamed
   Bool_t   ReadDigits(AliRawReader *rawReader);
   Bool_t   ReadDigits(TTree *digitsTree);
 
+  Bool_t   ReadTracklets();
+  Bool_t   ReadTracks();
+
   Bool_t   WriteClusters(Int_t det);
   void     ResetRecPoints();
   virtual TClonesArray    *RecPoints();
-  virtual TClonesArray    *TrackletsArray();
+  virtual TClonesArray    *TrackletsArray(const TString &trkltype = "");
   virtual TClonesArray    *TracksArray();
 
   Bool_t   Raw2Clusters(AliRawReader *rawReader);
index 1aee804f1129537f70557f9fd37f77ba013c990f..f8fddf60c1867e2dba953b990ec62c2a0c9105e9 100644 (file)
@@ -556,122 +556,247 @@ Int_t AliTRDrawStream::DecodeGTUtracks()
   // decode GTU track words
   // this depends on the hardware revision of the SMU
 
-  AliDebug(1, DumpRaw(Form("GTU tracks (hw rev %i)", fCurrHwRev),
+  Int_t sector = fCurrEquipmentId-kDDLOffset;
+
+  if ((sector < 0) || (sector > 17)) {
+    AliError(Form("Invalid sector %i for GTU tracks", sector));
+    return -1;
+  }
+
+  AliDebug(1, DumpRaw(Form("GTU tracks in sector %2i (hw rev %i)", sector, fCurrHwRev),
                      fPayloadCurr + 4, 10, 0xffe0ffff));
 
   if (fCurrHwRev < 1772) {
-    UInt_t trackWord[2] = { 0, 0 };
+    UInt_t    fastWord;                // fast trigger word
+    ULong64_t trackWord = 0;   // extended track word
     Int_t stack = 0;
     Int_t idx = 0;
     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
-      if (fPayloadCurr[iWord] == 0x10000000) {
+      if (fPayloadCurr[iWord] == 0x10000000) { // stack boundary marker
         stack++;
         idx = 0;
       }
       else {
         if ((idx == 0) &&
            ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
-         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fPayloadCurr[iWord]));
+         fastWord = fPayloadCurr[iWord];
+         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
          continue;
         }
         else if ((idx & 0x1) == 0x1) {
-         trackWord[1] = fPayloadCurr[iWord];
-         AliDebug(1,Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
-         // if (fTracks)
-         //   new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-kDDLOffset);
-        }
+         trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
+         AliDebug(1,Form("track debug word: 0x%016llx", trackWord));
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+
+           trk->SetSector(sector);
+           trk->SetStack((trackWord >> 60) & 0x7);
+           trk->SetA(0);
+           trk->SetB(0);
+           trk->SetPID(0);
+           trk->SetLayerMask((trackWord >> 16) & 0x3f);
+           trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
+           trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
+           trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
+           trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
+           trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
+           trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
+
+           trk->SetFlags(0);
+           trk->SetReserved(0);
+
+           Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
+           if (TMath::Abs(pt) > 0.1) {
+             trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
+           }
+         }
+       }
         else {
-         trackWord[0] = fPayloadCurr[iWord];
+         trackWord = fPayloadCurr[iWord];
         }
         idx++;
       }
     }
   }
   else if (fCurrHwRev < 1804) {
-    UInt_t trackWord[2] = { 0, 0 };
+    UInt_t    fastWord;                // fast trigger word
+    ULong64_t trackWord = 0;   // extended track word
     Int_t stack = 0;
     Int_t idx = 0;
     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
-      if (fPayloadCurr[iWord] == 0xffe0ffff) {
+      if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
         stack++;
         idx = 0;
       }
       else {
         if ((idx == 0) &&
            ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
-         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fPayloadCurr[iWord]));
+         fastWord = fPayloadCurr[iWord];
+         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
          continue;
         }
         else if ((idx & 0x1) == 0x1) {
-         trackWord[1] = fPayloadCurr[iWord];
-         AliDebug(1, Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
-         Float_t pt = (trackWord[0] & 0x8000) ? -1. * ((~(trackWord[0] & 0xffff)&0xffff) + 1)/128. : (trackWord[0] & 0xffff)/128.;
-         AliDebug(1, Form("pt = %f", pt));
-         // if (fTracks) {
-         //   AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-kDDLOffset);
-         //   if (TMath::Abs(pt) > 0.1) {
-         //     trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
-         //   }
-         //   trk->SetStack((trackWord[1] >> 28) & 0x7);
-         // }
+         trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
+         AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+
+           trk->SetSector(fCurrEquipmentId-kDDLOffset);
+           trk->SetStack((trackWord >> 60) & 0x7);
+           trk->SetA(0);
+           trk->SetB(0);
+           trk->SetPID(0);
+           trk->SetLayerMask((trackWord >> 16) & 0x3f);
+           trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
+           trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
+           trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
+           trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
+           trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
+           trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
+
+           trk->SetFlags(0);
+           trk->SetReserved(0);
+
+           Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
+           if (TMath::Abs(pt) > 0.1) {
+             trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
+           }
+         }
         }
         else {
-         trackWord[0] = fPayloadCurr[iWord];
+         trackWord = fPayloadCurr[iWord];
         }
         idx++;
       }
     }
   }
   else if (fCurrHwRev < 1819) {
-    UInt_t trackWord[2];
+    UInt_t    fastWord;                // fast trigger word
+    ULong64_t trackWord = 0;   // extended track word
     Int_t stack = 0;
     Int_t idx = 0;
     for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
-      if (fPayloadCurr[iWord] == 0xffe0ffff) {
+      if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
        stack++;
        idx = 0;
       }
       else {
        if ((idx == 0) &&
            ((fPayloadCurr[iWord] & 0xfffff0f0) == 0x13370000)) {
-         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fPayloadCurr[iWord]));
+         fastWord = fPayloadCurr[iWord];
+         AliDebug(1, Form("stack %i: fast trigger word: 0x%08x", stack, fastWord));
          continue;
        }
        else if ((idx & 0x1) == 0x1) {
-         trackWord[idx&0x1] = fPayloadCurr[iWord];
-         AliDebug(1, Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
-         printf("%4i %2i %i ",
-                fRawReader->GetEventIndex(),
-                fCurrEquipmentId-kDDLOffset, (trackWord[1] >> 28) & 0x7);
-         Float_t pt = (trackWord[0] & 0x8000) ? -1. * ((~(trackWord[0] & 0xffff)&0xffff) + 1)/128. : (trackWord[0] & 0xffff)/128.;
-         printf("%+7.2f ", pt);
-         printf("%i%i%i%i%i%i ", ((trackWord[0] >> 21) & 0x1),
-                ((trackWord[0] >> 20) & 0x1),
-                ((trackWord[0] >> 19) & 0x1),
-                ((trackWord[0] >> 18) & 0x1),
-                ((trackWord[0] >> 17) & 0x1),
-                ((trackWord[0] >> 16) & 0x1));
-         printf("0x%08x%08x\n", trackWord[1], trackWord[0]);
-         // if (fTracks) {
-         //   AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-kDDLOffset);
-         //   if (TMath::Abs(pt) > 0.1) {
-         //     trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
-         //   }
-         //   trk->SetStack((trackWord[1] >> 28) & 0x7);
-         // }
+         trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
+         AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
+
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+
+           trk->SetSector(fCurrEquipmentId-kDDLOffset);
+           trk->SetStack((trackWord >> 60) & 0x7);
+           trk->SetA(0);
+           trk->SetB(0);
+           // trk->SetPt(((trackWord & 0xffff) ^ 0x8000) - 0x8000);
+           trk->SetPID(0);
+           trk->SetLayerMask((trackWord >> 16) & 0x3f);
+           trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
+           trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
+           trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
+           trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
+           trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
+           trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
+
+           trk->SetFlags(0);
+           trk->SetReserved(0);
+
+           Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
+           if (TMath::Abs(pt) > 0.1) {
+             trk->SetA((Int_t) (0.15*51625./100./trk->Pt() / 160e-4 * 2));
+           }
+         }
        }
        else {
-         trackWord[idx&0x1] = fPayloadCurr[iWord];
+         trackWord = fPayloadCurr[iWord];
        }
        idx++;
       }
     }
   }
   else if (fCurrHwRev < 1860) {
-    AliError(Form("unsupported hardware rev %i", fCurrHwRev));
+    UInt_t    fastWord;                // fast trigger word
+    ULong64_t trackWord = 0;   // extended track word
+    Int_t stack = 0;
+    Int_t idx = 0;
+    Bool_t upperWord = kFALSE;
+    Int_t word = 0;
+    for (UInt_t iWord = 4; iWord < fCurrSmHeaderSize; iWord++) {
+      if (fPayloadCurr[iWord] == 0xffe0ffff) { // stack boundary marker
+        stack++;
+        idx = 0;
+       upperWord = kFALSE;
+      }
+      else {
+       // assemble the 32-bit words out of 16-bit blocks
+       if (upperWord) {
+         word |= (fPayloadCurr[iWord] & 0xffff0000);
+         upperWord = kFALSE;
+       }
+       else {
+         // lower word is read first
+         word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
+         upperWord = kTRUE;
+         continue;
+       }
+
+        if ((word & 0xffff0008) == 0x13370008) {
+         fastWord = word;
+         AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, fastWord));
+         continue;
+        }
+        else if ((idx & 0x1) == 0x1) {
+         trackWord |= ((ULong64_t) word) << 32;
+         AliDebug(1, Form("track debug word: 0x%016llx", trackWord));
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+
+           trk->SetSector(fCurrEquipmentId-kDDLOffset);
+           trk->SetStack((trackWord >> 60) & 0x7);
+           trk->SetA(0);
+           trk->SetB(0);
+           trk->SetPID(0);
+           trk->SetLayerMask((trackWord >> 16) & 0x3f);
+           trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
+           trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
+           trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
+           trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
+           trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
+           trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
+
+           trk->SetFlags(0);
+           trk->SetReserved(0);
+
+           Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
+           if (TMath::Abs(pt) > 0.1) {
+             trk->SetA((Int_t) (0.15*51625./100./pt / 160e-4 * 2));
+           }
+         }
+        }
+        else {
+         trackWord = word;
+        }
+        idx++;
+      }
+    }
+
   }
   else {
-    UInt_t trackWord[2] = { 0, 0 };
+    ULong64_t trackWord = 0; // this is the debug word
     Int_t stack = 0;
     Int_t idx = 0;
     Bool_t upperWord = kFALSE;
@@ -704,18 +829,35 @@ Int_t AliTRDrawStream::DecodeGTUtracks()
          continue;
        }
         else if ((idx & 0x1) == 0x1) {
-         trackWord[1] = word;
-         AliDebug(1, Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
-         // if (fTracks) {
-         //   AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-kDDLOffset);
-         //   if (TMath::Abs(trk->GetPt()) > 0.1) {
-         //     trk->SetA((Int_t) (0.15*51625./100./trk->GetPt() / 160e-4 * 2));
-         //   }
-         //   trk->SetStack((trackWord[1] >> 28) & 0x7);
-         // }
+         trackWord |= ((ULong64_t) word) << 32;
+         AliDebug(1, Form("track debug word: 0x%16llx", trackWord));
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+           trk->SetSector(fCurrEquipmentId-kDDLOffset);
+           trk->SetStack((trackWord >> 60) & 0x7);
+           trk->SetA(0);
+           trk->SetB(0);
+           trk->SetPID(0);
+           trk->SetLayerMask((trackWord >> 16) & 0x3f);
+           trk->SetTrackletIndex((trackWord >> 22) & 0x3f, 0);
+           trk->SetTrackletIndex((trackWord >> 28) & 0x3f, 1);
+           trk->SetTrackletIndex((trackWord >> 34) & 0x3f, 2);
+           trk->SetTrackletIndex((trackWord >> 40) & 0x3f, 3);
+           trk->SetTrackletIndex((trackWord >> 46) & 0x3f, 4);
+           trk->SetTrackletIndex((trackWord >> 52) & 0x3f, 5);
+
+           trk->SetFlags(0);
+           trk->SetReserved(0);
+
+           Float_t pt = (((Int_t) (trackWord & 0xffff) ^ 0x8000) - 0x8000)/128.;
+           if (TMath::Abs(pt) > 0.1) {
+             trk->SetA(-(Int_t) (0.15*51625./100./pt / 160e-4 * 2));
+           }
+         }
         }
         else {
-         trackWord[0] = word;
+         trackWord = word;
         }
         idx++;
       }
@@ -732,50 +874,87 @@ Int_t AliTRDrawStream::ReadTrackingHeader(Int_t stack)
 
   fCurrTrkHeaderIndexWord[stack] = *fPayloadCurr;
   fCurrTrkHeaderSize[stack]      = ((*fPayloadCurr) >> 16) & 0x3ff;
-  fPayloadCurr++;
 
-  AliDebug(1, Form("tracking header index word: 0x%08x, size: %i\n",
-                  fCurrTrkHeaderIndexWord[stack], fCurrTrkHeaderSize[stack]));
+  AliDebug(1, Form("tracking header index word: 0x%08x, size: %i (hw rev: %i)",
+                  fCurrTrkHeaderIndexWord[stack], fCurrTrkHeaderSize[stack], fCurrHwRev));
+  Int_t trackingTime = *fPayloadCurr & 0x3ff;
+
+  fPayloadCurr++;
 
   // data words
-  UInt_t trackWord[2] = { 0, 0 };
+  ULong64_t trackWord = 0;
   Int_t idx = 0;
-  Bool_t upperWord = kFALSE;
-  Int_t word = 0;
+  Int_t trackIndex = fTracks ? fTracks->GetEntriesFast() : -1;
+
   for (UInt_t iWord = 0; iWord < fCurrTrkHeaderSize[stack]; iWord++) {
-    // assemble the 32-bit words out of 16-bit blocks
-    if (upperWord) {
-      word |= (fPayloadCurr[iWord] & 0xffff0000);
-      upperWord = kFALSE;
-    }
-    else {
-      // lower word is read first
-      word = (fPayloadCurr[iWord] & 0xffff0000) >> 16;
-      upperWord = kTRUE;
-      continue;
-    }
 
-    if ((word & 0xffff0008) == 0x13370008) {
-      AliDebug(1, Form("stack %i: fast track word: 0x%08x", stack, word));
-      continue;
-    }
-    else if ((word & 0xffff0010) == 0x13370010) {
-      AliDebug(1, Form("stack %i: tracking done word: 0x%08x", stack, word));
-      continue;
-    }
-    else if ((idx & 0x1) == 0x1) {
-      trackWord[1] = word;
-      AliDebug(1, Form("track debug word: 0x%08x%08x", trackWord[1], trackWord[0]));
-      // if (fTracks) {
-      //       AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()]) AliESDTrdTrack(0, 0, trackWord[0], trackWord[1], fCurrEquipmentId-kDDLOffset);
-      //       if (TMath::Abs(trk->GetPt()) > 0.1) {
-      //         trk->SetA((Int_t) (0.15*51625./100./trk->GetPt() / 160e-4 * 2));
-      //       }
-      //       trk->SetStack((trackWord[1] >> 28) & 0x7);
-      // }
+    if (!(idx & 0x1)) {
+      // first part of 64-bit word
+      trackWord = fPayloadCurr[iWord];
     }
     else {
-      trackWord[0] = word;
+      trackWord |= ((ULong64_t) fPayloadCurr[iWord]) << 32;
+
+      if (trackWord & (1ul << 63)) {
+       if ((trackWord & (0x3ful << 56)) != 0) {
+         // track word
+         AliDebug(2, Form("track word: 0x%016llx", trackWord));
+
+         if (fTracks) {
+           AliESDTrdTrack *trk = new ((*fTracks)[fTracks->GetEntriesFast()])
+             AliESDTrdTrack();
+
+           trk->SetSector(fCurrEquipmentId-kDDLOffset);
+           trk->SetLayerMask((trackWord >> 56) & 0x3f);
+           trk->SetA( (((trackWord >> 38) & 0x3ffff) ^ 0x20000) - 0x20000);
+           trk->SetB( (((trackWord >> 20) & 0x3ffff) ^ 0x20000) - 0x20000);
+           trk->SetC( (((trackWord >> 8)  &  0xffff) ^  0x8000) -  0x8000);
+           trk->SetPID((trackWord >>  0) & 0xff);
+           trk->SetStack(stack);
+
+           // now compare the track word with the one generated from the ESD information
+           if (trackWord != trk->GetTrackWord(0)) {
+             AliError(Form("track word 0x%016llx does not match the read one 0x%016llx",
+                           trk->GetTrackWord(0), trackWord));
+           }
+         }
+       }
+       else {
+         // done marker (so far only used to set trigger flag)
+
+         AliDebug(2, Form("seg / stack / first / last / done / index : %i %i %lli %lli %lli %i",
+                          fCurrEquipmentId - kDDLOffset, stack,
+                          (trackWord >> 20) & 0x3ff,
+                          (trackWord >> 10) & 0x3ff,
+                          (trackWord >>  0) & 0x3ff,
+                          trackingTime));
+       }
+      }
+      else {
+       // extended track word
+       AliDebug(2, Form("extended track word: 0x%016llx", trackWord));
+
+       if (fTracks) {
+         AliESDTrdTrack *trk = (AliESDTrdTrack*) (*fTracks)[trackIndex];
+
+         trk->SetFlags((trackWord >> 52) & 0x7ff);
+         trk->SetReserved((trackWord >> 49) & 0x7);
+         trk->SetY((trackWord >> 36) & 0x1fff);
+         trk->SetTrackletIndex((trackWord >>  0) & 0x3f, 0);
+         trk->SetTrackletIndex((trackWord >>  6) & 0x3f, 1);
+         trk->SetTrackletIndex((trackWord >> 12) & 0x3f, 2);
+         trk->SetTrackletIndex((trackWord >> 18) & 0x3f, 3);
+         trk->SetTrackletIndex((trackWord >> 24) & 0x3f, 4);
+         trk->SetTrackletIndex((trackWord >> 30) & 0x3f, 5);
+
+         if (trackWord != trk->GetExtendedTrackWord(0)) {
+           AliError(Form("extended track word 0x%016llx does not match the read one 0x%016llx",
+                         trk->GetExtendedTrackWord(0), trackWord));
+           }
+
+         trackIndex++;
+       }
+      }
     }
     idx++;
   }
@@ -1307,16 +1486,17 @@ Int_t AliTRDrawStream::ReadZSData()
       }
 
       adcwc = 0;
-      AliDebug(3, Form("Now reading %i words for channel %2i", timebins / 3, channelno));
+      Int_t nADCwords = (timebins + 2) / 3;
+      AliDebug(3, Form("Now reading %i words for channel %2i", nADCwords, channelno));
       Int_t adccol = adccoloff - channelno;
       Int_t padcol = padcoloff - channelno;
 //      if (adccol < 3 || adccol > 165)
 //     AliInfo(Form("writing channel %i of det %3i %i:%2i to adcrow/-col: %i/%i padcol: %i",
 //                  channelno, fCurrHC/2, fCurrRobPos, fCurrMcmPos, row, adccol, padcol));
 
-      while (adcwc < timebins / 3 &&
-            *(fPayloadCurr) != fgkDataEndmarker &&
-            fPayloadCurr - fPayloadStart < fPayloadSize) {
+      while ((adcwc < nADCwords) &&
+            (*(fPayloadCurr) != fgkDataEndmarker) &&
+            (fPayloadCurr - fPayloadStart < fPayloadSize)) {
        int check = 0x3 & *fPayloadCurr;
        if (channelno % 2 != 0) { // odd channel
          if (check != 0x2 && channelno < 21) {
@@ -1348,7 +1528,7 @@ Int_t AliTRDrawStream::ReadZSData()
        fPayloadCurr++;
       }
 
-      if (adcwc != timebins / 3)
+      if (adcwc != nADCwords)
        MCMError(kAdcDataAbort);
 
       // adding index
@@ -1467,12 +1647,13 @@ Int_t AliTRDrawStream::ReadNonZSData()
       currentTimebin = 0;
 
       adcwc = 0;
-      AliDebug(2, Form("Now looking %i words", timebins / 3));
+      Int_t nADCwords = (timebins + 2) / 3;
+      AliDebug(2, Form("Now looking %i words", nADCwords));
       Int_t adccol = adccoloff - channelno;
       Int_t padcol = padcoloff - channelno;
-      while (adcwc < timebins / 3 &&
-            *(fPayloadCurr) != fgkDataEndmarker &&
-            fPayloadCurr - fPayloadStart < fPayloadSize) {
+      while ((adcwc < nADCwords) &&
+            (*(fPayloadCurr) != fgkDataEndmarker) &&
+            (fPayloadCurr - fPayloadStart < fPayloadSize)) {
        int check = 0x3 & *fPayloadCurr;
        if (channelno % 2 != 0) { // odd channel
          if (check != 0x2 && channelno < 21) {
@@ -1504,7 +1685,7 @@ Int_t AliTRDrawStream::ReadNonZSData()
        fPayloadCurr++;
       }
 
-      if (adcwc != timebins / 3)
+      if (adcwc != nADCwords)
        MCMError(kAdcDataAbort);
 
       // adding index
index f200305bf89cb622b9229ecfcbb08949310e0dd6..9bde53aa8b9d7cac5c5df3d44fe608569887ce31 100644 (file)
@@ -32,6 +32,7 @@ class AliTRDtrackletBase : public TObject {
     virtual Bool_t   CookPID() = 0;\r
 \r
     virtual Int_t    GetDetector() const = 0 ;\r
+    virtual Int_t    GetHCId() const { return 2 * GetDetector() + (GetYbin() > 0 ? 1 : 0); }\r
 \r
     virtual Float_t  GetX() const  = 0;\r
     virtual Float_t  GetY() const  = 0;\r
index ac23edc8c31cac0a14dea8453391480e316931db..ca0409329e8472770e9b4ef6001a0d687738d739 100644 (file)
@@ -34,7 +34,7 @@ class AliTRDtrackletWord : public AliTRDtrackletBase {
 
   // ----- Getters for offline corresponding values -----
   Bool_t CookPID() { return kFALSE; }
-  Double_t GetPID(Int_t /* is */) const { return 0; }
+  Double_t GetPID(Int_t /* is */) const { return (Double_t) GetPID()/256.; }
   Int_t GetDetector() const { return fHCId / 2; }
   Int_t GetHCId() const { return fHCId; }
   Float_t GetdYdX() const { return (GetdY() * 140e-4 / 3.); }
@@ -53,7 +53,7 @@ class AliTRDtrackletWord : public AliTRDtrackletBase {
  protected:
   Int_t fHCId;                  // half-chamber ID
   UInt_t fTrackletWord;                // tracklet word: PID | Z | deflection length | Y
-                               //          bits:  12   4            7          13
+                               //          bits:   8   4            7          13
   static AliTRDgeometry *fgGeo;  // pointer to TRD geometry for coordinate calculations
 
  private: