]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PMD/AliPMDClusterFinder.cxx
AliACORDEQAChecker::Check fixed
[u/mrichter/AliRoot.git] / PMD / AliPMDClusterFinder.cxx
index 298b7ce17dc1859ed67d8de79e8f2391b67f6975..7ca2996f0bbe61142e38d96da9fa5c0fcd7efeb5 100644 (file)
@@ -26,6 +26,7 @@
 #include <TTree.h>
 #include <TObjArray.h>
 #include <TClonesArray.h>
+#include <TSystem.h>
 
 #include "AliLog.h"
 #include "AliRunLoader.h"
 #include "AliPMDcluster.h"
 #include "AliPMDrecpoint1.h"
 #include "AliPMDrechit.h"
+#include "AliPMDisocell.h"
 #include "AliPMDRawStream.h"
 #include "AliPMDCalibData.h"
 #include "AliPMDPedestal.h"
 #include "AliPMDddldata.h"
+#include "AliPMDHotData.h"
+#include "AliPMDRecoParam.h"
+#include "AliPMDReconstructor.h"
 
 #include "AliDAQ.h"
 #include "AliCDBManager.h"
@@ -57,6 +62,8 @@ AliPMDClusterFinder::AliPMDClusterFinder():
   fPMDLoader(0),
   fCalibGain(GetCalibGain()),
   fCalibPed(GetCalibPed()),
+  fCalibHot(GetCalibHot()),
+  fRecoParam(0x0),
   fTreeD(0),
   fTreeR(0),
   fDigits(new TClonesArray("AliPMDdigit", 1000)),
@@ -77,6 +84,8 @@ AliPMDClusterFinder::AliPMDClusterFinder(AliRunLoader* runLoader):
   fPMDLoader(runLoader->GetLoader("PMDLoader")),
   fCalibGain(GetCalibGain()),
   fCalibPed(GetCalibPed()),
+  fCalibHot(GetCalibHot()),
+  fRecoParam(0x0),
   fTreeD(0),
   fTreeR(0),
   fDigits(new TClonesArray("AliPMDdigit", 1000)),
@@ -98,6 +107,8 @@ AliPMDClusterFinder::AliPMDClusterFinder(const AliPMDClusterFinder & finder):
   fPMDLoader(0),
   fCalibGain(GetCalibGain()),
   fCalibPed(GetCalibPed()),
+  fCalibHot(GetCalibHot()),
+  fRecoParam(0x0),
   fTreeD(0),
   fTreeR(0),
   fDigits(NULL),
@@ -124,22 +135,17 @@ AliPMDClusterFinder::~AliPMDClusterFinder()
   // Destructor
   if (fDigits)
     {
-      fDigits->Delete();
-      delete fDigits;
-      fDigits=0;
+      fDigits->Clear();
     }
   if (fRecpoints)
     {
-      fRecpoints->Delete();
-      delete fRecpoints;
-      fRecpoints=0;
+      fRecpoints->Clear();
     }
   if (fRechits)
     {
-      fRechits->Delete();
-      delete fRechits;
-      fRechits=0;
+      fRechits->Clear();
     }
+
 }
 // ------------------------------------------------------------------------- //
 
@@ -156,10 +162,21 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
   Int_t    idet;
   Float_t  clusdata[6];
 
+  AliPMDisocell *pmdiso = 0x0;
+
   TObjArray *pmdcont = new TObjArray();
+  TObjArray *pmdisocell = new TObjArray();
+
   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
 
-  pmdclust->SetEdepCut(fEcut);
+  // fetch the recoparam object from database
+  fRecoParam = AliPMDReconstructor::GetRecoParam();
+
+  if(fRecoParam == 0x0)
+    {
+       AliFatal("No Reco Param found for PMD!!!");
+    }
+
 
   fRunLoader->GetEvent(ievt);
 
@@ -203,7 +220,16 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
          xpos   = pmddigit->GetRow();
          ypos   = pmddigit->GetColumn();
          adc    = pmddigit->GetADC();
-         
+         if(xpos < 0 || xpos > 48 || ypos < 0 || ypos > 96)
+           {
+             AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
+                             xpos, ypos));
+             continue; 
+           }
+
+         // Hot cell - set the cell adc = 0
+         Float_t hotflag = fCalibHot->GetHotChannel(det,smn,xpos,ypos);
+         if (hotflag == 1) adc = 0;
          // CALIBRATION
          Float_t gain = fCalibGain->GetGainFact(det,smn,xpos,ypos);
          // printf("adc = %d gain = %f\n",adc,gain);
