]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
add optional restriction on leading hadron on opposite side of trigger
authorgconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 30 Sep 2013 12:30:38 +0000 (12:30 +0000)
committergconesab <gconesab@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 30 Sep 2013 12:30:38 +0000 (12:30 +0000)
PWGGA/CaloTrackCorrelations/AliAnaParticleHadronCorrelation.cxx
PWGGA/CaloTrackCorrelations/AliAnaParticleHadronCorrelation.h

index 287fe2f33d129a8a410ffc891f0799c6f8482f10..961cf62732bdcb220a3e2aa2a091dd5f085c7d20 100755 (executable)
@@ -74,7 +74,11 @@ ClassImp(AliAnaParticleHadronCorrelation)
     fListMixTrackEvents(),          fListMixCaloEvents(),
     fUseMixStoredInReader(0),       fFillNeutralEventMixPool(0),
     fM02MaxCut(0),                  fM02MinCut(0),  
-    fFillPileUpHistograms(0),       
+    fFillPileUpHistograms(0),
+    fSelectLeadingHadronAngle(0),
+    fMinLeadHadPhi(0),              fMaxLeadHadPhi(0),
+    fMinLeadHadPt(0),               fMaxLeadHadPt(0),
+
     //Histograms
     fhPtInput(0),                   fhPtFidCut(0),
     fhPtLeading(0),                 fhPtLeadingVtxBC0(0),
@@ -90,9 +94,9 @@ ClassImp(AliAnaParticleHadronCorrelation)
     fhDeltaPhiCharged(0),           fhDeltaEtaCharged(0), 
     fhDeltaPhiChargedPt(0),         fhDeltaPhiUeChargedPt(0), 
     fhUePart(0),
-    fhXECharged(0),                 fhXECharged_Cone2(0),            fhXEUeCharged(0),
+    fhXECharged(0),                 fhXECharged_Cone2(0),      fhXEUeCharged(0),
     fhXEPosCharged(0),              fhXENegCharged(0),
-    fhPtHbpXECharged(0),            fhPtHbpXECharged_Cone2(0),             fhPtHbpXEUeCharged(0),
+    fhPtHbpXECharged(0),            fhPtHbpXECharged_Cone2(0), fhPtHbpXEUeCharged(0),
     fhZTCharged(0),                 fhZTUeCharged(0),
     fhZTPosCharged(0),              fhZTNegCharged(0),
     fhPtHbpZTCharged(0),            fhPtHbpZTUeCharged(0),
@@ -2621,12 +2625,18 @@ void AliAnaParticleHadronCorrelation::InitParameters()
   fAssocPtBinLimit[18]  = 50.0 ;
   fAssocPtBinLimit[19]  = 200.0 ;
   
-  
   fUseMixStoredInReader = kTRUE;
   
   fM02MinCut   = -1 ;
   fM02MaxCut   = -1 ;
   
+  fSelectLeadingHadronAngle = kFALSE;
+  fMinLeadHadPhi = 150*TMath::DegToRad();
+  fMaxLeadHadPhi = 210*TMath::DegToRad();
+
+  fMinLeadHadPt  = 1;
+  fMaxLeadHadPt  = 100;
+  
 }
 
 //__________________________________________________________
