]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGLF/RESONANCES/AliRsnCutTrackQuality.cxx
Completed changes needed because of previous commit
[u/mrichter/AliRoot.git] / PWGLF / RESONANCES / AliRsnCutTrackQuality.cxx
index 0ffb024992e304f042dd74c86995aa78611e9e13..ad09b185cf1a4b33df7a5887ab40669166f18b37 100644 (file)
@@ -35,9 +35,12 @@ AliRsnCutTrackQuality::AliRsnCutTrackQuality(const char *name) :
    fFlagsOn(0x0),
    fFlagsOff(0x0),
    fRejectKinkDaughters(kTRUE),
-   fDCARfixed(kTRUE),
+   fDCARmaxfixed(kTRUE),
+   fDCARminfixed(kTRUE),
    fDCARptFormula(""),
+   fDCARptFormulaMin(""),
    fDCARmax(1E20),
+   fDCARmin(0),
    fDCAZfixed(kTRUE),
    fDCAZptFormula(""),
    fDCAZmax(1E20),
@@ -47,14 +50,20 @@ AliRsnCutTrackQuality::AliRsnCutTrackQuality(const char *name) :
    fTPCminNClusters(0),
    fTPCmaxChi2(1E20),
    fCutMaxChi2TPCConstrainedVsGlobal(1E20),
+   fTrackMaxChi2(1E20),
+   fIsUseCrossedRowsCut(kFALSE),
+   fTPCminNCrossedRows(0),
+   fTPCminCrossedRowsOverFindableCls(0),
+   fIsUseLengthActiveVolumeTPCCut(kFALSE),
+   fCutMinLengthActiveVolumeTPC(0),
    fAODTestFilterBit(-1),
+   fCheckOnlyFilterBit(kTRUE),
    fESDtrackCuts(0x0)
 {
 //
 // Default constructor.
 // Initializes all cuts in such a way that all of them are disabled.
 //
-
    SetPtRange(0.0, 1E20);
    SetEtaRange(-1E20, 1E20);
 }
@@ -65,9 +74,12 @@ AliRsnCutTrackQuality::AliRsnCutTrackQuality(const AliRsnCutTrackQuality &copy)
    fFlagsOn(copy.fFlagsOn),
    fFlagsOff(copy.fFlagsOff),
    fRejectKinkDaughters(copy.fRejectKinkDaughters),
-   fDCARfixed(copy.fDCARfixed),
+   fDCARmaxfixed(copy.fDCARmaxfixed),
+   fDCARminfixed(copy.fDCARminfixed),
    fDCARptFormula(copy.fDCARptFormula),
+   fDCARptFormulaMin(copy.fDCARptFormulaMin),
    fDCARmax(copy.fDCARmax),
+   fDCARmin(copy.fDCARmin),
    fDCAZfixed(copy.fDCAZfixed),
    fDCAZptFormula(copy.fDCAZptFormula),
    fDCAZmax(copy.fDCAZmax),
@@ -77,7 +89,14 @@ AliRsnCutTrackQuality::AliRsnCutTrackQuality(const AliRsnCutTrackQuality &copy)
    fTPCminNClusters(copy.fTPCminNClusters),
    fTPCmaxChi2(copy.fTPCmaxChi2),
    fCutMaxChi2TPCConstrainedVsGlobal(copy.fCutMaxChi2TPCConstrainedVsGlobal),
+   fTrackMaxChi2(copy.fTrackMaxChi2),
+   fIsUseCrossedRowsCut(copy.fIsUseCrossedRowsCut),
+   fTPCminNCrossedRows(copy.fTPCminNCrossedRows),
+   fTPCminCrossedRowsOverFindableCls(copy.fTPCminCrossedRowsOverFindableCls),
+   fIsUseLengthActiveVolumeTPCCut(copy.fIsUseLengthActiveVolumeTPCCut),
+   fCutMinLengthActiveVolumeTPC(copy.fCutMinLengthActiveVolumeTPC),
    fAODTestFilterBit(copy.fAODTestFilterBit),
