]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/KINK/AliAnalysisTaskKinkResonance.cxx
new AddTask macros for resonance analysis (P.Ganoti)
[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
92adf4f6 17// products.
f9afc48d 18//-----------------------------------------------------------------------------------------------------------------
f9afc48d 19
92adf4f6 20#include "AliESDEvent.h"
be1a7181 21
f9afc48d 22#include "AliAnalysisTaskKinkResonance.h"
be1a7181 23#include "AliResonanceKink.h"
f9afc48d 24
25ClassImp(AliAnalysisTaskKinkResonance)
26
f9afc48d 27//______________________________________________________________________________
92adf4f6 28AliAnalysisTaskKinkResonance::AliAnalysisTaskKinkResonance(const char *dname)
29 : AliAnalysisTaskSE(dname), fList(0), fKinkResonance(0)
f9afc48d 30
31{
32 // Constructor
33
34 // Define input and output slots here
35 // Input slot #0 works with a TChain
92adf4f6 36 //DefineInput(0, TChain::Class());
f9afc48d 37 DefineOutput(1, TList::Class());
38}
39
40//________________________________________________________________________
92adf4f6 41void AliAnalysisTaskKinkResonance::UserCreateOutputObjects()
f9afc48d 42{
43 // Create histograms
44 // Called once
45
f9afc48d 46 fList= fKinkResonance->GetHistogramList();
47
48}
49
50//________________________________________________________________________
92adf4f6 51void AliAnalysisTaskKinkResonance::UserExec(Option_t *)
f9afc48d 52{
53 // Main loop
54 // Called for each event
55
92adf4f6 56 AliVEvent *event = InputEvent();
57 if (!event) {
58 Printf("ERROR: Could not retrieve event");
59 return;
60 }
61
62 AliESDEvent* esd = dynamic_cast<AliESDEvent*>(event);
63 if (!esd) {
64 Printf("ERROR: Could not retrieve esd");
65 return;
66 }
67
68 AliMCEvent* mcEvent = MCEvent();
69 if (!mcEvent) {
70 Printf("ERROR: Could not retrieve MC event");
71 return;
72 }
73
0aa3a5b9 74 fKinkResonance->SetDebugLevel(fDebug);
75 fKinkResonance->Analyse(esd, mcEvent);
f9afc48d 76
0aa3a5b9 77 PostData(1, fList);
f9afc48d 78
79}
80
81//________________________________________________________________________
82void AliAnalysisTaskKinkResonance::Terminate(Option_t *)
83{
84 // Draw result to the screen
85 // Called once at the end of the query
86 fList = dynamic_cast<TList*> (GetOutputData(1));
87
88 if (!fList) {
89 Printf("ERROR: fList not available");
90 return;
91 }
92 //TH1D* h=(TH1D *) fList->At(15);
93 //TCanvas *c1 = new TCanvas("AliAnalysisTaskKinkResonance","Pt MC",10,10,510,510);
94 //c1->cd(1);
95 //if (h) h->DrawCopy("E");
96
97}
98