@@ -211,13 +237,28 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
          adc = adc*gain;
 
          //Int_t trno   = pmddigit->GetTrackNumber();
-         fCellADC[xpos][ypos] = (Double_t) adc;
+         fCellTrack[xpos][ypos] = pmddigit->GetTrackNumber();
+         fCellPid[xpos][ypos]   = pmddigit->GetTrackPid();
+         fCellADC[xpos][ypos]   = (Double_t) adc;
        }
 
       idet = det;
       ismn = smn;
-      pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
-      
+
+      // Set the minimum noise cut per module before clustering
+
+      Int_t imod = idet*24 + ismn;
+      fEcut = fRecoParam->GetNoiseCut(imod);   // default
+      // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
+      // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
+      // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
+
+      pmdclust->SetEdepCut(fEcut);
+
+
+      pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
+                       pmdisocell,pmdcont);
+
       Int_t nentries1 = pmdcont->GetEntries();
 
       AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
@@ -241,12 +282,35 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
            {
              Int_t celldataX = pmdcl->GetClusCellX(ihit);
              Int_t celldataY = pmdcl->GetClusCellY(ihit);
-             AddRecHit(celldataX, celldataY);
+             Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
+             Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
+             Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
+             AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
            }
          branch2->Fill();
          ResetRechit();
        }
-      pmdcont->Clear();
+      pmdcont->Delete();
+
+         // Added single isolated cell for offline gain calibration
+         nentries1 = pmdisocell->GetEntries();
+         AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
+         
+         for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
+           {
+             pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
+             idet = pmdiso->GetDetector();
+             ismn = pmdiso->GetSmn();
+             clusdata[0] = (Float_t) pmdiso->GetRow();
+             clusdata[1] = (Float_t) pmdiso->GetCol();
+             clusdata[2] = pmdiso->GetADC();
+             clusdata[3] = 1.;
+             clusdata[4] = -99.;
+             clusdata[5] = -99.;
+      
+             AddRecPoint(idet,ismn,clusdata);
+           }
+         pmdisocell->Delete();
       
       branch1->Fill();
       ResetRecpoint();
@@ -260,6 +324,7 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
   //   delete the pointers
   delete pmdclust;
   delete pmdcont;
+  delete pmdisocell;
     
 }
 // ------------------------------------------------------------------------- //
@@ -270,6 +335,7 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
   // Converts digits to recpoints after running clustering
   // algorithm on CPV plane and PREshower plane
   //
+  // This algorithm is called during the reconstruction from digits
 
   Int_t    det = 0,smn = 0;
   Int_t    xpos,ypos;
@@ -278,10 +344,21 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
   Int_t    idet;
   Float_t  clusdata[6];
 
+  AliPMDcluster *pmdcl = 0x0;
+  AliPMDisocell *pmdiso = 0x0;
+
   TObjArray *pmdcont = new TObjArray();
+  TObjArray *pmdisocell = new TObjArray();
   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
 
-  pmdclust->SetEdepCut(fEcut);
+  // Fetch the reco param object
+
+  fRecoParam = AliPMDReconstructor::GetRecoParam();
+  if(fRecoParam == 0x0)
+    {
+       AliFatal("No Reco Param found for PMD!!!");
+    }
+
 
   AliPMDdigit  *pmddigit;
   TBranch *branch = digitsTree->GetBranch("PMDDigit");
@@ -297,6 +374,8 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
 
   for (Int_t imodule = 0; imodule < nmodules; imodule++)
     {
+
+      Int_t totADCMod = 0;
       ResetCellADC();
       digitsTree->GetEntry(imodule); 
       Int_t nentries = fDigits->GetLast();
@@ -309,29 +388,58 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
          xpos   = pmddigit->GetRow();
          ypos   = pmddigit->GetColumn();
          adc    = pmddigit->GetADC();
-
+         if(xpos < 0 || xpos > 48 || ypos < 0 || ypos > 96)
+           {
+             AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
+                           xpos, ypos));
+             continue; 
+           }
+         
          // Pedestal Subtraction
          Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,xpos,ypos);
