]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AliAnalysisTaskProtons.cxx
Create the rec-point branch even in the case of no digits. Please review and fix...
[u/mrichter/AliRoot.git] / PWG2 / AliAnalysisTaskProtons.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TList.h"
4 #include "TH2F.h"
5 #include "TF1.h"
6 #include "TCanvas.h"
7 #include "TLorentzVector.h"
8
9 #include "AliAnalysisTask.h"
10 #include "AliAnalysisManager.h"
11
12 #include "AliESDEvent.h"
13 #include "AliESDInputHandler.h"
14 #include "AliAODEvent.h"
15 #include "AliAODInputHandler.h"
16
17 #include "PWG2spectra/SPECTRA/AliProtonAnalysis.h"
18 #include "AliAnalysisTaskProtons.h"
19
20 // Analysis task creating a the 2d y-p_t spectrum of p and antip
21 // Author: Panos Cristakoglou
22
23 ClassImp(AliAnalysisTaskProtons)
24
25 //________________________________________________________________________
26 AliAnalysisTaskProtons::AliAnalysisTaskProtons(const char *name) 
27 : AliAnalysisTask(name, ""), fESD(0), fAOD(0), fAnalysisType("ESD"), 
28   fList(0), fAnalysis(0), 
29   fElectronFunction(0), fMuonFunction(0),
30   fPionFunction(0), fKaonFunction(0), fProtonFunction(0),
31   fFunctionUsed(kFALSE) { 
32   // Constructor
33
34   // Define input and output slots here
35   // Input slot #0 works with a TChain
36   DefineInput(0, TChain::Class());
37   // Output slot #0 writes into a TList container
38   DefineOutput(0, TList::Class());
39 }
40
41 //________________________________________________________________________
42 void AliAnalysisTaskProtons::ConnectInputData(Option_t *) {
43   // Connect ESD or AOD here
44   // Called once
45
46   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
47   if (!tree) {
48     Printf("ERROR: Could not read chain from input slot 0");
49   } else {
50     // Disable all branches and enable only the needed ones
51     // The next two lines are different when data produced as AliESDEvent is read
52     if(fAnalysisType == "ESD") {
53       tree->SetBranchStatus("*", kFALSE);
54       tree->SetBranchStatus("fTracks.*", kTRUE);
55
56       AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
57       
58       if (!esdH) {
59         Printf("ERROR: Could not get ESDInputHandler");
60       } else
61         fESD = esdH->GetEvent();
62     }
63     else if(fAnalysisType == "AOD") {
64       AliAODInputHandler *aodH = dynamic_cast<AliAODInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
65       
66       if (!aodH) {
67         Printf("ERROR: Could not get AODInputHandler");
68       } else
69         fAOD = aodH->GetEvent();
70     }
71     else 
72       Printf("Wrong analysis type: Only ESD and AOD types are allowed!");
73   }
74 }
75
76 //________________________________________________________________________
77 void AliAnalysisTaskProtons::CreateOutputObjects() {
78   // Create histograms
79   // Called once
80   Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
81
82   fAnalysis = new AliProtonAnalysis();
83   fAnalysis->InitHistograms(10,-1.0,1.0,30,0.1,3.1);
84   if(fAnalysisType == "ESD") {
85     fAnalysis->SetMinTPCClusters(50);
86     fAnalysis->SetMinITSClusters(1);
87     fAnalysis->SetMaxChi2PerTPCCluster(3.5);
88     fAnalysis->SetMaxCov11(2.0);
89     fAnalysis->SetMaxCov22(2.0);
90     fAnalysis->SetMaxCov33(0.5);
91     fAnalysis->SetMaxCov44(0.5);
92     fAnalysis->SetMaxCov55(2.0);
93     fAnalysis->SetMaxSigmaToVertex(3.);
94     fAnalysis->SetITSRefit();
95     fAnalysis->SetTPCRefit();
96   }
97   if(fFunctionUsed)
98     fAnalysis->SetPriorProbabilityFunctions(fElectronFunction,
99                                             fMuonFunction,
100                                             fPionFunction,
101                                             fKaonFunction,
102                                             fProtonFunction);
103   else
104     fAnalysis->SetPriorProbabilities(partFrac);
105
106   fList = new TList();
107   fList->Add(fAnalysis->GetProtonYPtHistogram()); 
108   fList->Add(fAnalysis->GetAntiProtonYPtHistogram()); 
109 }
110
111 //________________________________________________________________________
112 void AliAnalysisTaskProtons::Exec(Option_t *) {
113   // Main loop
114   // Called for each event
115   
116   if(fAnalysisType == "ESD") {
117     if (!fESD) {
118       Printf("ERROR: fESD not available");
119       return;
120     }
121   
122     Printf("Proton analysis task: There are %d tracks in this event", fESD->GetNumberOfTracks());
123     fAnalysis->Analyze(fESD);
124   }
125   else if(fAnalysisType == "AOD") {
126     if (!fAOD) {
127       Printf("ERROR: fAOD not available");
128       return;
129     }
130     
131     Printf("Proton analysis task: There are %d tracks in this event", fAOD->GetNumberOfTracks());
132     fAnalysis->Analyze(fAOD);
133   }
134
135   // Post output data.
136   PostData(0, fList);
137 }      
138
139 //________________________________________________________________________
140 void AliAnalysisTaskProtons::Terminate(Option_t *) {
141   // Draw result to the screen
142   // Called once at the end of the query
143   
144   fList = dynamic_cast<TList*> (GetOutputData(0));
145   if (!fList) {
146     Printf("ERROR: fList not available");
147     return;
148   }
149    
150   TH2F *fHistYPtProtons = (TH2F *)fList->At(0);
151   TH2F *fHistYPtAntiProtons = (TH2F *)fList->At(1);
152     
153   TCanvas *c2 = new TCanvas("c2","p-\bar{p}",200,0,800,400);
154   c2->SetFillColor(10); c2->SetHighLightColor(10); c2->Divide(2,1);
155
156   c2->cd(1)->SetLeftMargin(0.15); c2->cd(1)->SetBottomMargin(0.15);  
157   if (fHistYPtProtons) fHistYPtProtons->DrawCopy("colz");
158   c2->cd(2)->SetLeftMargin(0.15); c2->cd(2)->SetBottomMargin(0.15);  
159   if (fHistYPtAntiProtons) fHistYPtAntiProtons->DrawCopy("colz");
160 }
161
162