]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Using gain factor from CDB (M.Ivanov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 3 Apr 2006 17:23:07 +0000 (17:23 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 3 Apr 2006 17:23:07 +0000 (17:23 +0000)
TPC/AliTPCDigitizer.cxx
TPC/AliTPCclustererMI.cxx
TPC/AliTPCclustererMI.h

index 2d61a4e135123641e9fe6f4cb1fb655d3d6057fc..fb7e1bb206e181edb1012e45093fbaa921460783 100644 (file)
 #include "AliSimDigits.h"
 #include "AliLog.h"
 
+#include "AliTPCcalibDB.h"
+#include "AliTPCCalPad.h"
+#include "AliTPCCalROC.h"
+
 ClassImp(AliTPCDigitizer)
 
 //___________________________________________
@@ -177,6 +181,7 @@ void AliTPCDigitizer::ExecFast(Option_t* option)
   param->SetZeroSup(2);
 
   Int_t zerosup = param->GetZeroSup(); 
+  AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor(); 
   //
   //Loop over segments of the TPC
     
@@ -188,7 +193,7 @@ void AliTPCDigitizer::ExecFast(Option_t* option)
       cerr<<"AliTPC warning: invalid segment ID ! "<<segmentID<<endl;
       continue;
      }
-
+    AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec);  // pad gains per given sector
     digrow->SetID(segmentID);
 
     Int_t nrows = 0;
@@ -260,6 +265,11 @@ void AliTPCDigitizer::ExecFast(Option_t* option)
            ptr[i]++;
          }
         q/=16.;  //conversion factor
+       Float_t gain = gainROC->GetValue(row,elem/nrows);  // get gain for given - pad-row pad
+       if (gain<0.5){
+         printf("problem\n");
+       }
+       q*= gain;
         //       Float_t noise  = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());  
         Float_t noise  = pTPC->GetNoise();
         q+=noise;
@@ -371,6 +381,7 @@ void AliTPCDigitizer::ExecSave(Option_t* option)
   Int_t zerosup = param->GetZeroSup();
   //Loop over segments of the TPC
     
+  AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
   for (Int_t n=0; n<nentries; n++) {
     rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
     gime = rl->GetLoader("TPCLoader");
@@ -378,7 +389,8 @@ void AliTPCDigitizer::ExecSave(Option_t* option)
 
     digarr[0]->ExpandBuffer();
     digarr[0]->ExpandTrackBuffer();
-           
+
+
     for (Int_t i=1;i<nInputs; i++){ 
 //      fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());      
       rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
@@ -397,6 +409,7 @@ void AliTPCDigitizer::ExecSave(Option_t* option)
       continue;
     }
 
+    AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec);  // pad gains per given sector
     digrow->SetID(digarr[0]->GetID());
 
     Int_t nrows = digarr[0]->GetNRows();
@@ -434,6 +447,8 @@ void AliTPCDigitizer::ExecSave(Option_t* option)
         }
        q/=16.;  //conversion factor
        //       Float_t noise  = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());  
+       Float_t gain = gainROC->GetValue(row,col);
+       q*= gain;
        Float_t noise  = pTPC->GetNoise();
        q+=noise;
         q=TMath::Nint(q);
index 07ba2b2039c700988dd0dd0860f8584a5a65f92d..8d222f2b629e281b0b78d954fdc679dd0a6ad8ff 100644 (file)
 #include "Riostream.h"
 #include <TTree.h>
 
+#include "AliTPCcalibDB.h"
+#include "AliTPCCalPad.h"
+#include "AliTPCCalROC.h"
+
+
 ClassImp(AliTPCclustererMI)
 
 
@@ -99,7 +104,7 @@ Float_t  AliTPCclustererMI::GetSigmaZ2(Int_t iz){
   return res;
 }
 
-void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Int_t *bins, UInt_t /*m*/,
+void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
 AliTPCclusterMI &c) 
 {
   Int_t i0=k/max;  //central pad
@@ -107,8 +112,8 @@ AliTPCclusterMI &c)
 
   // set pointers to data
   //Int_t dummy[5] ={0,0,0,0,0};