-         Int_t   pedrms1    = (Int_t) pedmeanrms%1000;
+         Int_t   pedrms1    = (Int_t) pedmeanrms%100;
          Float_t pedrms     = (Float_t)pedrms1/10.;
          Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
          //printf("%f %f\n",pedmean, pedrms);
 
          Float_t adc1 = adc - (pedmean + 3.0*pedrms);
 
+         // Hot cell - set the cell adc = 0
+         Float_t hotflag = fCalibHot->GetHotChannel(det,smn,xpos,ypos);
+         if (hotflag == 1) adc1 = 0;
+
          // CALIBRATION
          Float_t gain = fCalibGain->GetGainFact(det,smn,xpos,ypos);
          // printf("adc = %d gain = %f\n",adc,gain);
 
          adc = adc1*gain;
 
-         //Int_t trno   = pmddigit->GetTrackNumber();
+         fCellTrack[xpos][ypos] = pmddigit->GetTrackNumber();
+         fCellPid[xpos][ypos] = pmddigit->GetTrackPid();
          fCellADC[xpos][ypos] = (Double_t) adc;
+
+         totADCMod += (Int_t) adc;
+
        }
 
       idet = det;
       ismn = smn;
-      pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
+
+      if (totADCMod <= 0) continue;
+
+      // Set the minimum noise cut per module before clustering
+
+      Int_t imod = idet*24 + ismn;
+      fEcut = fRecoParam->GetNoiseCut(imod);       // default
+      // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
+      // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
+      // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
+
+
+      pmdclust->SetEdepCut(fEcut);
+
+      pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
+                       pmdisocell,pmdcont);
       
       Int_t nentries1 = pmdcont->GetEntries();
 
@@ -339,7 +447,7 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
 
       for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
        {
-         AliPMDcluster *pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
+         pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
          idet        = pmdcl->GetDetector();
          ismn        = pmdcl->GetSMN();
          clusdata[0] = pmdcl->GetClusX();
@@ -356,24 +464,48 @@ void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
            {
              Int_t celldataX = pmdcl->GetClusCellX(ihit);
              Int_t celldataY = pmdcl->GetClusCellY(ihit);
-             AddRecHit(celldataX, celldataY);
+             Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
+             Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
+             Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
+             AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
            }
          branch2->Fill();
          ResetRechit();
        }
-      pmdcont->Clear();
+      pmdcont->Delete();
+
+      // Added single isolated cell for offline gain calibration
+      nentries1 = pmdisocell->GetEntries();
+      AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
+      
+      for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
+       {
+         pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
+         idet = pmdiso->GetDetector();
+         ismn = pmdiso->GetSmn();
+         clusdata[0] = (Float_t) pmdiso->GetRow();
+         clusdata[1] = (Float_t) pmdiso->GetCol();
+         clusdata[2] = pmdiso->GetADC();
+         clusdata[3] = 1.;
+         clusdata[4] = -99.;
+         clusdata[5] = -99.;
+         
+         AddRecPoint(idet,ismn,clusdata);
+       }
+      pmdisocell->Delete();
       
       branch1->Fill();
       ResetRecpoint();
 
     } // modules
 
+
   ResetCellADC();
 
   //   delete the pointers
   delete pmdclust;
   delete pmdcont;
-    
+  delete pmdisocell;
 }
 // ------------------------------------------------------------------------- //
 
@@ -383,16 +515,58 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
   // Converts RAW data to recpoints after running clustering
   // algorithm on CPV and PREshower plane
   //
-  // This method is called at the time of reconstruction
+  // This method is called at the time of reconstruction from RAW data
+
+
+  AliPMDddldata *pmdddl = 0x0;
+  AliPMDcluster *pmdcl  = 0x0;
+  AliPMDisocell *pmdiso = 0x0;
 
 
   Float_t  clusdata[6];
   TObjArray pmdddlcont;
 
   TObjArray *pmdcont = new TObjArray();
+  TObjArray *pmdisocell = new TObjArray();
+
   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
 