+   fCheckOnlyFilterBit(copy.fCheckOnlyFilterBit),
    fESDtrackCuts(copy.fESDtrackCuts)
 {
 //
@@ -103,9 +122,12 @@ AliRsnCutTrackQuality &AliRsnCutTrackQuality::operator=(const AliRsnCutTrackQual
    fFlagsOn = copy.fFlagsOn;
    fFlagsOff = copy.fFlagsOff;
    fRejectKinkDaughters = copy.fRejectKinkDaughters;
-   fDCARfixed = copy.fDCARfixed;
+   fDCARmaxfixed = copy.fDCARmaxfixed;
+   fDCARminfixed = copy.fDCARminfixed;
    fDCARptFormula = copy.fDCARptFormula;
+   fDCARptFormulaMin = copy.fDCARptFormulaMin;
    fDCARmax = copy.fDCARmax;
+   fDCARmin = copy.fDCARmin;
    fDCAZfixed = copy.fDCAZfixed;
    fDCAZptFormula = copy.fDCAZptFormula;
    fDCAZmax = copy.fDCAZmax;
@@ -114,7 +136,16 @@ AliRsnCutTrackQuality &AliRsnCutTrackQuality::operator=(const AliRsnCutTrackQual
    fITSmaxChi2 = copy.fITSmaxChi2;
    fTPCminNClusters = copy.fTPCminNClusters;
    fTPCmaxChi2 = copy.fTPCmaxChi2;
+   fCutMaxChi2TPCConstrainedVsGlobal = copy.fCutMaxChi2TPCConstrainedVsGlobal;
+   fTrackMaxChi2 = copy.fTrackMaxChi2;
+   fIsUseCrossedRowsCut=copy.fIsUseCrossedRowsCut;
+   fTPCminNCrossedRows = copy.fTPCminNCrossedRows;
+   fTPCminCrossedRowsOverFindableCls = copy.fTPCminCrossedRowsOverFindableCls;
+   fIsUseLengthActiveVolumeTPCCut=copy.fIsUseLengthActiveVolumeTPCCut;
+   fCutMinLengthActiveVolumeTPC = copy.fCutMinLengthActiveVolumeTPC;
+   
    fAODTestFilterBit = copy.fAODTestFilterBit;
+   fCheckOnlyFilterBit = copy.fCheckOnlyFilterBit;
    fESDtrackCuts = copy.fESDtrackCuts;
    SetPtRange(copy.fPt[0], copy.fPt[1]);
    SetEtaRange(copy.fEta[0], copy.fEta[1]);
@@ -132,9 +163,12 @@ void AliRsnCutTrackQuality::DisableAll()
    fFlagsOn = 0x0;
    fFlagsOff = 0x0;
    fRejectKinkDaughters = kFALSE;
-   fDCARfixed = kTRUE;
+   fDCARmaxfixed = kTRUE;
+   fDCARminfixed = kTRUE;
    fDCARptFormula = "";
+   fDCARptFormulaMin = "";
    fDCARmax = 1E20;
+   fDCARmin = 0;
    fDCAZfixed = kTRUE;
    fDCAZptFormula = "";
    fDCAZmax = 1E20;
@@ -144,6 +178,14 @@ void AliRsnCutTrackQuality::DisableAll()
    fTPCminNClusters = 0;
    fTPCmaxChi2 = 1E20;
    fAODTestFilterBit = -1;
+   fCutMaxChi2TPCConstrainedVsGlobal = 1E20;
+   fTrackMaxChi2 = 1E20;
+   fIsUseCrossedRowsCut = 0;
+   fTPCminNCrossedRows = 0;
+   fTPCminCrossedRowsOverFindableCls = 0;
+   fIsUseLengthActiveVolumeTPCCut = 0;
+   fCutMinLengthActiveVolumeTPC = 0.0;
    if (fESDtrackCuts) {
       const char *cutsName = fESDtrackCuts->GetName();
       const char *cutsTitle = fESDtrackCuts->GetTitle();
@@ -154,6 +196,26 @@ void AliRsnCutTrackQuality::DisableAll()
    SetEtaRange(-1E20, 1E20);
 }
 
+//_________________________________________________________________________________________________
+void AliRsnCutTrackQuality::SetPtRange(Double_t a, Double_t b)
+{
+  //Set Pt range cut
+  fPt[0] = TMath::Min(a, b); 
+  fPt[1] = TMath::Max(a, b);
+  if (fESDtrackCuts) fESDtrackCuts->SetPtRange(fPt[0], fPt[1]);
+  return;
+}
+
+//_________________________________________________________________________________________________
+void AliRsnCutTrackQuality::SetEtaRange(Double_t a, Double_t b)   
+{
+  //Set Pt range cut
+  fEta[0] = TMath::Min(a, b);
+  fEta[1] = TMath::Max(a, b);
+  if (fESDtrackCuts) fESDtrackCuts->SetEtaRange(fEta[0], fEta[1]);
+  return;
+}
+
 //_________________________________________________________________________________________________
 Bool_t AliRsnCutTrackQuality::IsSelected(TObject *object)
 {
@@ -218,18 +280,23 @@ Bool_t AliRsnCutTrackQuality::CheckESD(AliESDtrack *track)
 // This is done using the default track checker for ESD.
 // It is declared static, not to recreate it every time.
 //
-
-   static AliESDtrackCuts cuts;
+  //static AliESDtrackCuts cuts;
+  AliESDtrackCuts cuts;
 
    // general acceptance/pt cuts
    cuts.SetPtRange(fPt[0], fPt[1]);
    cuts.SetEtaRange(fEta[0], fEta[1]);
 
    // transverse DCA cuts
-   if (fDCARfixed)
+   if (fDCARmaxfixed)
       cuts.SetMaxDCAToVertexXY(fDCARmax);
    else
       cuts.SetMaxDCAToVertexXYPtDep(fDCARptFormula.Data());
+      
+   if (fDCARminfixed)
+      cuts.SetMinDCAToVertexXY(fDCARmin);
+   else
+      cuts.SetMinDCAToVertexXYPtDep(fDCARptFormulaMin.Data());
 
    // longitudinal DCA cuts
    if (fDCAZfixed)
@@ -242,17 +309,29 @@ Bool_t AliRsnCutTrackQuality::CheckESD(AliESDtrack *track)
    cuts.SetRequireSigmaToVertex(kFALSE);
 
    // TPC related cuts for TPC+ITS tracks
-   cuts.SetMinNClustersTPC(fTPCminNClusters);
+   if (fIsUseCrossedRowsCut) {
+     cuts.SetMinNCrossedRowsTPC(fTPCminNCrossedRows);
+     cuts.SetMinRatioCrossedRowsOverFindableClustersTPC(fTPCminCrossedRowsOverFindableCls);
+   } else {
+     cuts.SetMinNClustersTPC(fTPCminNClusters);
+   }
    cuts.SetMaxChi2PerClusterTPC(fTPCmaxChi2);
    cuts.SetAcceptKinkDaughters(!fRejectKinkDaughters);
    cuts.SetMaxChi2TPCConstrainedGlobal(fCutMaxChi2TPCConstrainedVsGlobal);
 
+   if (fIsUseLengthActiveVolumeTPCCut)
+     cuts.SetMinLengthActiveVolumeTPC(fCutMinLengthActiveVolumeTPC);
+
    // ITS related cuts for TPC+ITS tracks
    if (fSPDminNClusters > 0)
       cuts.SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny);
    cuts.SetMaxChi2PerClusterITS(fITSmaxChi2);
 
    // now that all is initialized, do the check
+   if (!track) {
+     AliError("Invalid track object. Rejected.");
+     return kFALSE;
+   }
    return cuts.IsSelected(track);
 }
 
@@ -274,7 +353,7 @@ Bool_t AliRsnCutTrackQuality::CheckAOD(AliAODTrack *track)
       else {
          if (track->Pt() < fPt[0] || track->Pt() > fPt[1]) return kFALSE;
          if (track->Eta() < fEta[0] || track->Eta() > fEta[1]) return kFALSE;
-         return kTRUE;
+         if (fCheckOnlyFilterBit) return kTRUE;
       }
    }
 
@@ -296,27 +375,45 @@ Bool_t AliRsnCutTrackQuality::CheckAOD(AliAODTrack *track)
    }
 
 
