]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Raw data reconstruction (B.Vulpescu)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 8 Aug 2005 12:25:36 +0000 (12:25 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 8 Aug 2005 12:25:36 +0000 (12:25 +0000)
TRD/AliTRDDigits2Clusters.C [new file with mode: 0644]
TRD/AliTRDDigits2Raw.C [new file with mode: 0644]
TRD/AliTRDDigitsDDL2Clusters.C [new file with mode: 0644]
TRD/AliTRDRaw2Digits.C
TRD/AliTRDRawStream.cxx
TRD/AliTRDRawStream.h
TRD/AliTRDReconstructor.cxx
TRD/AliTRDReconstructor.h
TRD/AliTRDclusterizerV1.cxx
TRD/AliTRDclusterizerV1.h
TRD/AliTRDrawData.cxx

diff --git a/TRD/AliTRDDigits2Clusters.C b/TRD/AliTRDDigits2Clusters.C
new file mode 100644 (file)
index 0000000..f98a030
--- /dev/null
@@ -0,0 +1,54 @@
+void AliTRDDigits2Clusters() 
+{
+
+///////////////////////////////////////////////////////////////////////// 
+//
+// Creates cluster from the digit information. 
+//
+///////////////////////////////////////////////////////////////////////// 
+
+  // Dynamically link some shared libs
+  if (gClassTable->GetID("AliRun") < 0) {
+    gROOT->LoadMacro("loadlibs.C");
+    loadlibs();
+    cout << "Loaded shared libraries" << endl;
+  }
+
+  // Input and output file names
+  Char_t *infile  = "galice.root";
+
+  // Create the clusterizer
+  AliTRDclusterizerV1 *clusterizer = 
+    new AliTRDclusterizerV1("clusterizer","Clusterizer class"); 
+
+  // Read the parameter
+  TFile *parfile = TFile::Open(infile);
+  AliTRDparameter *par = (AliTRDparameter *) parfile->Get("TRDparameter"); 
+  par->ReInit();
+  clusterizer->SetParameter(par);
+
+  // Set the parameter
+  clusterizer->SetVerbose(1);
+
+  //Number of events
+  TTree * te = (TTree*)parfile->Get("TE");
+  Int_t nev = (Int_t)te->GetEntries();
+
+  for(Int_t iev=0;iev<nev;iev++) {
+
+    // Open the AliRoot file 
+    clusterizer->Open(infile,iev);
+
+    // Load the digits
+    clusterizer->ReadDigits();
+    clusterizer->Dump();
+    // Find the cluster
+    clusterizer->MakeClusters();
+
+    // Write the cluster tree into file AliTRDclusters.root
+    clusterizer->WriteClusters(-1);
+
+  }
+
+}
diff --git a/TRD/AliTRDDigits2Raw.C b/TRD/AliTRDDigits2Raw.C
new file mode 100644 (file)
index 0000000..92bfdb2
--- /dev/null
@@ -0,0 +1,25 @@
+void AliTRDDigits2Raw() {
+
+  const char * inFile = "galice.root";
+  AliRunLoader *rl = AliRunLoader::Open(inFile,"Event","read");
+
+  rl->GetEvent(0);
+  AliLoader *trdloader=rl->GetLoader("TRDLoader");
+  trdloader->LoadDigits();
+  TTree *digitsTree=trdloader->TreeD();
+
+  digitsTree->Print();
+
+  AliTRDrawData *raw = new AliTRDrawData();
+  raw->SetDebug(2);
+
+  if (raw->Digits2Raw(digitsTree)) {
+    printf("Digits2Raw done! \n");
+  } else {
+    printf("Digits2Raw returned FALSE! \n");
+  }
+
+  delete rl;
+
+}
+
diff --git a/TRD/AliTRDDigitsDDL2Clusters.C b/TRD/AliTRDDigitsDDL2Clusters.C
new file mode 100644 (file)
index 0000000..cc59194
--- /dev/null
@@ -0,0 +1,56 @@
+void AliTRDDigitsDDL2Clusters() 
+{
+
+///////////////////////////////////////////////////////////////////////// 
+//
+// Creates cluster from the digit information. 
+//
+///////////////////////////////////////////////////////////////////////// 
+
+  // Dynamically link some shared libs
+  if (gClassTable->GetID("AliRun") < 0) {
+    gROOT->LoadMacro("loadlibs.C");
+    loadlibs();
+    cout << "Loaded shared libraries" << endl;
+  }
+
+  // Input and output file names
+  Char_t *infile  = "galice.root";
+
+  // Create the clusterizer
+  AliTRDclusterizerV1 *clusterizer = 
+    new AliTRDclusterizerV1("clusterizer","Clusterizer class"); 
+
+  // Read the parameter
+  TFile *parfile = TFile::Open(infile);
+  AliTRDparameter *par = (AliTRDparameter *) parfile->Get("TRDparameter"); 
+  par->ReInit();
+  clusterizer->SetParameter(par);
+
+  // Set the parameter
+  clusterizer->SetVerbose(1);
+
+  //Number of events
+  TTree * te = (TTree*)parfile->Get("TE");
+  Int_t nev = (Int_t)te->GetEntries();
+
+  for(Int_t iev=0;iev<nev;iev++) {
+
+    AliRawReaderFile rawReader(iev);
+
+    // Open the AliRoot file 
+    clusterizer->Open(infile,iev);
+
+    // Load the digits
+    clusterizer->ReadDigits(&rawReader);
+    clusterizer->Dump();
+    // Find the cluster
+    clusterizer->MakeClusters();
+
+    // Write the cluster tree into file AliTRDclusters.root
+    clusterizer->WriteClusters(-1);
+
+  }
+
+}
index 0a950b48f7e3eced2f26f4fe477bf0dc54309e18..70974a2b6e7b177dd3aedf78ace179709d6d743c 100644 (file)
@@ -15,8 +15,9 @@ void AliTRDRaw2Digits(Int_t iEvent = 0, Int_t iDet = 0)
 {
 
   AliTRDrawData *raw = new AliTRDrawData();
-  raw->SetDebug(1);
+  raw->SetDebug(2);
   AliRawReaderFile rawReader(iEvent);
+
   AliTRDdigitsManager *digitsManagerRaw = raw->Raw2Digits(&rawReader);
 
   // The geometry object
@@ -37,10 +38,12 @@ void AliTRDRaw2Digits(Int_t iEvent = 0, Int_t iDet = 0)
   Int_t  rowMax = par->GetRowMax(iPla,iCha,iSec);
   Int_t  colMax = par->GetColMax(iPla);
   Int_t timeMax = par->GetTimeMax();
+  Int_t timeTotal = par->GetTimeTotal();
   cout << "Geometry: rowMax = "  <<  rowMax
                 << " colMax = "  <<  colMax
-                << " timeMax = " << timeMax << endl;
-  AliTRDmatrix *matrix = new AliTRDmatrix(rowMax,colMax,timeMax,iSec,iCha,iPla);
+                << " timeMax = " << timeMax
+                << " timeTotal = " << timeTotal << endl;
+  AliTRDmatrix *matrix = new AliTRDmatrix(rowMax,colMax,timeTotal,iSec,iCha,iPla);
 
   AliRunLoader* rl = AliRunLoader::Open("galice.root");
   AliLoader* loader = rl->GetLoader("TRDLoader");
@@ -58,7 +61,7 @@ void AliTRDRaw2Digits(Int_t iEvent = 0, Int_t iDet = 0)
   Int_t DigAmpRaw, DigAmp;
 
   // Loop through the detector pixel
-  for (Int_t time = 0; time < timeMax; time++) {
+  for (Int_t time = 0; time < timeTotal; time++) {
     for (Int_t  col = 0;  col <  colMax;  col++) {
       for (Int_t  row = 0;  row <  rowMax;  row++) {
 
@@ -89,6 +92,8 @@ void AliTRDRaw2Digits(Int_t iEvent = 0, Int_t iDet = 0)
   TCanvas *c1 = new TCanvas("c1","Canvas 1",10,10,600,500);
   DigDiff->Draw();
 
+  printf("Number of digits expected = %d \n",timeTotal*rowMax*colMax);
+
 }
 
 
index ef2031b3e9343558d15dd76140f4de2fc7ad248f..25406a291a56a55bf337a882bde36507e2b8da55 100644 (file)
 
 #include "AliTRDRawStream.h"
 #include "AliRawReader.h"
+#include "AliTRDparameter.h"
 
 ClassImp(AliTRDRawStream)
 
 
-AliTRDRawStream::AliTRDRawStream(AliRawReader* rawReader) :
+AliTRDRawStream::AliTRDRawStream(AliRawReader* rawReader, AliTRDparameter* parameter) :
   fRawReader(rawReader),
-  fTimeMax(15),
+  fTimeTotal(parameter->GetTimeTotal()),
   fCount(0),
   fDetector(-1),
   fPrevDetector(-1),
@@ -54,7 +55,7 @@ AliTRDRawStream::AliTRDRawStream(AliRawReader* rawReader) :
 AliTRDRawStream::AliTRDRawStream(const AliTRDRawStream& stream) :
   TObject(stream),
   fRawReader(NULL),
-  fTimeMax(15),
+  fTimeTotal(0),
   fCount(0),
   fDetector(-1),
   fPrevDetector(-1),
@@ -131,6 +132,12 @@ Bool_t AliTRDRawStream::Next()
        return kFALSE;
       }
       fCount += (UInt_t(data) << 8);
+      if (!fRawReader->ReadNextChar(data)) {
+        Error("Next", "could not read number of bytes");
+        fCount = -1;
+        return kFALSE;
+      }
+      fCount += (UInt_t(data) << 16);
 
       // read the number of active pads
       if (!fRawReader->ReadNextChar(data)) {
@@ -146,18 +153,12 @@ Bool_t AliTRDRawStream::Next()
       }
       fNPads += (UInt_t(data) << 8);
 
-      // read the empty byte
-      if (!fRawReader->ReadNextChar(data)) {
-       Error("Next", "could not read fill byte");
-       fCount = -1;
-       return kFALSE;
-      }
+      fTime = fTimeTotal;
 
-      fTime = fTimeMax;
     }
 
     // read the pad row and column number
-    if ((fTime >= fTimeMax) && (fCount > 2)) {
+    if ((fTime >= fTimeTotal) && (fCount > 2)) {
       if (!fRawReader->ReadNextChar(data)) {
        Error("Next", "could not read row number");
        fCount = -1;
index 1945db10c20e635bf97c92f4b3f6efa9b02b84f5..ad00848a28f91d8634731170732db9595a0252a5 100644 (file)
 #include <TObject.h>
 
 class AliRawReader;
-
+class AliTRDparameter;
 
 class AliTRDRawStream: public TObject {
   public :
-    AliTRDRawStream(AliRawReader* rawReader);
+    AliTRDRawStream(AliRawReader* rawReader, AliTRDparameter* parameter);
     virtual ~AliTRDRawStream();
 
     virtual Bool_t   Next();
@@ -44,7 +44,7 @@ class AliTRDRawStream: public TObject {
 
     AliRawReader*    fRawReader;    // object for reading the raw data
 
-    Int_t            fTimeMax;      // maximal time bin
+    Int_t            fTimeTotal;    // total time bins
     Int_t            fCount;        // counter of bytes to be read for current detector
 
     Int_t            fDetector;     // index of current detector
index 82a9e40a9ad267c1a3ff521dba4929e2de0e0a31..291ea9dc98d9eb45ec18bf0b84234560d8bc37a5 100644 (file)
@@ -29,7 +29,8 @@
 #include "AliTRDtracker.h"
 #include "AliTRDpidESD.h"
 #include <TFile.h>
-
+#include "AliRawReaderFile.h"
+#include "AliLog.h"
 
 ClassImp(AliTRDReconstructor)
 
@@ -63,6 +64,39 @@ void AliTRDReconstructor::Reconstruct(AliRunLoader* runLoader) const
   loader->UnloadRecPoints();
 }
 
+//_____________________________________________________________________________
+void AliTRDReconstructor::Reconstruct(AliRunLoader* runLoader,
+                                      AliRawReader* rawReader) const
+{
+// reconstruct clusters
+
+  AliInfo("Reconstruct TRD clusters from RAW data");
+
+  AliLoader *loader=runLoader->GetLoader("TRDLoader");
+  loader->LoadRecPoints("recreate");
+
+  AliTRDclusterizerV1 clusterer("clusterer", "TRD clusterizer");
+  runLoader->CdGAFile();
+  AliTRDparameter* trdParam = GetTRDparameter(runLoader); 
+  if (!trdParam) {
+    Error("Reconstruct", "no TRD parameters found");
+    return;
+  }
+  trdParam->ReInit();
+  clusterer.SetParameter(trdParam);
+  Int_t nEvents = runLoader->GetNumberOfEvents();
+
+  for (Int_t iEvent = 0; iEvent < nEvents; iEvent++) {
+    if (!rawReader->NextEvent()) break;
+    clusterer.Open(runLoader->GetFileName(), iEvent);
+    clusterer.ReadDigits(rawReader);
+    clusterer.MakeClusters();
+    clusterer.WriteClusters(-1);
+  }
+
+  loader->UnloadRecPoints();
+}
+
 //_____________________________________________________________________________
 AliTracker* AliTRDReconstructor::CreateTracker(AliRunLoader* runLoader) const
 {
index e40a1b2b70505e62a306be66eb4b08dfa59f7263..d691606994f64d625f582cf0a3e175f5b13ec93b 100644 (file)
@@ -20,7 +20,7 @@ public:
   AliTRDReconstructor(): AliReconstructor() {};
   virtual ~AliTRDReconstructor() {};
 
-  virtual void         Reconstruct(AliRunLoader*, AliRawReader*) const { };
+  virtual void         Reconstruct(AliRunLoader* runLoader, AliRawReader* rawReader) const;
   virtual void         Reconstruct(AliRawReader*, TTree*) const { };
   virtual void         Reconstruct(TTree*, TTree*) const { };
   virtual void         Reconstruct(AliRunLoader* runLoader) const;
index c5fff1f5e931a23ddb8536e8a19077d93ee39ac7..65167bcf54876d03232d1fc7053cb8c3735cbdf7 100644 (file)
@@ -29,6 +29,7 @@
 #include "AliRun.h"
 #include "AliRunLoader.h"
 #include "AliLoader.h"
+#include "AliRawReader.h"
 
 #include "AliTRDclusterizerV1.h"
 #include "AliTRDmatrix.h"
@@ -38,6 +39,7 @@
 #include "AliTRDdigitsManager.h"
 #include "AliTRDparameter.h"
 #include "AliTRDpadPlane.h"
+#include "AliTRDrawData.h"
 
 ClassImp(AliTRDclusterizerV1)
 
@@ -136,6 +138,22 @@ Bool_t AliTRDclusterizerV1::ReadDigits()
 
 }
 
+//_____________________________________________________________________________
+Bool_t AliTRDclusterizerV1::ReadDigits(AliRawReader* rawReader)
+{
+  //
+  // Reads the digits arrays from the ddl file
+  //
+
+  AliTRDrawData *raw = new AliTRDrawData();
+  raw->SetDebug(1);
+
+  fDigitsManager = raw->Raw2Digits(rawReader);
+
+  return kTRUE;
+
+}
+
 //_____________________________________________________________________________
 Bool_t AliTRDclusterizerV1::MakeClusters()
 {
index 52c7a14c1975ad649df0c18c252e04202a2dd7c9..9442584f553a9fe3b5cc5e6c1b9bcd51382f507c 100644 (file)
@@ -13,6 +13,7 @@
 
 class AliTRDdigitsManager;
 class AliTRDparameter;
+class AliRawReader;
 
 class AliTRDclusterizerV1 : public AliTRDclusterizer {
 
@@ -27,6 +28,7 @@ class AliTRDclusterizerV1 : public AliTRDclusterizer {
   virtual void     Copy(TObject &c) const;
   virtual Bool_t   MakeClusters();
   virtual Bool_t   ReadDigits();
+  virtual Bool_t   ReadDigits(AliRawReader* rawReader);
 
  protected:
 
index 573a4f455d4ec7ccf7f13af9488e487f0c37d477..0f10e8994f4b483098da4d31e65e0b9d67e5105a 100644 (file)
@@ -30,6 +30,7 @@
 #include "AliTRDdataArrayI.h"
 #include "AliTRDRawStream.h"
 #include "AliRawDataHeader.h"
+#include "AliRawReader.h"
 
 ClassImp(AliTRDrawData)
 
@@ -111,7 +112,7 @@ Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
   const Int_t kNumberOfDDLs         = 18;
   const Int_t kSubeventHeaderLength = 8;
   const Int_t kSubeventDummyFlag    = 0xBB;
-  int headerSubevent[2];
+  int headerSubevent[3];
 
   ofstream      *outputFile[kNumberOfDDLs];
   UInt_t         bHPosition[kNumberOfDDLs];
@@ -162,8 +163,8 @@ Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
     Int_t sect      = geo->GetSector(det);
     Int_t rowMax    = par->GetRowMax(plan,cham,sect);
     Int_t colMax    = par->GetColMax(plan);
-    Int_t timeMax   = par->GetTimeMax();
-    Int_t bufferMax = rowMax*colMax*timeMax;
+    Int_t timeTotal = par->GetTimeTotal();
+    Int_t bufferMax = rowMax*colMax*timeTotal;
     int  *buffer    = new int[bufferMax];
 
     npads   = 0;
@@ -182,7 +183,7 @@ Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
 
        // Check whether data exists for this pad
         Bool_t dataflag = kFALSE;
-        for (Int_t time = 0; time < timeMax; time++) {
+        for (Int_t time = 0; time < timeTotal; time++) {
           Int_t data = digits->GetDataUnchecked(row,col,time);
           if (data) {
             dataflag = kTRUE;
@@ -201,14 +202,14 @@ Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
           nbyte += 2;
 
           Int_t nzero = 0;
-          for (Int_t time = 0; time < timeMax; time++) {
+          for (Int_t time = 0; time < timeTotal; time++) {
 
             Int_t data = digits->GetDataUnchecked(row,col,time);
 
             if (!data) {
               nzero++;
               if ((nzero ==       256) || 
-                  (time  == timeMax-1)) {
+                  (time  == timeTotal-1)) {
                 *bytePtr++ = 0;
                 *bytePtr++ = nzero-1;
                 nbyte += 2;
@@ -255,9 +256,9 @@ Bool_t AliTRDrawData::Digits2Raw(TTree *digitsTree)
     *bytePtr++ = (det   >> 8);
     *bytePtr++ = (nbyte & 0xff);
     *bytePtr++ = (nbyte >> 8);
+    *bytePtr++ = (nbyte >> 16);
     *bytePtr++ = (npads & 0xff);
     *bytePtr++ = (npads >> 8);
-    *bytePtr++ = 0;
     outputFile[iDDL]->write((char*)headerPtr,kSubeventHeaderLength);
 
     // Write the buffer to the file
@@ -303,6 +304,9 @@ AliTRDdigitsManager* AliTRDrawData::Raw2Digits(AliRawReader* rawReader)
   //
 
   AliTRDdataArrayI *digits    = 0;
+  AliTRDdataArrayI *track0    = 0;
+  AliTRDdataArrayI *track1    = 0;
+  AliTRDdataArrayI *track2    = 0; 
 
   AliTRDgeometryFull *geo = new AliTRDgeometryFull();
   AliTRDparameter    *par = new AliTRDparameter("TRDparameter"
@@ -313,7 +317,7 @@ AliTRDdigitsManager* AliTRDrawData::Raw2Digits(AliRawReader* rawReader)
   digitsManager->SetDebug(fDebug);
   digitsManager->CreateArrays();
 
-  AliTRDRawStream input(rawReader);
+  AliTRDRawStream input(rawReader,par);
 
   // Loop through the digits
   while (input.Next()) {
@@ -324,6 +328,9 @@ AliTRDdigitsManager* AliTRDrawData::Raw2Digits(AliRawReader* rawReader)
     if (input.IsNewDetector()) {
 
       if (digits) digits->Compress(1,0);
+      if (track0) track0->Compress(1,0);
+      if (track1) track1->Compress(1,0);
+      if (track2) track2->Compress(1,0);
 
       if (fDebug > 2) {
        Info("Raw2Digits","Subevent header:");
@@ -337,21 +344,37 @@ AliTRDdigitsManager* AliTRDrawData::Raw2Digits(AliRawReader* rawReader)
       Int_t sect      = geo->GetSector(det);
       Int_t rowMax    = par->GetRowMax(plan,cham,sect);
       Int_t colMax    = par->GetColMax(plan);
-      Int_t timeMax   = par->GetTimeMax();
+      Int_t timeTotal = par->GetTimeTotal();
 
       // Add a container for the digits of this detector
       digits = digitsManager->GetDigits(det);
+      track0 = digitsManager->GetDictionary(det,0);
+      track1 = digitsManager->GetDictionary(det,1);
+      track2 = digitsManager->GetDictionary(det,2);
       // Allocate memory space for the digits buffer
       if (digits->GetNtime() == 0) {
-        digits->Allocate(rowMax,colMax,timeMax);
+        digits->Allocate(rowMax,colMax,timeTotal);
+        track0->Allocate(rowMax,colMax,timeTotal);
+        track1->Allocate(rowMax,colMax,timeTotal);
+        track2->Allocate(rowMax,colMax,timeTotal);
       }
 
     } 
+
     digits->SetDataUnchecked(input.GetRow(),input.GetColumn(),
                             input.GetTime(),input.GetSignal());
+    track0->SetDataUnchecked(input.GetRow(),input.GetColumn(),
+                             input.GetTime(),               -1);
+    track1->SetDataUnchecked(input.GetRow(),input.GetColumn(),
+                             input.GetTime(),               -1);
+    track2->SetDataUnchecked(input.GetRow(),input.GetColumn(),
+                             input.GetTime(),               -1);
   }
 
   if (digits) digits->Compress(1,0);
+  if (track0) track0->Compress(1,0);
+  if (track1) track1->Compress(1,0);
+  if (track2) track2->Compress(1,0);
 
   delete geo;
   delete par;