-  pmdclust->SetEdepCut(fEcut);
+  // open the ddl file info to know the module
+  TString ddlinfofileName(gSystem->Getenv("ALICE_ROOT"));
+  ddlinfofileName += "/PMD/PMD_ddl_info.dat";
+  
+  ifstream infileddl;
+  infileddl.open(ddlinfofileName.Data(), ios::in); // ascii file
+  if(!infileddl) AliError("Could not read the ddl info file");
+
+  Int_t ddlno;
+  Int_t modno;
+  Int_t modulePerDDL;
+  Int_t moduleddl[6];
+
+  for(Int_t jddl = 0; jddl < 6; jddl++)
+    {
+      if (infileddl.eof()) break;
+      infileddl >> ddlno >> modulePerDDL;
+      moduleddl[jddl] = modulePerDDL;
+
+      if (modulePerDDL == 0) continue;
+      for (Int_t im = 0; im < modulePerDDL; im++)
+       {
+         infileddl >> modno;
+       }
+    }
+
+  infileddl.close();
+
+  // Set the minimum noise cut per module before clustering
+
+  fRecoParam = AliPMDReconstructor::GetRecoParam();
+
+  if(fRecoParam == 0x0)
+    {
+       AliFatal("No Reco Param found for PMD!!!");
+    }
 
   ResetRecpoint();
 
@@ -401,24 +575,19 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
 
   TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize); 
 
-  const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
   const Int_t kRow = 48;
   const Int_t kCol = 96;
 
   Int_t idet = 0;
   Int_t iSMN = 0;
 
+  Int_t indexDDL = -1;
+  AliPMDRawStream pmdinput(rawReader);
 
-  for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
+  while ((indexDDL = pmdinput.DdlData(&pmdddlcont)) >=0)
     {
-      if (indexDDL < 4)
-       {
-         iSMN = 6;
-       }
-      else if (indexDDL >= 4)
-       {
-         iSMN = 12;
-       }
+      iSMN = moduleddl[indexDDL];
+
       Int_t ***precpvADC;
       precpvADC = new int **[iSMN];
       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
@@ -437,18 +606,12 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
            }
        }
       ResetCellADC();
-      rawReader->Reset();
-      AliPMDRawStream pmdinput(rawReader);
-
-      rawReader->Select("PMD", indexDDL, indexDDL);
-
-      pmdinput.DdlData(indexDDL,&pmdddlcont);
 
       Int_t indexsmn = 0;
       Int_t ientries = pmdddlcont.GetEntries();
       for (Int_t ient = 0; ient < ientries; ient++)
        {
-         AliPMDddldata *pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
+         pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
          
          Int_t det = pmdddl->GetDetector();
          Int_t smn = pmdddl->GetSMN();
@@ -458,9 +621,22 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
          Int_t col = pmdddl->GetColumn();
          Int_t sig = pmdddl->GetSignal();
 
+         if(smn == -1)
+           {
+             AliError(Form("*MODULE NUMBER WRONG %d *",smn));
+             continue; 
+           }
+         if(row < 0 || row > 48 || col < 0 || col > 96)
+           {
+             AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
+                           row, col));
+
+             continue; 
+           }
+
          // Pedestal Subtraction
          Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,row,col);
-         Int_t   pedrms1    = (Int_t) pedmeanrms%1000;
+         Int_t   pedrms1    = (Int_t) pedmeanrms%100;
          Float_t pedrms     = (Float_t)pedrms1/10.;
          Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
 
@@ -469,12 +645,33 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
          // Float_t sig1 = (Float_t) sig;
          Float_t sig1 = (Float_t) sig - (pedmean + 3.0*pedrms);
 
+         // Hot cell - set the cell adc = 0
+         Float_t hotflag = fCalibHot->GetHotChannel(det,smn,row,col);
+         if (hotflag == 1) sig1 = 0;
+
          // CALIBRATION
          Float_t gain = fCalibGain->GetGainFact(det,smn,row,col);
          //printf("sig = %d gain = %f\n",sig,gain);
          sig = (Int_t) (sig1*gain);
 
-         if (indexDDL < 4)
+         if (indexDDL == 0)
+           {
+             if (det != 0)
+               AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
+                             indexDDL, det));
+             if (iSMN == 6)
+               {
+                 indexsmn = smn;
+               }
+             else if (iSMN == 12)
+               {
+                 if (smn < 6)
+                   indexsmn = smn;
+                 else if (smn >= 18 && smn < 24)
+                   indexsmn = smn-12;
+               }
+           }
+         else if (indexDDL >= 1 && indexDDL < 4)
            {
              if (det != 0)
                AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
@@ -500,34 +697,54 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
              if (det != 1)
                AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
                              indexDDL, det));
