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