]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - PWG2/FLOW/AliFlowCommon/AliFlowCommonHist.cxx
FMD and ITS additions
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowCommonHist.cxx
index a57f2e63e2e3e408287dc131722edbcd9b7073ab..6f91cf92100e9d7ccac22aa5b3d29e5cda59b9d0 100644 (file)
@@ -23,6 +23,7 @@
 #include "TProfile.h"
 #include "TMath.h"   //needed as include
 #include "TList.h"
+#include "TH2F.h"
 #include "AliFlowVector.h"
 #include "TBrowser.h"
 
@@ -40,16 +41,25 @@ ClassImp(AliFlowCommonHist)
 
 //-----------------------------------------------------------------------
 
-AliFlowCommonHist::AliFlowCommonHist():TNamed(),
+AliFlowCommonHist::AliFlowCommonHist():
+  TNamed(),
   fHistMultOrig(NULL),
   fHistMultRP(NULL),
   fHistMultPOI(NULL),
   fHistPtRP(NULL),
   fHistPtPOI(NULL),
+  fHistPtSub0(NULL),
+  fHistPtSub1(NULL),
   fHistPhiRP(NULL),
   fHistPhiPOI(NULL),
+  fHistPhiSub0(NULL),
+  fHistPhiSub1(NULL),
   fHistEtaRP(NULL),
   fHistEtaPOI(NULL),
+  fHistEtaSub0(NULL),
+  fHistEtaSub1(NULL),
+  fHistPhiEtaRP(NULL),
+  fHistPhiEtaPOI(NULL),
   fHistProMeanPtperBin(NULL),
   fHistQ(NULL),
   fHarmonic(NULL),
@@ -67,10 +77,18 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistMultPOI(new TH1F(*a.fHistMultPOI)),
   fHistPtRP(new TH1F(*a.fHistPtRP)),
   fHistPtPOI(new TH1F(*a.fHistPtPOI)),
+  fHistPtSub0(new TH1F(*a.fHistPtSub0)),
+  fHistPtSub1(new TH1F(*a.fHistPtSub1)),
   fHistPhiRP(new TH1F(*a.fHistPhiRP)),
   fHistPhiPOI(new TH1F(*a.fHistPhiPOI)),
+  fHistPhiSub0(new TH1F(*a.fHistPhiSub0)),
+  fHistPhiSub1(new TH1F(*a.fHistPhiSub1)),
   fHistEtaRP(new TH1F(*a.fHistEtaRP)),
   fHistEtaPOI(new TH1F(*a.fHistEtaPOI)),
+  fHistEtaSub0(new TH1F(*a.fHistEtaSub0)),
+  fHistEtaSub1(new TH1F(*a.fHistEtaSub1)),
+  fHistPhiEtaRP(new TH2F(*a.fHistPhiEtaRP)),
+  fHistPhiEtaPOI(new TH2F(*a.fHistPhiEtaPOI)),
   fHistProMeanPtperBin(new TProfile(*a.fHistProMeanPtperBin)),
   fHistQ(new TH1F(*a.fHistQ)),
   fHarmonic(new TProfile(*a.fHarmonic)),
@@ -83,11 +101,19 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistList-> Add(fHistMultRP);        
   fHistList-> Add(fHistMultPOI);       
   fHistList-> Add(fHistPtRP);          
-  fHistList-> Add(fHistPtPOI);         
+  fHistList-> Add(fHistPtPOI);
+  fHistList-> Add(fHistPtSub0);
+  fHistList-> Add(fHistPtSub1);
   fHistList-> Add(fHistPhiRP);          
-  fHistList-> Add(fHistPhiPOI);         
+  fHistList-> Add(fHistPhiPOI);
+  fHistList-> Add(fHistPhiSub0);
+  fHistList-> Add(fHistPhiSub1);    
   fHistList-> Add(fHistEtaRP);          
-  fHistList-> Add(fHistEtaPOI);         
+  fHistList-> Add(fHistEtaPOI); 
+  fHistList-> Add(fHistEtaSub0);
+  fHistList-> Add(fHistEtaSub1);
+  fHistList-> Add(fHistPhiEtaRP);
+  fHistList-> Add(fHistPhiEtaPOI);
   fHistList-> Add(fHistProMeanPtperBin); 
   fHistList-> Add(fHarmonic);  
   fHistList-> Add(fHistQ);           
@@ -116,21 +142,30 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
 
 //-----------------------------------------------------------------------
 
