]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
from Alex Kalweit:
authormkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 14 Nov 2013 21:16:48 +0000 (21:16 +0000)
committermkrzewic <mkrzewic@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 14 Nov 2013 21:16:48 +0000 (21:16 +0000)
  patch for the mip momentum range and the pad equalization

TPC/Calib/AliTPCcalibGainMult.cxx
TPC/Calib/AliTPCcalibGainMult.h
TPC/Calib/AliTPCcalibTimeGain.cxx
TPC/Calib/AliTPCcalibTimeGain.h

index d566def32251339f65fdda64039d9ddd5cb17ad0..2d47564282e327a8270469301e028184ddfab8cb 100644 (file)
@@ -83,6 +83,9 @@ AliTPCcalibGainMult::AliTPCcalibGainMult()
    fCutRequireITSrefit(0),
    fCutMaxDcaXY(0),
    fCutMaxDcaZ(0),
+   fMinMomentumMIP(0),
+   fMaxMomentumMIP(0),
+   fAlephParameters(),
    fHistNTracks(0),
    fHistClusterShape(0),
    fHistQA(0),
@@ -114,6 +117,9 @@ AliTPCcalibGainMult::AliTPCcalibGainMult(const Text_t *name, const Text_t *title
    fCutRequireITSrefit(0),
    fCutMaxDcaXY(0),
    fCutMaxDcaZ(0),
+   fMinMomentumMIP(0),
+   fMaxMomentumMIP(0),
+   fAlephParameters(),
    fHistNTracks(0),
    fHistClusterShape(0),
    fHistQA(0),
@@ -145,6 +151,14 @@ AliTPCcalibGainMult::AliTPCcalibGainMult(const Text_t *name, const Text_t *title
   fCutRequireITSrefit = kFALSE;
   fCutMaxDcaXY = 3.5;
   fCutMaxDcaZ  = 25.;
+  // default values for MIP window selection
+  fMinMomentumMIP = 0.4;
+  fMaxMomentumMIP = 0.6;
+  fAlephParameters[0] = 0.07657; // the following parameters work for most of the periods and are therefore default
+  fAlephParameters[1] = 10.6654; // but they can be overwritten in the train setup of cpass0
+  fAlephParameters[2] = 2.51466e-14;
+  fAlephParameters[3] = 2.05379;
+  fAlephParameters[4] = 1.84288;
 
   //
   fHistNTracks = new TH1F("ntracks","Number of Tracks per Event; number of tracks per event; number of tracks",1001,-0.5,1000.5);
@@ -274,7 +288,7 @@ void AliTPCcalibGainMult::Process(AliESDEvent *event) {
   //const Double_t kMaxEta=0.8; // maximal eta fo the track to be accepted
   //const Double_t kMaxDCAR=10; // maximal DCA R of the track
   //const Double_t kMaxDCAZ=5;  // maximal DCA Z of the track
-  const Double_t kMIPPt=0.525; // MIP pt
+  //  const Double_t kMIPPt=0.525; // MIP pt
   
   if (!event) {
     Printf("ERROR: ESD not available");
@@ -319,6 +333,15 @@ void AliTPCcalibGainMult::Process(AliESDEvent *event) {
     Int_t ncls = track->GetTPCNcls();
     Int_t nCrossedRows = track->GetTPCCrossedRows();
 
+    // correction factor of dE/dx in MIP window
+    Float_t corrFactorMip = AliExternalTrackParam::BetheBlochAleph(meanP/0.13957, 
+                                                                  fAlephParameters[0], 
+                                                                  fAlephParameters[1], 
+                                                                  fAlephParameters[2], 
+                                                                  fAlephParameters[3],
+                                                                  fAlephParameters[4]);
+    if (TMath::Abs(corrFactorMip) < 1e-10) continue;
+    
     if (nCrossedRows < fCutCrossRows) continue;     
     // exclude tracks which do not look like primaries or are simply too short or on wrong sectors
     if (TMath::Abs(trackIn->Eta()) > fCutEtaWindow) continue;
@@ -393,36 +416,39 @@ void AliTPCcalibGainMult::Process(AliESDEvent *event) {
       fHistQA->Fill(meanP, signal, 3);
       fHistQA->Fill(meanP, mipSignalOroc, 4);
       //
-      // "dEdxRatioMax","dEdxRatioTot","padType","mult","driftlength", "1_pt"
-      Float_t meanMax = (1/3.)*(signalArrayMax[0] + signalArrayMax[1] + signalArrayMax[2]);
-      Float_t meanTot = (1/3.)*(signalArrayTot[0] + signalArrayTot[1] + signalArrayTot[2]); 
+      // normalize pad regions to their common mean
+      //
+      Float_t meanMax = (63./159)*signalArrayMax[0] + (64./159)*signalArrayMax[1] + (32./159)*signalArrayMax[2];
+      Float_t meanTot = (63./159)*signalArrayTot[0] + (64./159)*signalArrayTot[1] + (32./159)*signalArrayTot[2]; 
       if (meanMax < 1e-5 || meanTot < 1e-5) continue;
+      //
       for(Int_t ipad = 0; ipad < 4; ipad ++) {
+       // "dEdxRatioMax","dEdxRatioTot","padType","mult","driftlength", "1_pt"
        Double_t vecPadEqual[5] = {signalArrayMax[ipad]/meanMax, signalArrayTot[ipad]/meanTot, ipad, nContributors, meanDrift};
-       if ( TMath::Abs(meanP-kMIPPt)<0.025 ) fHistPadEqual->Fill(vecPadEqual);
+       if (fMinMomentumMIP > meanP && meanP < fMaxMomentumMIP) fHistPadEqual->Fill(vecPadEqual);
       }
       //
-      //      if (meanP > 0.4 && meanP < 0.55) {
-      if ( TMath::Abs(meanP-kMIPPt)<0.025 ) {
-       Double_t vecMult[6] = {seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,row0,row1),
-                              seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,row0,row1),
+      //
+      if (meanP < fMaxMomentumMIP && meanP > fMinMomentumMIP) {
+       Double_t vecMult[6] = {seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,1,row0,row1)/corrFactorMip,
+                              seed->CookdEdxAnalytical(fLowerTrunc,fUpperTrunc,0,row0,row1)/corrFactorMip,
                               meanDrift,
                               3,
                               nContributors,
                               ncls};
        //
        fHistGainMult->Fill(vecMult);
-       vecMult[0]=mipSignalShort; vecMult[1]=mipSignalShort; vecMult[3]=0;
+       vecMult[0]=mipSignalShort/corrFactorMip; vecMult[1]=mipSignalShort/corrFactorMip; vecMult[3]=0;
        fHistGainMult->Fill(vecMult);
-       vecMult[0]=mipSignalMed; vecMult[1]=mipSignalMed; vecMult[3]=1;
+       vecMult[0]=mipSignalMed/corrFactorMip; vecMult[1]=mipSignalMed/corrFactorMip; vecMult[3]=1;
        fHistGainMult->Fill(vecMult);
-       vecMult[0]=mipSignalLong; vecMult[1]=mipSignalLong; vecMult[3]=2;
+       vecMult[0]=mipSignalLong/corrFactorMip; vecMult[1]=mipSignalLong/corrFactorMip; vecMult[3]=2;
        fHistGainMult->Fill(vecMult);
        //
       }
       //
       //
-      if ( TMath::Abs(meanP-kMIPPt)>0.025 ) continue;  // only MIP pions
+      if (fMinMomentumMIP < meanP || meanP > fMaxMomentumMIP) continue;  // only MIP pions
       //if (meanP > 0.5 || meanP < 0.4) continue; // only MIP pions
       //
       // for each track, we look at the three different pad regions, split it into tracklets, check if the sector does not change and fill the histogram
@@ -450,7 +476,7 @@ void AliTPCcalibGainMult::Process(AliESDEvent *event) {
       //
       //                        MIP, sect,  pad,   run
       //
-      Double_t vecMip[5] = {mipSignalShort, mipSignalMed, mipSignalLong, signal, mipSignalOroc};
+      Double_t vecMip[5] = {mipSignalShort/corrFactorMip, mipSignalMed/corrFactorMip, mipSignalLong/corrFactorMip, signal/corrFactorMip, mipSignalOroc/corrFactorMip};
       //
       for(Int_t ipad = 0; ipad < 3; ipad++) {
        // AK. -  run Number To be removed - not needed 
index fc30797004113276d4ea5723efaa4a3062c95764..4e0eefc73b61278a02189b054e23c2edbc992a01 100644 (file)
@@ -63,6 +63,10 @@ public:
   void SetCutMaxDcaXY(Float_t maxXY){fCutMaxDcaXY = maxXY;}; 
   void SetCutMaxDcaZ(Float_t maxZ){fCutMaxDcaZ = maxZ;}; 
   //
+  void SetMinMomentumMIP(Float_t minMom = 0.4){fMinMomentumMIP = minMom;};
+  void SetMaxMomentumMIP(Float_t maxMom = 0.6){fMaxMomentumMIP = maxMom;};
+  void SetAlephParameters(Float_t * parameters){for(Int_t j=0;j<5;j++) fAlephParameters[j] = parameters[j];};
+  //
   //
   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
   void     Process(AliTPCseed *track){return AliTPCcalibBase::Process(track);}
@@ -96,6 +100,12 @@ private:
   Float_t fCutMaxDcaXY;                 // max dca_xy (only TPConly resolution is guaranteed!)
   Float_t fCutMaxDcaZ;                  // max dca_z  (dangerous if vDrift is not calibrated)
   //
+  // definition of MIP window
+  //
+  Float_t fMinMomentumMIP;              // minimum momentum of MIP region, e.g. 400 MeV
+  Float_t fMaxMomentumMIP;              // maximum momentum of MIP region, e.g. 600 MeV
+  Float_t fAlephParameters[5];          // parameters for equalization in MIP window, parameter set should be =1 at MIP
+  //
   // histograms
   //
   TH1F  *fHistNTracks;            //  histogram showing number of ESD tracks per event
@@ -117,7 +127,7 @@ private:
   AliTPCcalibGainMult(const AliTPCcalibGainMult&); 
   AliTPCcalibGainMult& operator=(const AliTPCcalibGainMult&); 
 
-  ClassDef(AliTPCcalibGainMult, 2); 
+  ClassDef(AliTPCcalibGainMult, 3); 
 };
 
 #endif
index e8bcc6bb197e222216f46d6e72f4cf4995e30ba6..a06b99adfaecf5b556a8b41d26251f8b639dee4d 100644 (file)
@@ -194,6 +194,9 @@ AliTPCcalibTimeGain::AliTPCcalibTimeGain()
    fCutRequireITSrefit(0),
    fCutMaxDcaXY(0),
    fCutMaxDcaZ(0),
+   fMinMomentumMIP(0),
+   fMaxMomentumMIP(0),
+   fAlephParameters(),
    fUseMax(0),
    fLowerTrunc(0),
    fUpperTrunc(0),
@@ -223,6 +226,9 @@ AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title
    fCutRequireITSrefit(0),
    fCutMaxDcaXY(0),
    fCutMaxDcaZ(0),
+   fMinMomentumMIP(0),
+   fMaxMomentumMIP(0),
+   fAlephParameters(),
    fUseMax(0),
    fLowerTrunc(0),
    fUpperTrunc(0),
@@ -265,6 +271,15 @@ AliTPCcalibTimeGain::AliTPCcalibTimeGain(const Text_t *name, const Text_t *title
   fCutMaxDcaXY = 3.5;
   fCutMaxDcaZ  = 25.;
 
+  // default values for MIP window selection
+  fMinMomentumMIP = 0.4;
+  fMaxMomentumMIP = 0.6;
+  fAlephParameters[0] = 0.07657; // the following parameters work for most of the periods and are therefore default
+  fAlephParameters[1] = 10.6654; // but they can be overwritten in the train setup of cpass0
+  fAlephParameters[2] = 2.51466e-14;
+  fAlephParameters[3] = 2.05379;
+  fAlephParameters[4] = 1.84288;
+
   // default values for dE/dx
   fMIP = 50.;
   fUseMax = kTRUE;
@@ -442,19 +457,25 @@ void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) {
 
     if (seed) {
       Int_t particleCase = 0;
-      if (meanP < 0.55 && meanP > 0.5)  particleCase = 2; // MIP pions
+      if (meanP < fMaxMomentumMIP && meanP > fMinMomentumMIP)  particleCase = 2; // MIP pions
       if (meanP < 0.57 && meanP > 0.56) particleCase = 3; // protons 1
       if (meanP < 0.66 && meanP > 0.65) particleCase = 4; // protons 2
       //
       if (fLowMemoryConsumption && particleCase == 0) continue;
       //
       Double_t tpcSignal = GetTPCdEdx(seed);
-      /*
-      if (particleCalse == 0) {
-       Float_t corrFactor = AliExternalTrackParam::BetheBlochAleph(x/0.13957);
-         tpcSignal /= corrFactor; 
-      }
-      */
+      //
+      // flattens signal in MIP window
+      //
+      if (particleCase == 0) {
+       Float_t corrFactor = AliExternalTrackParam::BetheBlochAleph(meanP/0.13957, 
+                                                                   fAlephParameters[0], 
+                                                                   fAlephParameters[1], 
+                                                                   fAlephParameters[2], 
+                                                                   fAlephParameters[3],
+                                                                   fAlephParameters[4]);
+       tpcSignal /= corrFactor; 
+      }        
       fHistDeDxTotal->Fill(meanP, tpcSignal);
       //
       //dE/dx, time, type (1-muon cosmic,2-pion beam data, 3&4 protons), momenta, runNumner, eta
@@ -500,7 +521,7 @@ void AliTPCcalibTimeGain::ProcessBeamEvent(AliESDEvent *event) {
       }    
       if (seed) { 
        if (fLowMemoryConsumption) {
-         if (meanP > 0.5 || meanP < 0.4) continue;
+         if (meanP > fMaxMomentumMIP || meanP < fMinMomentumMIP) continue;
          meanP = 0.45; // set all momenta to one in order to save memory
       }
        Double_t tpcSignal = GetTPCdEdx(seed);
index efe0de40abfade38be4a38cfbb76f9490eb44835..385417b51ad61e9f172b63e70bf065d84f199bf0 100644 (file)
@@ -66,7 +66,9 @@ public:
   void SetCutMaxDcaXY(Float_t maxXY){fCutMaxDcaXY = maxXY;}; 
   void SetCutMaxDcaZ(Float_t maxZ){fCutMaxDcaZ = maxZ;}; 
   //
-
+  void SetMinMomentumMIP(Float_t minMom = 0.4){fMinMomentumMIP = minMom;};
+  void SetMaxMomentumMIP(Float_t maxMom = 0.6){fMaxMomentumMIP = maxMom;};
+  void SetAlephParameters(Float_t * parameters){for(Int_t j=0;j<5;j++) fAlephParameters[j] = parameters[j];};
 
   static void SetMergeEntriesCut(Double_t entriesCut){fgMergeEntriesCut = entriesCut;}
 
@@ -93,6 +95,12 @@ private:
   Float_t fCutMaxDcaXY;                 // max dca_xy (only TPConly resolution is guaranteed!)
   Float_t fCutMaxDcaZ;                  // max dca_z  (dangerous if vDrift is not calibrated)
   //
+  // definition of MIP window
+  //
+  Float_t fMinMomentumMIP;              // minimum momentum of MIP region, e.g. 400 MeV
+  Float_t fMaxMomentumMIP;              // maximum momentum of MIP region, e.g. 600 MeV
+  Float_t fAlephParameters[5];          // parameters for equalization in MIP window, parameter set should be =1 at MIP
+  //
   //
   Bool_t  fUseMax;                      // true: use max charge for dE/dx calculation, false: use total charge for dE/dx calculation
   Float_t fLowerTrunc;                  // lower truncation of dE/dx ; at most 5%
@@ -110,7 +118,7 @@ private:
   void     Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
   void     Process(AliTPCseed *track){return AliTPCcalibBase::Process(track);}
 
-  ClassDef(AliTPCcalibTimeGain, 1); 
+  ClassDef(AliTPCcalibTimeGain, 2);
 };
 
 #endif