]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/KINK/AliAnalysisTaskKinkResonance.cxx
Hardcoded module number removed; fixed errors of previous commit.
[u/mrichter/AliRoot.git] / PWG2 / KINK / AliAnalysisTaskKinkResonance.cxx
CommitLineData
f9afc48d 1/**************************************************************************
2 * Author: Paraskevi Ganoti, University of Athens (pganoti@phys.uoa.gr) *
3 * Contributors are mentioned in the code where appropriate. *
4 * *
5 * Permission to use, copy, modify and distribute this software and its *
6 * documentation strictly for non-commercial purposes is hereby granted *
7 * without fee, provided that the above copyright notice appears in all *
8 * copies and that both the copyright notice and this permission notice *
9 * appear in the supporting documentation. The authors make no claims *
10 * about the suitability of this software for any purpose. It is *
11 * provided "as is" without express or implied warranty. *
12 **************************************************************************/
13
14//----------------------------------------------------------------------------------------------------------------
15// class AliAnalysisTaskKinkResonance
16// Example of an analysis task for reconstructing resonances having at least one kaon-kink in their decay
17// products. It provides basic plots as well as plots helping to calculate the corrections.
18// Usage: To analyse a resonance having a kaon in its decay products, one has to modify the integer
19// variables resonancePDG, daughter1 and daughter2 accordingly as well as daughter1Mass and daughter2Mass.
20// Also, depending on the analysis mode (ESD or MC), fAnalysisType in the constructor must also be changed
21//-----------------------------------------------------------------------------------------------------------------
22#include "TCanvas.h"
23#include "TChain.h"
24#include "TTree.h"
25#include "TH2D.h"
26
27#include "AliAnalysisTaskSE.h"
28#include "AliAnalysisManager.h"
29#include "AliESDInputHandler.h"
30#include "AliMCEventHandler.h"
31#include "AliMCEvent.h"
32#include "AliResonanceKink.h"
33#include "AliAnalysisTaskKinkResonance.h"
34
35ClassImp(AliAnalysisTaskKinkResonance)
36
37//_____________________________________________________________________________
38AliAnalysisTaskKinkResonance::AliAnalysisTaskKinkResonance()
39 : AliAnalysisTaskSE(), fESD(0), fmcEventH(0), fList(0), fKinkResonance(0)
40
41{
42 // Constructor
43}
44
45//______________________________________________________________________________
46AliAnalysisTaskKinkResonance::AliAnalysisTaskKinkResonance(const char *name)
47 : AliAnalysisTaskSE(name), fESD(0), fmcEventH(0), fList(0), fKinkResonance(0)
48
49{
50 // Constructor
51
52 // Define input and output slots here
53 // Input slot #0 works with a TChain
54 DefineInput(0, TChain::Class());
55 DefineOutput(1, TList::Class());
56}
57
58//________________________________________________________________________
59void AliAnalysisTaskKinkResonance::ConnectInputData(Option_t *)
60{
61 // Connect ESD or AOD here
62 // Called once
63
64 TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
65 if (!tree) {
66 Printf("ERROR: Could not read chain from input slot 0");
67 } else {
68
69 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
70
71 if (!esdH) {
72 Printf("ERROR: Could not get ESDInputHandler");
73 } else
74 fESD = esdH->GetEvent();
75
76 AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
77
78 if (!eventHandler) {
79 Printf("ERROR: Could not retrieve MC event handler");
80 return;
81 }
82
83 fmcEventH = eventHandler->MCEvent();
84 }
85}
86
87//________________________________________________________________________
88void AliAnalysisTaskKinkResonance::CreateOutputObjects()
89{
90 // Create histograms
91 // Called once
92
93 fList=new TList();
94 fList= fKinkResonance->GetHistogramList();
95
96}
97
98//________________________________________________________________________
99void AliAnalysisTaskKinkResonance::Exec(Option_t *)
100{
101 // Main loop
102 // Called for each event
103
104 fKinkResonance->Analyse(fESD, fmcEventH);
105
106 PostData(1, fList);
107
108}
109
110//________________________________________________________________________
111void AliAnalysisTaskKinkResonance::Terminate(Option_t *)
112{
113 // Draw result to the screen
114 // Called once at the end of the query
115 fList = dynamic_cast<TList*> (GetOutputData(1));
116
117 if (!fList) {
118 Printf("ERROR: fList not available");
119 return;
120 }
121 //TH1D* h=(TH1D *) fList->At(15);
122 //TCanvas *c1 = new TCanvas("AliAnalysisTaskKinkResonance","Pt MC",10,10,510,510);
123 //c1->cd(1);
124 //if (h) h->DrawCopy("E");
125
126}
127