-             if (smn >= 6 && smn < 12)
-               {
-                 indexsmn = smn - 6;
-               }
-             else if (smn >= 12 && smn < 18)
+             if (smn >= 6 && smn < 18)
                {
                  indexsmn = smn - 6;
                }
            }         
+
          precpvADC[indexsmn][row][col] = sig;
        }
       
-      pmdddlcont.Clear();
+      pmdddlcont.Delete();
+
+      Int_t totAdcMod = 0;
 
       Int_t ismn = 0;
-      for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
+      for (indexsmn = 0; indexsmn < iSMN; indexsmn++)
        {
          ResetCellADC();
+         totAdcMod = 0;
          for (Int_t irow = 0; irow < kRow; irow++)
            {
              for (Int_t icol = 0; icol < kCol; icol++)
                {
+                 fCellTrack[irow][icol] = -1;
+                 fCellPid[irow][icol]   = -1;
+
                  fCellADC[irow][icol] = 
                    (Double_t) precpvADC[indexsmn][irow][icol];
+                 totAdcMod += precpvADC[indexsmn][irow][icol];
                } // row
            }     // col
          
-         if (indexDDL < 4)
+         if (indexDDL == 0)
+           {
+             if (iSMN == 6)
+               {
+                 ismn = indexsmn;
+               }
+             else if (iSMN == 12)
+               {
+                 
+                 if (indexsmn < 6)
+                   ismn = indexsmn;
+                 else if (indexsmn >= 6 && indexsmn < 12)
+                   ismn = indexsmn + 12;
+               }
+             idet = 0;
+           }
+         else if (indexDDL >= 1 && indexDDL < 4)
            {
              ismn = indexsmn + indexDDL * 6;
              idet = 0;
@@ -546,26 +763,32 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
            }
          else if (indexDDL == 5)
            {
-             if (indexsmn < 6)
-               {
-                 ismn = indexsmn + 6;
-               }
-             else if (indexsmn >= 6 && indexsmn < 12)
-               {
-                 ismn = indexsmn + 6;
-               }
+             ismn = indexsmn + 6;
              idet = 1;
            }
 
-         pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
+         if (totAdcMod <= 0) continue;
+
+         Int_t imod = idet*24 + ismn;
+
+         fEcut = fRecoParam->GetNoiseCut(imod);       // default
+         // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
+         // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
+         // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
+
+
+         pmdclust->SetEdepCut(fEcut);
+
+         pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
+                           pmdisocell,pmdcont);
+
          Int_t nentries1 = pmdcont->GetEntries();
 
          AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
 
          for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
            {
-             AliPMDcluster *pmdcl = 
-               (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
+             pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
              idet        = pmdcl->GetDetector();
              ismn        = pmdcl->GetSMN();
              clusdata[0] = pmdcl->GetClusX();
@@ -582,13 +805,37 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
                {
                  Int_t celldataX = pmdcl->GetClusCellX(ihit);
                  Int_t celldataY = pmdcl->GetClusCellY(ihit);
-                 AddRecHit(celldataX, celldataY);
+                 Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
+                 Int_t celldataPid   = pmdcl->GetClusCellPid(ihit);
+                 Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
+                 AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
                }
              branch2->Fill();
              ResetRechit();
 
            }
-         pmdcont->Clear();
+         pmdcont->Delete();
+
+
+         // Added single isolated cell for offline gain calibration
+         nentries1 = pmdisocell->GetEntries();
+         AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
+         
+         for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
+           {
+             pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
+             idet = pmdiso->GetDetector();
+             ismn = pmdiso->GetSmn();
+             clusdata[0] = (Float_t) pmdiso->GetRow();
+             clusdata[1] = (Float_t) pmdiso->GetCol();
+             clusdata[2] = pmdiso->GetADC();
+             clusdata[3] = 1.;
+             clusdata[4] = -99.;
+             clusdata[5] = -99.;
+      
+             AddRecPoint(idet,ismn,clusdata);
+           }
+         pmdisocell->Delete();
          
          branch1->Fill();
          ResetRecpoint();
@@ -601,14 +848,17 @@ void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
          for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
        }
       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