-  Int_t * matrix[5]; //5x5 matrix with digits  - indexing i = 0 ..4  j = -2..2
-  Int_t * resmatrix[5];
+  Float_t * matrix[5]; //5x5 matrix with digits  - indexing i = 0 ..4  j = -2..2
+  Float_t * resmatrix[5];
   for (Int_t di=-2;di<=2;di++){
     matrix[di+2]  =  &bins[k+di*max];
     resmatrix[di+2]  =  &fResBins[k+di*max];
@@ -205,7 +210,7 @@ AliTPCclusterMI &c)
     //remove cluster data from data
     for (Int_t di=-2;di<=2;di++)
       for (Int_t dj=-2;dj<=2;dj++){
-       resmatrix[di+2][dj] -= Int_t(vmatrix[di+2][dj+2]);
+       resmatrix[di+2][dj] -= vmatrix[di+2][dj+2];
        if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
       }
     resmatrix[2][0] =0;
@@ -215,8 +220,8 @@ AliTPCclusterMI &c)
   //unfolding when neccessary  
   //
   
-  Int_t * matrix2[7]; //7x7 matrix with digits  - indexing i = 0 ..6  j = -3..3
-  Int_t dummy[7]={0,0,0,0,0,0};
+  Float_t * matrix2[7]; //7x7 matrix with digits  - indexing i = 0 ..6  j = -3..3
+  Float_t dummy[7]={0,0,0,0,0,0};
   for (Int_t di=-3;di<=3;di++){
     matrix2[di+3] =  &bins[k+di*max];
     if ((k+di*max)<3)  matrix2[di+3] = &dummy[3];
@@ -248,7 +253,7 @@ AliTPCclusterMI &c)
 
 
 
-void AliTPCclustererMI::UnfoldCluster(Int_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj, 
+void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj, 
                                      Float_t & sumu, Float_t & overlap )
 {
   //
@@ -461,6 +466,8 @@ void AliTPCclustererMI::Digits2Clusters()
     return;
   }
 
+  AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
+
   AliSimDigits digarr, *dummy=&digarr;
   fRowDig = dummy;
   fInput->GetBranch("Segment")->SetAddress(&dummy);
@@ -469,7 +476,7 @@ void AliTPCclustererMI::Digits2Clusters()
   fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3   after
     
   Int_t nclusters  = 0;
-  
+
   for (Int_t n=0; n<nentries; n++) {
     fInput->GetEvent(n);
     Int_t row;
@@ -477,7 +484,9 @@ void AliTPCclustererMI::Digits2Clusters()
       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
       continue;
     }
-        
+    AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector);  // pad gains per given sector
+    
+    //
     AliTPCClustersRow *clrow= new AliTPCClustersRow();
     fRowCl = clrow;
     clrow->SetClass("AliTPCclusterMI");
@@ -504,16 +513,18 @@ void AliTPCclustererMI::Digits2Clusters()
     
     
     fMaxBin=fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
-    fBins    =new Int_t[fMaxBin];
-    fResBins =new Int_t[fMaxBin];  //fBins with residuals after 1 finder loop 
-    memset(fBins,0,sizeof(Int_t)*fMaxBin);
+    fBins    =new Float_t[fMaxBin];
+    fResBins =new Float_t[fMaxBin];  //fBins with residuals after 1 finder loop 
+    memset(fBins,0,sizeof(Float_t)*fMaxBin);
+    memset(fResBins,0,sizeof(Float_t)*fMaxBin);
     
     if (digarr.First()) //MI change
       do {
-       Short_t dig=digarr.CurrentDigit();
+       Float_t dig=digarr.CurrentDigit();
        if (dig<=fParam->GetZeroSup()) continue;
        Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
-       fBins[i*fMaxTime+j]=dig;
+       Float_t gain = gainROC->GetValue(row,digarr.CurrentRow()/fParam->GetMaxTBin());
+       fBins[i*fMaxTime+j]=dig/gain;
       } while (digarr.Next());
     digarr.ExpandTrackBuffer();
 
@@ -555,8 +566,8 @@ void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
   Int_t zeroSup = fParam->GetZeroSup();
 
   fBins = NULL;
-  Int_t** splitRows = new Int_t* [kNS*2];
-  Int_t** splitRowsRes = new Int_t* [kNS*2];
+  Float_t** splitRows = new Float_t* [kNS*2];
+  Float_t** splitRowsRes = new Float_t* [kNS*2];
   for (Int_t iSector = 0; iSector < kNS*2; iSector++)
     splitRows[iSector] = NULL;
   Int_t iSplitRow = -1;
