]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWGCF/Correlations/Base/AliTwoPlusOneContainer.cxx
2+1 analysis update from Markus
[u/mrichter/AliRoot.git] / PWGCF / Correlations / Base / AliTwoPlusOneContainer.cxx
index beac5ac123699dd65c66d542723e1b378b3b8db4..f92153709c39015ac50a540a08d21509f35aeb0e 100644 (file)
@@ -35,6 +35,7 @@
 #include "AliVParticle.h"
 
 #include "TH1F.h"
+#include "TH2F.h"
 #include "TMath.h"
 
 ClassImp(AliTwoPlusOneContainer)
@@ -42,6 +43,9 @@ ClassImp(AliTwoPlusOneContainer)
 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* binning, Double_t alpha) : 
   TNamed(name, name),
   fTwoPlusOne(0),
+  fAsymmetry(0),
+  fAsymmetryMixed(0),
+  fTriggerPt(0),
   fTriggerPt1Min(0),
   fTriggerPt1Max(0),
   fTriggerPt2Min(0),
@@ -49,6 +53,7 @@ AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* bin
   fPtAssocMin(0),
   fPtAssocMax(0),
   fAlpha(alpha),
+  fUseLeadingPt(1),
   fMergeCount(0)
 {
   // Constructor
@@ -77,7 +82,13 @@ AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* bin
   fTriggerPt2Max = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(6)->GetXmax();
   fPtAssocMin = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(1)->GetXmin();
   fPtAssocMax = fTwoPlusOne->GetTrackHist(AliUEHist::kToward)->GetGrid((AliUEHist::CFStep) AliTwoPlusOneContainer::kSameNS)->GetAxis(1)->GetXmax();
-
+  
+  fAsymmetry  = new TH1F("fAsymmetry", ";A;dN", 50, 0, 1);
+  fAsymmetryMixed  = new TH1F("fAsymmetryMixed", ";A;dN", 50, 0, 1);
+  Int_t pt1_bins = (fTriggerPt1Max - fTriggerPt1Min)/0.1;
+  Int_t pt2_bins = (fTriggerPt2Max - fTriggerPt2Min)/0.1;
+  fTriggerPt = new TH2F("fTriggerPt", ";p_{T,1};p_{T,2}", pt1_bins, fTriggerPt1Min, fTriggerPt1Max, pt2_bins, fTriggerPt2Min, fTriggerPt2Max);
+  
   TH1::AddDirectory();
 }
 
@@ -85,6 +96,9 @@ AliTwoPlusOneContainer::AliTwoPlusOneContainer(const char* name, const char* bin
 AliTwoPlusOneContainer::AliTwoPlusOneContainer(const AliTwoPlusOneContainer &c) :
   TNamed(fName, fTitle),
   fTwoPlusOne(0),
+  fAsymmetry(0),
+  fAsymmetryMixed(0),
+  fTriggerPt(0),
   fTriggerPt1Min(0),
   fTriggerPt1Max(0),
   fTriggerPt2Min(0),
@@ -125,6 +139,24 @@ void AliTwoPlusOneContainer::DeleteContainers()
     delete fTwoPlusOne;
     fTwoPlusOne = 0;
   }
+  
+  if(fAsymmetry)
+  {
+    delete fAsymmetry;
+    fAsymmetry = 0;
+  }
+
+  if(fAsymmetryMixed)
+  {
+    delete fAsymmetryMixed;
+    fAsymmetryMixed = 0;
+  }
+
+  if(fTriggerPt)
+  {
+    delete fTriggerPt;
+    fTriggerPt = 0;
+  }
 }
 
 
@@ -161,6 +193,34 @@ void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx,
     Int_t triggerAway_entries = triggerAway->GetEntriesFast();
     AliVParticle* found_particle[triggerAway_entries];
 
