]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/muon/AliAnalysisTaskLUT.cxx
unused variable removed
[u/mrichter/AliRoot.git] / PWG3 / muon / AliAnalysisTaskLUT.cxx
CommitLineData
27de2dfb 1
2/* $Id$ */
3
39280342 4//===================================================================
5// Class AliAnalysisTaskLUT
6//
7// Extract ESD muon tracks information and store in ntuple.
8// Used for offline calculation/adjustment of Look-up-Tables for the
9// trigger chambers.
10//
11// Clermont-Ferrand 2008
12//===================================================================
13
14#define AliAnalysisTaskLUT_cxx
15
16#include "Riostream.h"
17#include "TChain.h"
18#include "TTree.h"
19#include "TNtuple.h"
20
21#include "AliESDEvent.h"
22#include "AliESDVertex.h"
23#include "AliESDMuonTrack.h"
24
25#include "AliAnalysisTask.h"
26#include "AliAnalysisTaskLUT.h"
27
28ClassImp(AliAnalysisTaskLUT)
29
30//________________________________________________________________________
31AliAnalysisTaskLUT::AliAnalysisTaskLUT(const char *name) :
32 AliAnalysisTask(name,""),
33 fESDEvent(0),
34 fTracksLUT(0)
35{
36 // Constructor.
37 // Input slot #0 works with a TChain
38 DefineInput(0, TChain::Class());
39
40 // Output slot #0 writes into a TNtuple container
41 DefineOutput(0, TNtuple::Class());
42
43}
44
45//___________________________________________________________________________
46void AliAnalysisTaskLUT::ConnectInputData(Option_t *)
47{
48 // Input ESD tree
49 printf(" ConnectInputData of task %s\n", GetName());
50
51 TTree *tinput = (TTree*)GetInputData(0);
52 fESDEvent = new AliESDEvent();
53 fESDEvent->ReadFromTree(tinput);
54
55}
56
57//___________________________________________________________________________
58void AliAnalysisTaskLUT::CreateOutputObjects()
59{
60 // Output object (TNtuple)
61 printf(" CreateOutputObjects of task %s\n", GetName());
62
63 fTracksLUT = new TNtuple("ntTracksLUT","ntTracksLUT","TrigMask:VertZ:NHit:Circ:StripX:StripY:Dev:Lpt:Hpt:Chi2:Chi2trg:Z:P:Pt:Phi:The:TrackZ:PU:PtU:PhiU:TheU:TrackZU");
64
65}
66
67//________________________________________________________________________
68void AliAnalysisTaskLUT::Exec(Option_t *)
69{
70 // Execution
71 Float_t ntvar[22];
72
73 Float_t px, py, pz, pt;
74 Float_t zVertex;
75 Int_t dev;
76
77 if (!fESDEvent) return;
78
79 Int_t nTracks = fESDEvent->GetNumberOfMuonTracks();
80
81 AliESDVertex* vertex = (AliESDVertex*) fESDEvent->GetVertex();
82 zVertex = 0.0;
83 if (vertex) {
84 zVertex = vertex->GetZv();
85 }
86 ULong64_t mask = fESDEvent->GetTriggerMask();
87
88 ntvar[ 0] = mask;
89 ntvar[ 1] = zVertex;
90
91 for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
92
93 AliESDMuonTrack* muonTrack = fESDEvent->GetMuonTrack(iTracks);
94 if (muonTrack == 0x0) continue;
95
96 if (!muonTrack->GetMatchTrigger()) continue;
97
98 dev = muonTrack->LoDev();
99 dev *= (Int_t)(TMath::Sign(1.0,muonTrack->GetInverseBendingMomentum()));
100 dev += 15;
101
102 ntvar[ 2] = muonTrack->GetNHit();
103 ntvar[ 3] = muonTrack->LoCircuit();
104 ntvar[ 4] = muonTrack->LoStripX();
105 ntvar[ 5] = muonTrack->LoStripY();
106 ntvar[ 6] = dev;
107 ntvar[ 7] = muonTrack->LoLpt();
108 ntvar[ 8] = muonTrack->LoHpt();
109 ntvar[ 9] = muonTrack->GetChi2() / (2.0 * muonTrack->GetNHit() - 5);
110 ntvar[10] = muonTrack->GetChi2MatchTrigger();
111 ntvar[11] = TMath::Sign(1.,muonTrack->GetInverseBendingMomentum());
112
113 // corrected to vertex
114
115 px = muonTrack->Px();
116 py = muonTrack->Py();
117 pz = muonTrack->Pz();
118 pt = TMath::Sqrt(px*px+py*py);
119 ntvar[12] = muonTrack->P();
120 ntvar[13] = TMath::Sqrt(px*px+py*py);
121 ntvar[14] = TMath::ATan2(py,px)*180./TMath::Pi();
122 ntvar[15] = TMath::ATan2(pt,pz)*180./TMath::Pi();
123 ntvar[16] = muonTrack->GetZ();
124
125 // uncorrected
126
127 px = muonTrack->PxUncorrected();
128 py = muonTrack->PyUncorrected();
129 pz = muonTrack->PzUncorrected();
130 pt = TMath::Sqrt(px*px+py*py);
131 ntvar[17] = muonTrack->PUncorrected();
132 ntvar[18] = TMath::Sqrt(px*px+py*py);
133 ntvar[19] = TMath::ATan2(py,px)*180./TMath::Pi();
134 ntvar[20] = TMath::ATan2(pt,pz)*180./TMath::Pi();
135 ntvar[21] = muonTrack->GetZUncorrected();
136
137 fTracksLUT->Fill(ntvar);
138
139 } // end ESD track loop
140
141 // Post final data. It will be written to a file with option "RECREATE"
142 PostData(0, fTracksLUT);
143
144}
145
146//________________________________________________________________________
147void AliAnalysisTaskLUT::Terminate(Option_t *) {
148 // End function, empty
149}
150