-      delete precpvADC;
+      delete [] precpvADC;
+
     } // DDL Loop
+
   
   ResetCellADC();
   
   //   delete the pointers
   delete pmdclust;
   delete pmdcont;
+  delete pmdisocell;
 
 }
 // ------------------------------------------------------------------------- //
@@ -620,12 +870,55 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
   //
 
   Float_t  clusdata[6];
+
   TObjArray pmdddlcont;
+
+  AliPMDcluster *pmdcl  = 0x0;
+  AliPMDisocell *pmdiso  = 0x0;
+
+
   TObjArray *pmdcont = new TObjArray();
+  TObjArray *pmdisocell = new TObjArray();
 
   AliPMDClustering *pmdclust = new AliPMDClusteringV1();
 
-  pmdclust->SetEdepCut(fEcut);
+  // open the ddl file info to know the module
+  TString ddlinfofileName(gSystem->Getenv("ALICE_ROOT"));
+  ddlinfofileName += "/PMD/PMD_ddl_info.dat";
+  
+  ifstream infileddl;
+  infileddl.open(ddlinfofileName.Data(), ios::in); // ascii file
+  if(!infileddl) AliError("Could not read the ddl info file");
+
+  Int_t ddlno;
+  Int_t modno;
+  Int_t modulePerDDL;
+  Int_t moduleddl[6];
+
+  for(Int_t jddl = 0; jddl < 6; jddl++)
+    {
+      if (infileddl.eof()) break;
+      infileddl >> ddlno >> modulePerDDL;
+      moduleddl[jddl] = modulePerDDL;
+
+      if (modulePerDDL == 0) continue;
+      for (Int_t im = 0; im < modulePerDDL; im++)
+       {
+         infileddl >> modno;
+       }
+    }
+
+  infileddl.close();
+
+  // Set the minimum noise cut per module before clustering
+
+  fRecoParam = AliPMDReconstructor::GetRecoParam();
+
+  if(fRecoParam == 0x0)
+    {
+       AliFatal("No Reco Param found for PMD!!!");
+    }
+
 
   fRunLoader->GetEvent(ievt);
 
@@ -641,23 +934,20 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
   TBranch *branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize); 
   TBranch *branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize); 
 
-  const Int_t kDDL = AliDAQ::NumberOfDdls("PMD");
   const Int_t kRow = 48;
   const Int_t kCol = 96;
 
   Int_t idet = 0;
   Int_t iSMN = 0;
-  
-  for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
+
+  AliPMDRawStream pmdinput(rawReader);
+  Int_t indexDDL = -1;
+
+  while ((indexDDL = pmdinput.DdlData(&pmdddlcont)) >=0)
     {
-      if (indexDDL < 4)
-       {
-         iSMN = 6;
-       }
-      else if (indexDDL >= 4)
-       {
-         iSMN = 12;
-       }
+      
+      iSMN = moduleddl[indexDDL];
+
       Int_t ***precpvADC;
       precpvADC = new int **[iSMN];
       for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
@@ -676,11 +966,7 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
            }
        }
       ResetCellADC();
-      rawReader->Reset();
-      rawReader->Select("PMD", indexDDL, indexDDL);
 
-      AliPMDRawStream pmdinput(rawReader);
-      pmdinput.DdlData(indexDDL,&pmdddlcont);
     
       Int_t indexsmn = 0;
       Int_t ientries = pmdddlcont.GetEntries();
@@ -695,25 +981,48 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
          Int_t row = pmdddl->GetRow();
          Int_t col = pmdddl->GetColumn();
          Int_t sig = pmdddl->GetSignal();
-
+         if(row < 0 || row > 48 || col < 0 || col > 96)
+           {
+             AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
+                             row, col));
+             continue; 
+           }
          // Pedestal Subtraction
          Int_t   pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,row,col);
-         Int_t   pedrms1    = (Int_t) pedmeanrms%1000;
+         Int_t   pedrms1    = (Int_t) pedmeanrms%100;
          Float_t pedrms     = (Float_t)pedrms1/10.;
          Float_t pedmean    = (Float_t) (pedmeanrms - pedrms1)/1000.0;
-
          //printf("%f %f\n",pedmean, pedrms);
-
-         //Float_t sig1 = (Float_t) sig;
          Float_t sig1 = (Float_t) sig - (pedmean + 3.0*pedrms);