+    Bool_t do_not_use_T1 = false;
+
+    //in case only the leading pt of a jet should be used, check every particle on the trigger near side if it's closer than alpha and if it has a higher pt than trigger 1
+    if(fUseLeadingPt){
+      for (Int_t i2=0; i2<triggerNear->GetEntriesFast(); i2++){
+       if(i==i2)
+         continue;
+
+       AliVParticle* part_i2 = (AliVParticle*) triggerNear->UncheckedAt(i2);
+
+       if(part_i2->Pt()<=part_pt)
+         continue;
+    
+       Double_t dphi_check = part_phi-part_i2->Phi(); 
+       if(dphi_check>1.5*TMath::Pi()) dphi_check -= TMath::TwoPi();
+       else if(dphi_check<-0.5*TMath::Pi()) dphi_check += TMath::TwoPi();
+
+       if(TMath::Abs(dphi_check)<fAlpha){
+         do_not_use_T1 = true;
+         break;
+       }
+      }
+    }
+
+    //if there is a particle with higher energy than T1 closer than alpha to T1, do not use this T1
+     if(do_not_use_T1)
+      continue;
+
     //have to fake the away side triggers for the 1+1 analysis
     if(is1plus1){
       found_particle[ind_found] = part;//in 1plus1 use first trigger particle also as pseudo second trigger particle
@@ -172,18 +232,17 @@ void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx,
        AliVParticle* part2 = (AliVParticle*) triggerAway->UncheckedAt(j);
 
        Double_t part2_pt = part2->Pt();
-       //check if pT of trigger 2 is within the trigger range
-       //also pt of trigger 2 needs to be smaller than the pt of trigger 1 (to have an ordering if both pt are close to each other)
-       if(part2_pt<fTriggerPt2Min || part2_pt>fTriggerPt2Max || part2_pt>part_pt)
+       //check if pT of trigger 2 has enough energy to be a trigger
+       //maximum energy is checked later after checking this particle may have to much energy for trigger 1 or 2
+       if(part2_pt<fTriggerPt2Min)
          continue;
-      
+
        // don't use the same particle (is in any case impossible because the Delta phi angle will be 0)
        if(part==part2){
          continue;
        }
        
        Double_t dphi_triggers = part_phi-part2->Phi();
-       
        if(dphi_triggers>1.5*TMath::Pi()) dphi_triggers -= TMath::TwoPi();
        else if(dphi_triggers<-0.5*TMath::Pi()) dphi_triggers += TMath::TwoPi();
       
@@ -203,37 +262,75 @@ void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx,
        if(TMath::Abs(dphi_triggers)>fAlpha)
          continue;
 
+       //check if pT of trigger 2 is too high
+       if(part2_pt>fTriggerPt2Max || part2_pt>part_pt){
+         //pt of trigger 2 needs to be smaller than the pt of trigger 1 (to have an ordering if both pt are close to each other)
+         if(fUseLeadingPt){
+           do_not_use_T1 = true;
+           break;
+         }else
+           continue;
+       }
+
        found_particle[ind_found] = part2;
        if(ind_max_found_pt==-1 || part2_pt>found_particle[ind_max_found_pt]->Pt()) ind_max_found_pt = ind_found;
        ind_found++;
       }//end loop to search for the second trigger particle
     }
 
+    //if there is a particle with higher energy than T1 or max(T2) within Delta phi = pi +/- alpha to T1, do not use this T1
+    if(do_not_use_T1)
+      continue;
+
     //if no second trigger particle was found continue to search for the next first trigger particle
     if(ind_found==0)
       continue;
 
