]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/AliAnalysisTaskProtons.cxx
Removing not used parameter and putting the correct flag for the ITS clusters
[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 "TCanvas.h"
6 #include "TLorentzVector.h"
7
8 #include "AliAnalysisTask.h"
9 #include "AliAnalysisManager.h"
10
11 #include "AliESDEvent.h"
12 #include "AliESDInputHandler.h"
13
14 #include "PWG2spectra/SPECTRA/AliProtonAnalysis.h"
15 #include "AliAnalysisTaskProtons.h"
16
17 // Analysis task creating a the 2d y-p_t spectrum of p and antip
18 // Author: Panos Cristakoglou
19
20 ClassImp(AliAnalysisTaskProtons)
21
22 //________________________________________________________________________
23 AliAnalysisTaskProtons::AliAnalysisTaskProtons(const char *name) 
24   : AliAnalysisTask(name, ""), fESD(0), fList(0), fAnalysis(0) { 
25   // Constructor
26
27   // Define input and output slots here
28   // Input slot #0 works with a TChain
29   DefineInput(0, TChain::Class());
30   // Output slot #0 writes into a TList container
31   DefineOutput(0, TList::Class());
32 }
33
34 //________________________________________________________________________
35 void AliAnalysisTaskProtons::ConnectInputData(Option_t *) {
36   // Connect ESD or AOD here
37   // Called once
38
39   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
40   if (!tree) {
41     Printf("ERROR: Could not read chain from input slot 0");
42   } else {
43     // Disable all branches and enable only the needed ones
44     // The next two lines are different when data produced as AliESDEvent is read
45     tree->SetBranchStatus("*", kFALSE);
46     tree->SetBranchStatus("fTracks.*", kTRUE);
47
48     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
49
50     if (!esdH) {
51       Printf("ERROR: Could not get ESDInputHandler");
52     } else
53       fESD = esdH->GetEvent();
54   }
55 }
56
57 //________________________________________________________________________
58 void AliAnalysisTaskProtons::CreateOutputObjects() {
59   // Create histograms
60   // Called once
61   Double_t partFrac[5] = {0.01, 0.01, 0.85, 0.10, 0.05};
62
63   fAnalysis = new AliProtonAnalysis();
64   fAnalysis->InitHistograms(10,-1.0,1.0,30,0.1,3.1);
65   fAnalysis->SetMinTPCClusters(50);
66   fAnalysis->SetMinITSClusters(1);
67   fAnalysis->SetMaxChi2PerTPCCluster(3.5);
68   fAnalysis->SetMaxCov11(2.0);
69   fAnalysis->SetMaxCov22(2.0);
70   fAnalysis->SetMaxCov33(0.5);
71   fAnalysis->SetMaxCov44(0.5);
72   fAnalysis->SetMaxCov55(2.0);
73   fAnalysis->SetMaxSigmaToVertex(3.);
74   fAnalysis->SetITSRefit();
75   fAnalysis->SetTPCRefit();
76   fAnalysis->SetPriorProbabilities(partFrac);
77 }
78
79 //________________________________________________________________________
80 void AliAnalysisTaskProtons::Exec(Option_t *) {
81   // Main loop
82   // Called for each event
83   
84   if (!fESD) {
85     Printf("ERROR: fESD not available");
86     return;
87   }
88   
89   Printf("Task2: There are %d tracks in this event", fESD->GetNumberOfTracks());
90   fAnalysis->Analyze(fESD);
91   fList = new TList();
92   fList->Add(fAnalysis->GetProtonYPtHistogram()); 
93   fList->Add(fAnalysis->GetAntiProtonYPtHistogram()); 
94
95   // Post output data.
96   PostData(0, fList);
97 }      
98
99 //________________________________________________________________________
100 void AliAnalysisTaskProtons::Terminate(Option_t *) {
101   // Draw result to the screen
102   // Called once at the end of the query
103   
104   fList = dynamic_cast<TList*> (GetOutputData(0));
105   if (!fList) {
106     Printf("ERROR: fList not available");
107     return;
108   }
109    
110   TH2F *fHistYPtProtons = (TH2F *)fList->At(0);
111   TH2F *fHistYPtAntiProtons = (TH2F *)fList->At(1);
112     
113   TCanvas *c2 = new TCanvas("c2","p-\bar{p}",200,0,800,400);
114   c2->SetFillColor(10); c2->SetHighLightColor(10); c2->Divide(2,1);
115
116   c2->cd(1)->SetLeftMargin(0.15); c2->cd(1)->SetBottomMargin(0.15);  
117   if (fHistYPtProtons) fHistYPtProtons->DrawCopy("colz");
118   c2->cd(2)->SetLeftMargin(0.15); c2->cd(2)->SetBottomMargin(0.15);  
119   if (fHistYPtAntiProtons) fHistYPtAntiProtons->DrawCopy("colz");
120 }
121
122