+
+         // Hot cell - set the cell adc = 0
+         Float_t hotflag = fCalibHot->GetHotChannel(det,smn,row,col);
+         if (hotflag == 1) sig1 = 0;
+
          // CALIBRATION
          Float_t gain = fCalibGain->GetGainFact(det,smn,row,col);
 
          //printf("sig = %d gain = %f\n",sig,gain);
          sig = (Int_t) (sig1*gain);
 
-
-         if (indexDDL < 4)
+         if (indexDDL == 0)
+           {
+             if (det != 0)
+               AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
+                             indexDDL, det));
+             if (iSMN == 6)
+               {
+                 indexsmn = smn;
+               }
+             else if (iSMN == 12)
+               {
+                 if (smn < 6)
+                   indexsmn = smn;
+                 else if (smn >= 18 && smn < 24)
+                   indexsmn = smn-12;
+               }
+           }
+         else if (indexDDL >= 1 && indexDDL < 4)
            {
              if (det != 0)
                AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
@@ -739,36 +1048,51 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
              if (det != 1)
                AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
                              indexDDL, det));
-             if (smn >= 6 && smn < 12)
-               {
-                 indexsmn = smn - 6;
-               }
-             else if (smn >= 12 && smn < 18)
+             if (smn >= 6 && smn < 18)
                {
                  indexsmn = smn - 6;
                }
            }         
+         
          precpvADC[indexsmn][row][col] = sig;
 
        }
       
-      pmdddlcont.Clear();
+      pmdddlcont.Delete();
 
       Int_t ismn = 0;
-      for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
+      for (indexsmn = 0; indexsmn < iSMN; indexsmn++)
        {
          ResetCellADC();
          for (Int_t irow = 0; irow < kRow; irow++)
            {
              for (Int_t icol = 0; icol < kCol; icol++)
                {
+                 fCellTrack[irow][icol] = -1;
+                 fCellPid[irow][icol]   = -1;
                  fCellADC[irow][icol] = 
                    (Double_t) precpvADC[indexsmn][irow][icol];
                } // row
            }     // col
 
-         
-         if (indexDDL < 4)
+
+         if (indexDDL == 0)
+           {
+             if (iSMN == 6)
+               {
+                 ismn = indexsmn;
+               }
+             else if (iSMN == 12)
+               {
+                 
+                 if (indexsmn < 6)
+                   ismn = indexsmn;
+                 else if (indexsmn >= 6 && indexsmn < 12)
+                   ismn = indexsmn + 12;
+               }
+             idet = 0;
+           }
+         else if (indexDDL >= 1 && indexDDL < 4)
            {
              ismn = indexsmn + indexDDL * 6;
              idet = 0;
@@ -787,26 +1111,28 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
            }
          else if (indexDDL == 5)
            {
-             if (indexsmn < 6)
-               {
-                 ismn = indexsmn + 6;
-               }
-             else if (indexsmn >= 6 && indexsmn < 12)
-               {
-                 ismn = indexsmn + 6;
-               }
+             ismn = indexsmn + 6;
              idet = 1;
            }
 
-         pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
+         Int_t imod = idet*24 + ismn;
+         fEcut = fRecoParam->GetNoiseCut(imod);       // default
+         // fEcut = fRecoParam->GetPbPbParam()->GetNoiseCut(imod);
+         // fEcut = fRecoParam->GetPPParam()->GetNoiseCut(imod);
+         // fEcut = fRecoParam->GetCosmicParam()->GetNoiseCut(imod);
+
+         pmdclust->SetEdepCut(fEcut);
+
+         pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,
+                           pmdisocell,pmdcont);
+
          Int_t nentries1 = pmdcont->GetEntries();
 
          AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
 
          for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
            {
-             AliPMDcluster *pmdcl = 
-               (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
+             pmdcl       = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
              idet        = pmdcl->GetDetector();
              ismn        = pmdcl->GetSMN();
              clusdata[0] = pmdcl->GetClusX();
@@ -823,13 +1149,36 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
                {
                  Int_t celldataX = pmdcl->GetClusCellX(ihit);
                  Int_t celldataY = pmdcl->GetClusCellY(ihit);
-                 AddRecHit(celldataX, celldataY);
+                 Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
+                 Int_t celldataPid = pmdcl->GetClusCellPid(ihit);
+                 Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
+                 AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
                }
              branch2->Fill();
              ResetRechit();
 
            }
-         pmdcont->Clear();
+         pmdcont->Delete();
+
+         // Added single isolated cell for offline gain calibration
+         nentries1 = pmdisocell->GetEntries();
+         AliDebug(1,Form("Total number of isolated single cell clusters = %d",nentries1));
+         
+         for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
+           {
+             pmdiso = (AliPMDisocell*)pmdisocell->UncheckedAt(ient1);
+             idet = pmdiso->GetDetector();
+             ismn = pmdiso->GetSmn();
+             clusdata[0] = (Float_t) pmdiso->GetRow();
+             clusdata[1] = (Float_t) pmdiso->GetCol();
+             clusdata[2] = pmdiso->GetADC();
+             clusdata[3] = 1.;
+             clusdata[4] = -99.;
+             clusdata[5] = -99.;
+      
+             AddRecPoint(idet,ismn,clusdata);
+           }
+         pmdisocell->Delete();
          
          branch1->Fill();
          ResetRecpoint();
@@ -844,7 +1193,8 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
       for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
       delete precpvADC;
     } // DDL Loop