+    //use only the highest energetic particle on the away side, if there is only 1 away side trigger this is already the case
+    if(fUseLeadingPt && ind_found>1){
+      found_particle[0] = found_particle[ind_max_found_pt];
+      ind_found=1;
+      ind_max_found_pt = 0;
+    }
+
     //the energy of the second trigger particle is set for the near side to the maximum energy of all trigger 2 particles on the away side
     // this leads to the fact that the number of accepted trigger combinations can be artificial smaller than the real number if there is a cut on the pT 2 energy from the top; cutting away the smallest energy of pT 2 is still save; this is the reason why it is not allowed to use a cut on the top pt of trigger particle 2
     //fill trigger particles
     if(ind_found>0){
-       Double_t vars[4];
-       vars[0] = part_pt;
-       vars[1] = centrality;
-       vars[2] = zVtx;
-       vars[3] = found_particle[ind_max_found_pt]->Pt();
-       if(is1plus1)
-         vars[3] = (fTriggerPt2Max+fTriggerPt2Min)/2;
-
-       event_hist->Fill(vars, stepUEHist, weight);//near side
-
-       if(!is1plus1)
-         for(Int_t k=0; k< ind_found; k++){
-           vars[3] = found_particle[k]->Pt();
-           event_hist->Fill(vars, stepUEHist+1, weight);//away side
+      Double_t vars[4];
+      vars[0] = part_pt;
+      vars[1] = centrality;
+      vars[2] = zVtx;
+      vars[3] = found_particle[ind_max_found_pt]->Pt();
+      if(is1plus1)
+       vars[3] = (fTriggerPt2Max+fTriggerPt2Min)/2;
+      
+      event_hist->Fill(vars, stepUEHist, weight);//near side
+      
+      if(!is1plus1)
+       for(Int_t k=0; k< ind_found; k++){
+         vars[3] = found_particle[k]->Pt();
+         event_hist->Fill(vars, stepUEHist+1, weight);//away side
+       }
+       
+      //fill fTriggerPt only once, choosed kSameNS
+      if(step==AliTwoPlusOneContainer::kSameNS)
+       for(Int_t k=0; k< ind_found; k++)
+         fTriggerPt->Fill(part_pt, found_particle[k]->Pt());
+
+      //fill asymmetry only for kSameNS and kMixedNS
+      if(step==AliTwoPlusOneContainer::kSameNS||step==AliTwoPlusOneContainer::kMixedNS){
+       for(Int_t k=0; k< ind_found; k++){
+         Float_t asymmetry = (part_pt-found_particle[k]->Pt())/(part_pt+found_particle[k]->Pt());
+         if(step==AliTwoPlusOneContainer::kSameNS){
+           fAsymmetry->Fill(asymmetry);
+         }else{
+           fAsymmetryMixed->Fill(asymmetry);
          }
+       }
+      }
     }
+    
     //add correlated particles on the near side
     for (Int_t k=0; k<assocNear->GetEntriesFast(); k++){
       AliVParticle* part3 = (AliVParticle*) assocNear->UncheckedAt(k);
@@ -264,6 +361,8 @@ void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx,
       vars[4] = dphi_near;
       vars[5] = zVtx;
       vars[6] = found_particle[ind_max_found_pt]->Pt();
+      if(is1plus1)
+       vars[6] = (fTriggerPt2Max+fTriggerPt2Min)/2;
 
       track_hist->Fill(vars, stepUEHist, weight);
     }
@@ -310,8 +409,6 @@ void AliTwoPlusOneContainer::FillCorrelations(Double_t centrality, Float_t zVtx,
   }//end loop to search for the first trigger particle
 }
 
-
-
 //____________________________________________________________________
 AliTwoPlusOneContainer &AliTwoPlusOneContainer::operator=(const AliTwoPlusOneContainer &c)
 {
@@ -334,7 +431,16 @@ void AliTwoPlusOneContainer::Copy(TObject& c) const
 
   if (fTwoPlusOne)
     target.fTwoPlusOne = dynamic_cast<AliUEHist*> (fTwoPlusOne->Clone());
+  
+  if (fAsymmetry)
+    target.fAsymmetry = dynamic_cast<TH1F*> (fAsymmetry->Clone());
 
+  if (fAsymmetryMixed)
+    target.fAsymmetryMixed = dynamic_cast<TH1F*> (fAsymmetryMixed->Clone());
+
+  if (fTriggerPt)
+    target.fTriggerPt = dynamic_cast<TH2F*> (fTriggerPt->Clone());
+  
 }
 
 //____________________________________________________________________
@@ -351,9 +457,13 @@ Long64_t AliTwoPlusOneContainer::Merge(TCollection* list)
 
   TIterator* iter = list->MakeIterator();
   TObject* obj;
-
+  
   // collections of objects
-  TList* lists = new TList;
+  const Int_t kMaxLists = 4;
+  TList* lists[kMaxLists];
+
+  for (Int_t i=0; i<kMaxLists; i++)
+    lists[i] = new TList;
 
   Int_t count = 0;
   while ((obj = iter->Next())) {
@@ -362,15 +472,24 @@ Long64_t AliTwoPlusOneContainer::Merge(TCollection* list)
     if (entry == 0) 
       continue;
 
-    lists->Add(entry->fTwoPlusOne);
+    lists[0]->Add(entry->fTwoPlusOne);
+    lists[1]->Add(entry->fAsymmetry);
+    lists[2]->Add(entry->fAsymmetryMixed);
+    lists[3]->Add(entry->fTriggerPt);
 
     fMergeCount += entry->fMergeCount;
     count++;
   }
   
-  fTwoPlusOne->Merge(lists);
+  fTwoPlusOne->Merge(lists[0]);
+  fAsymmetry->Merge(lists[1]);
+  fAsymmetryMixed->Merge(lists[2]);
+  fTriggerPt->Merge(lists[3]);
+
+  for (Int_t i=0; i<kMaxLists; i++)
+  delete lists[i];
 
-  delete lists;
+  //  delete lists;
   return count+1;
 
 }