]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Example for mixing analysis.
authormorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 12 Aug 2008 11:25:07 +0000 (11:25 +0000)
committermorsch <morsch@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 12 Aug 2008 11:25:07 +0000 (11:25 +0000)
ANALYSIS/AliAnalysisTaskPhiCorr.cxx [new file with mode: 0644]
ANALYSIS/AliAnalysisTaskPhiCorr.h [new file with mode: 0644]

diff --git a/ANALYSIS/AliAnalysisTaskPhiCorr.cxx b/ANALYSIS/AliAnalysisTaskPhiCorr.cxx
new file mode 100644 (file)
index 0000000..3e23877
--- /dev/null
@@ -0,0 +1,108 @@
+
+#include "TH1F.h"
+#include "TList.h"
+#include "TCanvas.h"
+#include "AliAODEvent.h"
+#include "AliAODTrack.h"
+#include "AliAnalysisTaskPhiCorr.h"
+#include "AliMultiAODInputHandler.h"
+
+ClassImp(AliAnalysisTaskPhiCorr)
+
+//________________________________________________________________________
+AliAnalysisTaskPhiCorr::AliAnalysisTaskPhiCorr(const char *name) 
+    : AliAnalysisTaskME(name), fHists(0), fHistDphiCO(0),  fHistDphiUC(0)
+{
+  // Constructor
+  // Define input and output slots here
+  DefineOutput(1, TList::Class());
+}
+
+
+//________________________________________________________________________
+void AliAnalysisTaskPhiCorr::UserCreateOutputObjects()
+{
+  // Create histograms
+  // Called once
+
+  fHistDphiCO = new TH1F("fHistDPhiCO", "#Delta #Phi distribution", 64, 0., 3.2);
+  fHistDphiCO->GetXaxis()->SetTitle("#Delta#Phi [rad]");
+  fHistDphiCO->GetYaxis()->SetTitle("dN/d#Phi");
+  fHistDphiCO->SetMarkerStyle(kFullCircle);
+
+  fHistDphiUC = new TH1F("fHistDPhiUC", "#Delta #Phi distribution", 64, 0., 3.2);
+  fHistDphiUC->GetXaxis()->SetTitle("#Delta#Phi [rad]");
+  fHistDphiUC->GetYaxis()->SetTitle("dN/d#Phi");
+  fHistDphiUC->SetMarkerStyle(kOpenCircle);
+  
+  fHists = new TList();
+  fHists->Add(fHistDphiCO);
+  fHists->Add(fHistDphiUC);  
+}
+
+//________________________________________________________________________
+void AliAnalysisTaskPhiCorr::UserExec(Option_t *) 
+{
+  // Main loop
+  // Called for each event
+    // Uncorrelated tracks
+    Int_t nev = fInputHandler->GetBufferSize();
+    Float_t wgt = 1./(nev*(nev-1));
+    
+    for (Int_t iev = 0; iev < nev; iev++) {
+       for (Int_t jev = (iev + 1); jev < nev; jev++) {
+           AliAODEvent* aod1 = GetEvent(iev);
+           AliAODEvent* aod2 = GetEvent(jev);
+           
+           Int_t ntracks1 =  aod1->GetNumberOfTracks();
+           Int_t ntracks2 =  aod2->GetNumberOfTracks();
+           
+           printf("Number of tracks %5d:%5d %5d:%5d\n", iev, ntracks1, jev, ntracks2);
+           
+           for (Int_t iTracks = 0; iTracks < ntracks1; iTracks++) {
+               for (Int_t jTracks = 0; jTracks < ntracks2; jTracks++) {
+                   
+                   
+                   AliAODTrack* track1 = aod1->GetTrack(iTracks);
+                   AliAODTrack* track2 = aod2->GetTrack(jTracks);
+                   
+                   Float_t phi1 = track1->Phi();
+                   Float_t phi2 = track2->Phi();
+                   Float_t dphi = TMath::Abs(phi1 - phi2);
+                   if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
+                   fHistDphiUC->Fill(dphi, wgt);
+               } // tracks
+           } // tracks 
+       } // event loop
+    } // event loop
+    
+//    Correlated 
+    AliAODEvent* aod = GetEvent(0);
+    Int_t ntracks = aod->GetNumberOfTracks();
+    for (Int_t iTracks = 0; iTracks < ntracks; iTracks++) {
+       for (Int_t jTracks = (iTracks+1); jTracks < ntracks; jTracks++) {
+             AliAODTrack* track1 = aod->GetTrack(iTracks);
+             AliAODTrack* track2 = aod->GetTrack(jTracks);
+             Float_t phi1 = track1->Phi();
+             Float_t phi2 = track2->Phi();
+             Float_t dphi = TMath::Abs(phi1 - phi2);
+             if (dphi > TMath::Pi()) dphi = 2. * TMath::Pi() - dphi;
+             fHistDphiCO->Fill(dphi);
+       } // tracks
+    }
+    
+  // Post output data.
+  PostData(1, fHists);
+}      
+
+//________________________________________________________________________
+void AliAnalysisTaskPhiCorr::Terminate(Option_t *) 
+{
+  // Draw result to the screen
+  // Called once at the end of the query
+
+  TCanvas *c1 = new TCanvas("AliAnalysisTaskPt","Pt",10,10,510,510);
+  c1->cd(1)->SetLogy();
+  fHistDphiUC->DrawCopy("E");
+  fHistDphiCO->DrawCopy("Esame");
+}
diff --git a/ANALYSIS/AliAnalysisTaskPhiCorr.h b/ANALYSIS/AliAnalysisTaskPhiCorr.h
new file mode 100644 (file)
index 0000000..c95f5b2
--- /dev/null
@@ -0,0 +1,34 @@
+// gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include");\r
+#ifndef AliAnalysisTaskPt_cxx\r
+#define AliAnalysisTaskPt_cxx\r
+\r
+// example of an analysis task creating a p_t spectrum\r
+// Authors: Panos Cristakoglou, Jan Fiete Grosse-Oetringhaus, Christian Klein-Boesing\r
+\r
+class TH1F;\r
+class TList;\r
+class AliESDEvent;\r
+\r
+#include "AliAnalysisTaskME.h"\r
+\r
+class AliAnalysisTaskPhiCorr : public AliAnalysisTaskME {\r
+ public:\r
+  AliAnalysisTaskPhiCorr(const char *name = "AliAnalysisTaskPt");\r
+  virtual ~AliAnalysisTaskPhiCorr() {}\r
+  \r
+  virtual void   UserCreateOutputObjects();\r
+  virtual void   UserExec(Option_t *option);\r
+  virtual void   Terminate(Option_t *);\r
+  \r
+ private:\r
+  TList       *fHists;      // List with histos\r
+  TH1F        *fHistDphiCO; // Pt spectrum\r
+  TH1F        *fHistDphiUC; // Pt spectrum\r
+   \r
+  AliAnalysisTaskPhiCorr(const AliAnalysisTaskPhiCorr&); // not implemented\r
+  AliAnalysisTaskPhiCorr& operator=(const AliAnalysisTaskPhiCorr&); // not implemented\r
+  \r
+  ClassDef(AliAnalysisTaskPhiCorr, 1); // example of analysis\r
+};\r
+\r
+#endif\r