-   // step #1: check number of clusters in TPC
-   if (track->GetTPCNcls() < fTPCminNClusters) {
+   //step #1: check number of clusters 
+   if ((!fIsUseCrossedRowsCut) && (track->GetTPCNcls() < fTPCminNClusters)) {
       AliDebug(AliLog::kDebug + 2, "Too few TPC clusters. Rejected");
       return kFALSE;
    }
    if (track->GetITSNcls() < fITSminNClusters) {
       AliDebug(AliLog::kDebug + 2, "Too few ITS clusters. Rejected");
       return kFALSE;
    }
 
-   // step #2: check chi square
-   if (track->Chi2perNDF() > fTPCmaxChi2) {
-      AliDebug(AliLog::kDebug + 2, "Bad chi2. Rejected");
-      return kFALSE;
-   }
-   if (track->Chi2perNDF() > fITSmaxChi2) {
+   //check track chi square
+   if (track->Chi2perNDF() > fTrackMaxChi2) {
       AliDebug(AliLog::kDebug + 2, "Bad chi2. Rejected");
       return kFALSE;
    }
 
-   // step #3: reject kink daughters
+   //step #2a: check number of crossed rows in TPC
+   if (fIsUseCrossedRowsCut) {
+     Float_t nCrossedRowsTPC = track->GetTPCNCrossedRows();
+     if (nCrossedRowsTPC < fTPCminNCrossedRows) {
+       AliDebug(AliLog::kDebug + 2, "Too few TPC crossed rows. Rejected");
+       return kFALSE;
+     }
+     if (track->GetTPCNclsF()>0) {
+       Float_t ratioCrossedRowsOverFindableClustersTPC = nCrossedRowsTPC / track->GetTPCNclsF();
+       if (ratioCrossedRowsOverFindableClustersTPC < fTPCminCrossedRowsOverFindableCls){
+        AliDebug(AliLog::kDebug + 2, "Too few TPC crossed rows/findable clusters. Rejected");
+        return kFALSE;
+       }
+     } else {
+       AliDebug(AliLog::kDebug + 2, "Negative value for TPC crossed rows/findable clusters. Rejected");
+       return kFALSE;
+     }
+   }
+   //step #2b: check on track length in active volume of TPC implemented only for ESD tracks
+   //if (fIsUseLengthActiveVolumeTPCCut) { // not yet implemented in AODs}
+   //step #3: reject kink daughters
    AliAODVertex *vertex = track->GetProdVertex();
    if (vertex && fRejectKinkDaughters) {
       if (vertex->GetType() == AliAODVertex::kKink) {
@@ -342,25 +439,36 @@ Bool_t AliRsnCutTrackQuality::CheckAOD(AliAODTrack *track)
       return kFALSE;
    }
    // if the DCA cut is not fixed, compute current value
-   if (!fDCARfixed) {
-      static TString str(fDCARptFormula);
+   if (!fDCARmaxfixed) {
+      TString str(fDCARptFormula);
       str.ReplaceAll("pt", "x");
-      static const TFormula dcaXY(Form("%s_dcaXY", GetName()), str.Data());
+      TFormula dcaXY(Form("%s_dcaXY", GetName()), str.Data());
       fDCARmax = dcaXY.Eval(track->Pt());
    }
+   if (!fDCARminfixed) {   
+      TString str2(fDCARptFormulaMin);
+      str2.ReplaceAll("pt", "x");
+      TFormula dcaXY_2(Form("%s_dcaXY_2", GetName()), str2.Data());
+      fDCARmin = dcaXY_2.Eval(track->Pt());
+   }
    // check the cut
    if (TMath::Abs(b[0]) > fDCARmax) {
       AliDebug(AliLog::kDebug + 2, "Too large transverse DCA");
       return kFALSE;
    }
+   
+   if (TMath::Abs(b[0]) < fDCARmin) {
+      AliDebug(AliLog::kDebug + 2, "Too short transverse DCA");
+      return kFALSE;
+   }
 
    // step #5: DCA cut (longitudinal)
    // the DCA has already been computed above
    // if the DCA cut is not fixed, compute current value
    if (!fDCAZfixed) {
-      static TString str(fDCAZptFormula);
+      TString str(fDCAZptFormula);
       str.ReplaceAll("pt", "x");
-      static const TFormula dcaZ(Form("%s_dcaXY", GetName()), str.Data());
+      TFormula dcaZ(Form("%s_dcaXY", GetName()), str.Data());
       fDCAZmax = dcaZ.Eval(track->Pt());
    }
    // check the cut
@@ -395,13 +503,22 @@ void AliRsnCutTrackQuality::Print(const Option_t *) const
    AliInfo(Form("Required flags (off, on): %lx %lx", fFlagsOn, fFlagsOff));
    AliInfo(Form("Ranges in eta, pt       : %.2f - %.2f, %.2f - %.2f", fEta[0], fEta[1], fPt[0], fPt[1]));
    AliInfo(Form("Kink daughters are      : %s", (fRejectKinkDaughters ? "rejected" : "accepted")));
-   AliInfo(Form("TPC requirements        : min. cluster = %d, max chi2 = %f", fTPCminNClusters, fTPCmaxChi2));
+   AliInfo(Form("TPC requirements (clusters)       : min. cluster = %i, max chi2 = %f", fTPCminNClusters, fTPCmaxChi2));
+   AliInfo(Form("TPC requirements (crossed rows)   : min. crossed rows = %f, min. crossed rows/findable clusters = %f", fTPCminNCrossedRows, fTPCminCrossedRowsOverFindableCls));
+   AliInfo(Form("TPC requirements (track length)   : min. track length in active volume TPC = %f", fCutMinLengthActiveVolumeTPC));
+
    AliInfo(Form("ITS requirements        : min. cluster = %d (all), %d (SPD), max chi2 = %f", fITSminNClusters, fSPDminNClusters, fITSmaxChi2));
 
-   if (fDCARfixed) {
-      AliInfo(Form("DCA r cut               : fixed to %f cm", fDCARmax));
+   if (fDCARmaxfixed) {
+      AliInfo(Form("Max DCA r cut               : fixed to %f cm", fDCARmax));
    } else {
-      AliInfo(Form("DCA r cut formula       : %s", fDCARptFormula.Data()));
+      AliInfo(Form("Max DCA r cut formula       : %s", fDCARptFormula.Data()));
+   }
+   
+   if (fDCARminfixed) {
+      AliInfo(Form("Min DCA r cut               : fixed to %f cm", fDCARmin));
+   } else {
+      AliInfo(Form("Min DCA r cut formula       : %s", fDCARptFormulaMin.Data()));
    }
 
    if (fDCAZfixed) {
@@ -409,30 +526,60 @@ void AliRsnCutTrackQuality::Print(const Option_t *) const
    } else {
       AliInfo(Form("DCA z cut formula       : %s", fDCAZptFormula.Data()));
    }
+
+   AliInfo(Form("fAODTestFilterBit       : filter bit %i",fAODTestFilterBit));
+   AliInfo(Form("fCheckOnlyFilterBit     : %i",((int) fCheckOnlyFilterBit)));
 }
 //__________________________________________________________________________________________________
-void AliRsnCutTrackQuality::SetDefaults2010()
+void AliRsnCutTrackQuality::SetDefaults2010(Bool_t useTPCCrossedRows, Bool_t useDefaultKinematicCuts)
 {
 //
 // Default settings for cuts used in 2010
 //
-   AddStatusFlag(AliESDtrack::kTPCin   , kTRUE);
-   AddStatusFlag(AliESDtrack::kTPCrefit, kTRUE);
-   AddStatusFlag(AliESDtrack::kITSrefit, kTRUE);
-   SetPtRange(0.15, 1E+20);
-   SetEtaRange(-0.8, 0.8);
-   SetDCARPtFormula("0.0182+0.0350/pt^1.01");
-   SetDCAZmax(2.0);
-   SetSPDminNClusters(1);
-   SetITSminNClusters(0);
-   SetITSmaxChi2(36);
-   SetMaxChi2TPCConstrainedGlobal(36);
-   SetTPCminNClusters(70);
-   SetTPCmaxChi2(4.0);
-   SetRejectKinkDaughters();
-   SetAODTestFilterBit(5);
+  
+  fIsUseCrossedRowsCut=useTPCCrossedRows;
+  fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE, fIsUseCrossedRowsCut);
+  if (useDefaultKinematicCuts) {
+    SetPtRange(0.15, 1E+20);
+    SetEtaRange(-0.8, 0.8);
+  } 
+  SetAODTestFilterBit(5);
+  return;
+}
+
+//__________________________________________________________________________________________________
+void AliRsnCutTrackQuality::SetDefaultsHighPt2011(Bool_t useTPCCrossedRows, Bool_t useDefaultKinematicCuts)
+{
+//
+// Default settings for cuts used in 2011 (for high-pT)
+//
+  fIsUseCrossedRowsCut=useTPCCrossedRows;
+  fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kTRUE, fIsUseCrossedRowsCut);
+  fESDtrackCuts->SetMinNCrossedRowsTPC(120); //default is min 70 crossed rows -> use 120 to go to higher pt
+  fESDtrackCuts->SetMaxFractionSharedTPCClusters(0.4);//default is not set --> use to go to higher pt
+  if (useDefaultKinematicCuts) {
+    SetPtRange(0.15, 1E+20);
+    SetEtaRange(-0.8, 0.8);
+  } 
+  SetAODTestFilterBit(10);
+  return;
 }
 
+//__________________________________________________________________________________________________
+void AliRsnCutTrackQuality::SetDefaults2011(Bool_t useTPCCrossedRows, Bool_t useDefaultKinematicCuts)
+{
+//
+// Default std cuts 2011 with crossed rows (=70)
+//
+  fIsUseCrossedRowsCut=useTPCCrossedRows;
+  fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011(kTRUE,fIsUseCrossedRowsCut);
+  if (useDefaultKinematicCuts) {
+    SetPtRange(0.15, 1E+20);
+    SetEtaRange(-0.8, 0.8);
+  } 
+  SetAODTestFilterBit(5);
+  return;
+}
 //__________________________________________________________________________________________________
 const char *AliRsnCutTrackQuality::Binary(UInt_t number)
 {