-  AliFlowCommonHist::AliFlowCommonHist(const char *anInput,const char *title):TNamed(anInput,title),
-   fHistMultOrig(NULL),
-   fHistMultRP(NULL),
-   fHistMultPOI(NULL),
-   fHistPtRP(NULL),
-   fHistPtPOI(NULL),
-   fHistPhiRP(NULL),
-   fHistPhiPOI(NULL),
-   fHistEtaRP(NULL),
-   fHistEtaPOI(NULL),
-   fHistProMeanPtperBin(NULL),
-   fHistQ(NULL),
-   fHarmonic(NULL),
-   fHistList(NULL)
- {
+  AliFlowCommonHist::AliFlowCommonHist(const char *anInput,const char *title):
+    TNamed(anInput,title),
+    fHistMultOrig(NULL),
+    fHistMultRP(NULL),
+    fHistMultPOI(NULL),
+    fHistPtRP(NULL),
+    fHistPtPOI(NULL),
+    fHistPtSub0(NULL),
+    fHistPtSub1(NULL),
+    fHistPhiRP(NULL),
+    fHistPhiPOI(NULL),
+    fHistPhiSub0(NULL),
+    fHistPhiSub1(NULL),
+    fHistEtaRP(NULL),
+    fHistEtaPOI(NULL),
+    fHistEtaSub0(NULL),
+    fHistEtaSub1(NULL),
+    fHistPhiEtaRP(NULL),
+    fHistPhiEtaPOI(NULL),
+    fHistProMeanPtperBin(NULL),
+    fHistQ(NULL),
+    fHarmonic(NULL),
+    fHistList(NULL)
+{
 
   //constructor creating histograms 
   Int_t iNbinsMult = AliFlowCommonConstants::GetMaster()->GetNbinsMult();
@@ -191,6 +226,18 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistPtPOI ->SetXTitle("P_{t} (GeV/c) for POI selection");
   fHistPtPOI ->SetYTitle("Counts");
 
+  sName = "Control_Flow_PtSub0_";
+  sName +=anInput;
+  fHistPtSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax); 
+  fHistPtSub0 ->SetXTitle("P_{t} (GeV/c) for Subevent 0 selection");
+  fHistPtSub0 ->SetYTitle("Counts");
+
+  sName = "Control_Flow_PtSub1_";
+  sName +=anInput;
+  fHistPtSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax); 
+  fHistPtSub1 ->SetXTitle("P_{t} (GeV/c) for Subevent 1 selection");
+  fHistPtSub1 ->SetYTitle("Counts");
+
   //Phi
   sName = "Control_Flow_PhiRP_";
   sName +=anInput;
@@ -204,6 +251,18 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistPhiPOI ->SetXTitle("#phi for POI selection");
   fHistPhiPOI ->SetYTitle("Counts");
 
+  sName = "Control_Flow_PhiSub0_";
+  sName +=anInput;
+  fHistPhiSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
+  fHistPhiSub0 ->SetXTitle("#phi for Subevent 0 selection");
+  fHistPhiSub0 ->SetYTitle("Counts");
+
+  sName = "Control_Flow_PhiSub1_";
+  sName +=anInput;
+  fHistPhiSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
+  fHistPhiSub1 ->SetXTitle("#phi for Subevent 1 selection");
+  fHistPhiSub1 ->SetYTitle("Counts");
+
   //Eta
   sName = "Control_Flow_EtaRP_";
   sName +=anInput;
@@ -217,6 +276,31 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistEtaPOI ->SetXTitle("#eta for POI selection");
   fHistEtaPOI ->SetYTitle("Counts");
 
+  sName = "Control_Flow_EtaSub0_";
+  sName +=anInput;
+  fHistEtaSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
+  fHistEtaSub0 ->SetXTitle("#eta for Subevent 0 selection");
+  fHistEtaSub0 ->SetYTitle("Counts");
+
+  sName = "Control_Flow_EtaSub1_";
+  sName +=anInput;
+  fHistEtaSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
+  fHistEtaSub1 ->SetXTitle("#eta for Subevent 1 selection");
+  fHistEtaSub1 ->SetYTitle("Counts");
+
+  //Phi vs Eta
+  sName = "Control_Flow_PhiEtaRP_";
+  sName +=anInput;
+  fHistPhiEtaRP = new TH2F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax, iNbinsPhi, dPhiMin, dPhiMax);
+  fHistPhiEtaRP ->SetXTitle("#eta");
+  fHistPhiEtaRP ->SetYTitle("#phi");
+
+  sName = "Control_Flow_PhiEtaPOI_";
+  sName +=anInput;
+  fHistPhiEtaPOI = new TH2F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax, iNbinsPhi, dPhiMin, dPhiMax);
+  fHistPhiEtaPOI ->SetXTitle("#eta");
+  fHistPhiEtaPOI ->SetYTitle("#phi");
+
   //Mean Pt per pt bin 
   sName = "Control_FlowPro_MeanPtperBin_";
   sName +=anInput;