-  
+
+
   ResetCellADC();
   
   fPMDLoader = fRunLoader->GetLoader("PMDLoader");  
@@ -853,7 +1203,7 @@ void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
   //   delete the pointers
   delete pmdclust;
   delete pmdcont;
-
+  delete pmdisocell;
 }
 // ------------------------------------------------------------------------- //
 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
@@ -872,13 +1222,15 @@ void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
   delete newrecpoint;
 }
 // ------------------------------------------------------------------------- //
-void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY)
+void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY,
+                                   Int_t celldataTr, Int_t celldataPid,
+                                   Float_t celldataAdc)
 {
   // Add associated cell hits to the Reconstructed points
   //
   TClonesArray &lrechits = *fRechits;
   AliPMDrechit *newrechit;
-  newrechit = new AliPMDrechit(celldataX, celldataY);
+  newrechit = new AliPMDrechit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
   new(lrechits[fNhit++]) AliPMDrechit(newrechit);
   delete newrechit;
 }
@@ -891,12 +1243,13 @@ void AliPMDClusterFinder::ResetCellADC()
     {
       for(Int_t icol = 0; icol < fgkCol; icol++)
        {
-         fCellADC[irow][icol] = 0.;
+         fCellTrack[irow][icol] = -1;
+         fCellPid[irow][icol]   = -1;
+         fCellADC[irow][icol]   = 0.;
        }
     }
 }
 // ------------------------------------------------------------------------- //
-
 void AliPMDClusterFinder::ResetRecpoint()
 {
   // Clear the list of reconstructed points
@@ -941,7 +1294,6 @@ void AliPMDClusterFinder::UnLoadClusters()
   fPMDLoader->UnloadRecPoints();
 }
 // ------------------------------------------------------------------------- //
-
 AliPMDCalibData* AliPMDClusterFinder::GetCalibGain() const
 {
   // The run number will be centralized in AliCDBManager,
@@ -958,9 +1310,7 @@ AliPMDCalibData* AliPMDClusterFinder::GetCalibGain() const
   
   return calibdata;
 }
-
 // ------------------------------------------------------------------------- //
-
 AliPMDPedestal* AliPMDClusterFinder::GetCalibPed() const
 {
   // The run number will be centralized in AliCDBManager,
@@ -976,3 +1326,20 @@ AliPMDPedestal* AliPMDClusterFinder::GetCalibPed() const
   
   return pedestal;
 }
+//--------------------------------------------------------------------//
+AliPMDHotData* AliPMDClusterFinder::GetCalibHot() const
+{
+  // The run number will be centralized in AliCDBManager,
+  // you don't need to set it here!
+  AliCDBEntry  *entry = AliCDBManager::Instance()->Get("PMD/Calib/Hot");
+  
+  if(!entry) AliFatal("HotData object retrieval failed!");
+  
+  AliPMDHotData *hot = 0;
+  if (entry) hot = (AliPMDHotData*) entry->GetObject();
+  
+  if (!hot)  AliFatal("No hot data from  database !");
+  
+  return hot;
+}
+