@@ -2897,6 +2907,8 @@ Bool_t  AliAnaParticleHadronCorrelation::MakeChargedCorrelation(AliAODPWG4Partic
   Float_t eta      = -100. ;
   Float_t pout     = -100. ;
   Float_t deltaPhi = -100. ;
+  Float_t ptLeadHad  = -100 ;
+  Float_t phiLeadHad = -100 ;
   
   TVector3 p3;  
   TLorentzVector photonMom ;   
@@ -2927,9 +2939,38 @@ Bool_t  AliAnaParticleHadronCorrelation::MakeChargedCorrelation(AliAODPWG4Partic
   if(fPi0Trigger && bFillHisto) decayFound = GetDecayPhotonMomentum(aodParticle,decayMom1, decayMom2);
 
   //-----------------------------------------------------------------------
-  //Track loop, select tracks with good pt, phi and fill AODs or histograms
+  // Track loop, select tracks with good pt, phi and fill AODs or histograms
   //-----------------------------------------------------------------------
 
+  // select events where the leading particle in the opposite hemisphere to the trigger particle is
+  // is in a window centered at 180 from the trigger
+  // Find the leading hadron
+  if(fSelectLeadingHadronAngle)
+  {
+    for(Int_t ipr = 0;ipr < pl->GetEntriesFast() ; ipr ++ )
+    {
+      AliVTrack * track = (AliVTrack *) (pl->At(ipr)) ;
+      Double_t mom[3] = {track->Px(),track->Py(),track->Pz()};
+      p3.SetXYZ(mom[0],mom[1],mom[2]);
+      pt   = p3.Pt();
+      phi  = p3.Phi() ;
+
+      if(pt > ptLeadHad && TMath::Abs(phi-phiTrig) > TMath::PiOver2())
+      {
+        ptLeadHad = pt ;
+        phiLeadHad = phi;
+      }
+      
+    }// track loop
+    
+    if( ptLeadHad < fMinLeadHadPt ||
+        ptLeadHad > fMaxLeadHadPt ) return kFALSE;
+    
+    if( TMath::Abs(phiLeadHad-phiTrig) < fMinLeadHadPhi ||
+        TMath::Abs(phiLeadHad-phiTrig) > fMaxLeadHadPhi ) return kFALSE;
+    
+  }// select leading hadron
+  
   for(Int_t ipr = 0;ipr < pl->GetEntriesFast() ; ipr ++ )
   {
     AliVTrack * track = (AliVTrack *) (pl->At(ipr)) ;
@@ -2959,7 +3000,7 @@ Bool_t  AliAnaParticleHadronCorrelation::MakeChargedCorrelation(AliAODPWG4Partic
     {
       if(pt > ptTrig)  return kFALSE;
     }
-    
+        
     //Only for mixed event
     Int_t evtIndex2 = 0 ; 
     if (GetMixedEvent()) 
index ed66bb3242bb45492bc8bc4fe14222d5d479a8a4..f4e5d5b00acf48b6ecfd36cd0be5ca0f4b2ce41a 100755 (executable)
@@ -106,23 +106,37 @@ class AliAnaParticleHadronCorrelation : public AliAnaCaloTrackCorrBaseClass {
   // Parameter setter and getter
   
   Float_t      GetMinimumTriggerPt()       const { return fMinTriggerPt          ; }
-  
+  void         SetMinimumTriggerPt(Float_t min)  { fMinTriggerPt = min           ; }
+
   Float_t      GetMaximumAssociatedPt()    const { return fMaxAssocPt            ; }
   Float_t      GetMinimumAssociatedPt()    const { return fMinAssocPt            ; }
-  
+  void         SetAssociatedPtRange(Float_t min, Float_t max)
+               { fMaxAssocPt   = max ;           fMinAssocPt  = min              ; }
+
   Double_t     GetDeltaPhiMaxCut()         const { return fDeltaPhiMaxCut        ; }
   Double_t     GetDeltaPhiMinCut()         const { return fDeltaPhiMinCut        ; }
+  void         SetDeltaPhiCutRange(Double_t phimin, Double_t phimax)
+               { fDeltaPhiMaxCut   = phimax ;    fDeltaPhiMinCut   = phimin      ; }
   
-  Double_t     GetUeDeltaPhiMaxCut()       const { return fUeDeltaPhiMaxCut      ; }
-  Double_t     GetUeDeltaPhiMinCut()       const { return fUeDeltaPhiMinCut      ; }
+  // Leading Hadron
+  Double_t     GetLeadHadronPhiMaxCut()    const { return fMaxLeadHadPhi         ; }
+  Double_t     GetLeadHadronPhiMinCut()    const { return fMinLeadHadPhi         ; }
+  void         SetLeadHadronPhiCut(Float_t min, Float_t max)
+               { fMaxLeadHadPhi = max ;          fMinLeadHadPhi  = min           ; }
+
+  Double_t     GetLeadHadronPtMinCut()     const { return fMinLeadHadPt          ; }
+  Double_t     GetLeadHadronPtMaxCut()     const { return fMaxLeadHadPt          ; }
+  void         SetLeadHadronPtCut(Float_t min, Float_t max)
+               { fMaxLeadHadPt  = max ;           fMinLeadHadPt  = min           ; }
   
-  void         SetMinimumTriggerPt(Float_t min)  { fMinTriggerPt = min           ; }
+  Bool_t       IsLeadHadronCutOn()        const { return fSelectLeadingHadronAngle ; }
+  void         SwitchOnLeadHadronSelection()    { fSelectLeadingHadronAngle = kTRUE  ; }
+  void         SwitchOffLeadHadronSelection()   { fSelectLeadingHadronAngle = kFALSE ; }
   
-  void         SetAssociatedPtRange(Float_t min, Float_t max)
-                  { fMaxAssocPt   = max ;           fMinAssocPt  = min           ; }
+  // UE
   
-  void         SetDeltaPhiCutRange(Double_t phimin, Double_t phimax)
-                  { fDeltaPhiMaxCut   = phimax ;    fDeltaPhiMinCut   = phimin   ; }
+  Double_t     GetUeDeltaPhiMaxCut()       const { return fUeDeltaPhiMaxCut      ; }
+  Double_t     GetUeDeltaPhiMinCut()       const { return fUeDeltaPhiMinCut      ; }
   
   void         SetUeDeltaPhiCutRange(Double_t uephimin, Double_t uephimax)
                   { fUeDeltaPhiMaxCut = uephimax ;  fUeDeltaPhiMinCut = uephimin ; }
@@ -220,6 +234,12 @@ class AliAnaParticleHadronCorrelation : public AliAnaCaloTrackCorrBaseClass {
   
   Bool_t       fFillPileUpHistograms;          // Fill pile-up related histograms
 
+  Bool_t       fSelectLeadingHadronAngle;      // Select events with leading particle within a range
+  Float_t      fMinLeadHadPhi;                 // Minimum angle between the trigger and leading hadron
+  Float_t      fMaxLeadHadPhi;                 // Maximum ange between the trigger and leading hadron
+  Float_t      fMinLeadHadPt;                  // Minimum pT of leading hadron
+  Float_t      fMaxLeadHadPt;                  // Maximum pT of leading hadron
+
   //Histograms
 
   //leading particles 
@@ -439,7 +459,7 @@ class AliAnaParticleHadronCorrelation : public AliAnaCaloTrackCorrBaseClass {
   AliAnaParticleHadronCorrelation(              const AliAnaParticleHadronCorrelation & ph) ; // cpy ctor
   AliAnaParticleHadronCorrelation & operator = (const AliAnaParticleHadronCorrelation & ph) ; // cpy assignment
        
-  ClassDef(AliAnaParticleHadronCorrelation,29)
+  ClassDef(AliAnaParticleHadronCorrelation,30)
 } ;