@@ -243,16 +327,23 @@ AliFlowCommonHist::AliFlowCommonHist(const AliFlowCommonHist& a):
   fHistList-> Add(fHistMultRP);        
   fHistList-> Add(fHistMultPOI);       
   fHistList-> Add(fHistPtRP);          
-  fHistList-> Add(fHistPtPOI);         
+  fHistList-> Add(fHistPtPOI); 
+  fHistList-> Add(fHistPtSub0);
+  fHistList-> Add(fHistPtSub1);
   fHistList-> Add(fHistPhiRP);          
-  fHistList-> Add(fHistPhiPOI);         
+  fHistList-> Add(fHistPhiPOI);
+  fHistList-> Add(fHistPhiSub0);
+  fHistList-> Add(fHistPhiSub1);
   fHistList-> Add(fHistEtaRP);          
-  fHistList-> Add(fHistEtaPOI);         
+  fHistList-> Add(fHistEtaPOI); 
+  fHistList-> Add(fHistEtaSub0); 
+  fHistList-> Add(fHistEtaSub1);
+  fHistList-> Add(fHistPhiEtaRP);  
+  fHistList-> Add(fHistPhiEtaPOI);
   fHistList-> Add(fHistProMeanPtperBin);
   fHistList-> Add(fHarmonic);  
   fHistList-> Add(fHistQ);           
 
-
 }
 
 
@@ -265,11 +356,19 @@ AliFlowCommonHist::~AliFlowCommonHist()
   delete fHistMultRP;       
   delete fHistMultPOI;      
   delete fHistPtRP;         
-  delete fHistPtPOI;       
+  delete fHistPtPOI;
+  delete fHistPtSub0;
+  delete fHistPtSub1;
   delete fHistPhiRP;        
-  delete fHistPhiPOI;       
+  delete fHistPhiPOI;
+  delete fHistPhiSub0;
+  delete fHistPhiSub1;
   delete fHistEtaRP;        
   delete fHistEtaPOI;
+  delete fHistEtaSub0;
+  delete fHistEtaSub1;
+  delete fHistPhiEtaRP;
+  delete fHistPhiEtaPOI;
   delete fHistProMeanPtperBin;
   delete fHistQ;
   delete fHarmonic;
@@ -314,24 +413,49 @@ Bool_t AliFlowCommonHist::FillControlHistograms(AliFlowEventSimple* anEvent)
     pTrack = anEvent->GetTrack(i);
     if (pTrack ) {
       if (pTrack->InRPSelection()){
+       //pt
        dPt = pTrack->Pt();
        fHistPtRP->Fill(dPt);
+       //phi
        dPhi = pTrack->Phi();
        if (dPhi<0.) dPhi+=2*TMath::Pi();
        fHistPhiRP->Fill(dPhi);
+       //eta
        dEta = pTrack->Eta();
        fHistEtaRP->Fill(dEta);
+       //eta vs phi
+       fHistPhiEtaRP->Fill(dEta,dPhi);
+       //count
        iMultRP++;
+       if (pTrack->InSubevent(0)){
+         //Fill distributions for the subevent
+         fHistPtSub0 -> Fill(dPt);
+         fHistPhiSub0 -> Fill(dPhi);
+         fHistEtaSub0 -> Fill(dEta);
+       } 
+       else if (pTrack->InSubevent(1)){
+         //Fill distributions for the subevent
+         fHistPtSub1 -> Fill(dPt);
+         fHistPhiSub1 -> Fill(dPhi);
+         fHistEtaSub1 -> Fill(dEta);
+       } 
       }
       if (pTrack->InPOISelection()){
+       //pt
        dPt = pTrack->Pt();
        fHistPtPOI->Fill(dPt);
+       //phi
        dPhi = pTrack->Phi();
        if (dPhi<0.) dPhi+=2*TMath::Pi();
        fHistPhiPOI->Fill(dPhi);
+       //eta
        dEta = pTrack->Eta();
        fHistEtaPOI->Fill(dEta);
+       //eta vs phi
+       fHistPhiEtaPOI->Fill(dEta,dPhi);
+       //mean pt
        fHistProMeanPtperBin->Fill(dPt,dPt);
+       //count
        iMultPOI++;
       }
     } //track