@@ -615,9 +626,10 @@ void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
     
       fMaxBin = fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
       if ((iSplitRow < 0) || !splitRows[fSector + kNS*iSplitRow]) {
-       fBins    = new Int_t[fMaxBin];
-       fResBins = new Int_t[fMaxBin];  //fBins with residuals after 1 finder loop 
-       memset(fBins, 0, sizeof(Int_t)*fMaxBin);
+       fBins    = new Float_t[fMaxBin];
+       fResBins = new Float_t[fMaxBin];  //fBins with residuals after 1 finder loop 
+       memset(fBins, 0, sizeof(Float_t)*fMaxBin);
+       memset(fResBins, 0, sizeof(Float_t)*fMaxBin);
       } else {
        fBins    = splitRows[fSector + kNS*iSplitRow];
        fResBins = splitRowsRes[fSector + kNS*iSplitRow];
@@ -688,7 +700,7 @@ void AliTPCclustererMI::FindClusters()
       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);  
       if (gDebug>4) printf("\n%f\n",amp0);
     }  
-    fBins[i+2*fMaxTime] = Int_t(amp0);    
+    fBins[i+2*fMaxTime] = amp0;    
     amp0 = 0;
     amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
     if (amp1>0){
@@ -698,7 +710,7 @@ void AliTPCclustererMI::FindClusters()
       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);          
       if (gDebug>4) printf("\n%f\n",amp0);
     }        
-    fBins[(fMaxPad+3)*fMaxTime+i] = Int_t(amp0);           
+    fBins[(fMaxPad+3)*fMaxTime+i] = amp0;           
   }
 
 //  memcpy(fResBins,fBins, fMaxBin*2);
@@ -707,7 +719,7 @@ void AliTPCclustererMI::FindClusters()
   fNcluster=0;
   //first loop - for "gold cluster" 
   fLoop=1;
-  Int_t *b=&fBins[-1]+2*fMaxTime;
+  Float_t *b=&fBins[-1]+2*fMaxTime;
   Int_t crtime = Int_t((fParam->GetZLength()-AliTPCReconstructor::GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
 
   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
index 70b5e8d5ba3bd64fcf05e0c5b5c7ab1da512733f..648579198f0259b21b142d92e8ec7e914b7bfa2d 100644 (file)
@@ -33,23 +33,23 @@ public:
   virtual void SetInput(TTree * tree);  // set input tree with digits    
   virtual void SetOutput(TTree * tree); //set output tree with 
 private:
-  Bool_t IsMaximum(Int_t k, Int_t max, const Int_t *bins); 
-  void MakeCluster2(Int_t k,Int_t max,Int_t *bins,UInt_t m,
+  Bool_t IsMaximum(Float_t k, Int_t max, const Float_t *bins); 
+  void MakeCluster2(Int_t k,Int_t max,Float_t *bins,UInt_t m,
    AliTPCclusterMI &c);  
-  void MakeCluster(Int_t k,Int_t max,Int_t *bins,UInt_t m,
+  void MakeCluster(Int_t k,Int_t max,Float_t *bins,UInt_t m,
    AliTPCclusterMI &c); 
   Float_t  GetSigmaY2(Int_t iz);
   Float_t  GetSigmaZ2(Int_t iz);
   Float_t  FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz);
   void AddCluster(AliTPCclusterMI &c);  // add the cluster to the array
-  void UnfoldCluster(Int_t * matrix[7], Float_t recmatrix[5][5], 
+  void UnfoldCluster(Float_t * matrix[7], Float_t recmatrix[5][5], 
                     Float_t & meani, Float_t & meanj, Float_t & sum, Float_t &overlap );
   void FindClusters();
 
 
 
-  Int_t * fBins;       //!digits array
-  Int_t * fResBins;    //!digits array with res. after 1 finder
+  Float_t * fBins;       //!digits array
+  Float_t * fResBins;    //!digits array with res. after 1 finder
   Int_t fLoop;         //loop - cf in 2 loops
   Int_t fMaxBin;
   Int_t fMaxTime;
@@ -70,7 +70,7 @@ private:
   ClassDef(AliTPCclustererMI,1)  // Time Projection Chamber digits
 };
 
-inline Bool_t AliTPCclustererMI::IsMaximum(Int_t q,Int_t max,const Int_t *bins){
+inline Bool_t AliTPCclustererMI::IsMaximum(Float_t q,Int_t max,const Float_t *bins){
   //is this a local maximum ?
   if (bins[-max] >= q) return kFALSE;
   if (bins[-1  ] >= q) return kFALSE;