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