]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ESDCheck/AliMUONQATask.cxx
New method to clone current raw-data event and create a single-event raw-reader....
[u/mrichter/AliRoot.git] / ESDCheck / AliMUONQATask.cxx
... / ...
CommitLineData
1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
18/// An analysis task to check the MUON data in simulated data
19/// This class checks out the ESD tree, providing the matching with
20/// the trigger,trigger responses for Low and High Pt cuts
21/// (in Single, Unlike Sign and Like Sign) and gives Pt, Y, ITS vertex
22/// and multiplicity distributions. All results are in histogram form.
23/// The output is a root file and eps files in MUON.tar.gz.
24
25//*-- Frederic Yermia, yermia@to.infn.it
26//////////////////////////////////////////////////////////////////////////////
27//////////////////////////////////////////////////////////////////////////////
28
29#include <TCanvas.h>
30#include <TChain.h>
31#include <TFile.h>
32#include <TH1F.h>
33#include <TROOT.h>
34#include <TLorentzVector.h>
35#include <TString.h>
36
37#include "AliMUONQATask.h"
38#include "AliESD.h"
39#include "AliLog.h"
40#include "AliESDVertex.h"
41#include "AliESDMuonTrack.h"
42#include <TLorentzVector.h>
43//______________________________________________________________________________
44AliMUONQATask::AliMUONQATask(const char *name) :
45 AliAnalysisTask(name,""),
46 fChain(0),
47 fESD(0),
48 fOutputContainer(0),
49 fV1(),
50 fnTrackTrig(0),
51 ftracktot(0),
52 fnevents(0),
53 fSLowpt(0),
54 fUSLowpt(0),
55 fUSHighpt(0),
56 fLSLowpt(0),
57 fLSHighpt(0),
58 fmuonMass(0.105658389),
59 fthetaX(0),
60 fthetaY(0),
61 fpYZ(0),
62 fPxRec1(0),
63 fPyRec1(0),
64 fPzRec1(0),
65 fE1(0),
66 fZ1(0),
67 fhMUONVertex(0),
68 fhMUONMult(0),
69 fhPt(0),
70 fhY(0),
71 fheffMatchT(0),
72 fhSLowpt(0),
73 fhUSLowpt(0),
74 fhUSHighpt(0),
75 fhLSLowpt(0),
76 fhLSHighpt(0),
77 fhChi2(0),
78 fhChi2match(0)
79{
80 // Constructor.
81 // Input slot #0 works with an Ntuple
82 DefineInput(0, TChain::Class());
83 // Output slot #0 writes into a TH1 container
84 DefineOutput(0, TObjArray::Class()) ;
85
86}
87
88
89//______________________________________________________________________________
90AliMUONQATask::~AliMUONQATask()
91{
92 // dtor
93 fOutputContainer->Clear() ;
94 delete fOutputContainer ;
95
96 delete fhMUONVertex ;
97 delete fhMUONMult ;
98 delete fhPt ;
99 delete fhY ;
100 delete fheffMatchT ;
101 delete fhSLowpt ;
102 delete fhUSLowpt ;
103 delete fhUSHighpt;
104 delete fhLSLowpt ;
105 delete fhLSHighpt;
106 delete fhChi2 ;
107 delete fhChi2match ;
108}
109
110//______________________________________________________________________________
111void AliMUONQATask::ConnectInputData(const Option_t*)
112{
113 // Initialisation of branch container and histograms
114
115 AliInfo(Form("*** Initialization of %s", GetName())) ;
116
117 // Get input data
118 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
119 if (!fChain) {
120 AliError(Form("Input 0 for %s not found\n", GetName()));
121 return ;
122 }
123
124 // One should first check if the branch address was taken by some other task
125 char ** address = (char **)GetBranchAddress(0, "ESD");
126 if (address) {
127 fESD = (AliESD*)(*address);
128 } else {
129 fESD = new AliESD();
130 SetBranchAddress(0, "ESD", &fESD);
131 }
132}
133
134//________________________________________________________________________
135void AliMUONQATask::CreateOutputObjects()
136{
137 // create histograms
138
139 OpenFile(0) ;
140
141 fhMUONVertex = new TH1F("hMUONVertex","ITS Vertex" ,100, -25., 25.);
142 fhMUONMult = new TH1F("hMUONMult" ,"Multiplicity of ESD tracks",10, -0.5, 9.5);
143 fhPt = new TH1F("hPt","Pt",100, 0.,20.);
144 fhY = new TH1F("hY","Rapidity",100,-5.,-1.);
145 fheffMatchT = new TH1F("heff_matchT","Trigger Matching Efficiency",100, 0.,100.);
146 fhSLowpt = new TH1F("hSLowpt","Single Low Pt Response (%)",101, 0.,101.);
147 fhUSLowpt = new TH1F("hUSLowpt","Unlike Sign Low Pt Response (%)",101, 0.,101.);
148 fhUSHighpt = new TH1F("hUSHighpt","Unlike Sign High Pt Response (%)",101, 0.,101.);
149 fhLSLowpt = new TH1F("hLSLowpt","Like Sign Low Pt Response (%)",101, 0.,101.);
150 fhLSHighpt = new TH1F("hLSHighpt","Like Sign High Pt Response (%)",101, 0.,101.);
151 fhChi2 = new TH1F("hChi2","Chi2 by d.o.f.",100, 0.,20.);
152 fhChi2match = new TH1F("hChi2match","Chi2 of trig/track matching",100, 0.,20.);
153 // create output container
154
155 fOutputContainer = new TObjArray(12) ;
156 fOutputContainer->SetName(GetName()) ;
157 fOutputContainer->AddAt(fhMUONVertex, 0) ;
158 fOutputContainer->AddAt(fhMUONMult, 1) ;
159 fOutputContainer->AddAt(fhPt, 2) ;
160 fOutputContainer->AddAt(fhY, 3) ;
161 fOutputContainer->AddAt(fheffMatchT, 4) ;
162 fOutputContainer->AddAt(fhSLowpt, 5) ;
163 fOutputContainer->AddAt(fhUSLowpt, 6) ;
164 fOutputContainer->AddAt(fhUSHighpt, 7) ;
165 fOutputContainer->AddAt(fhLSLowpt, 8) ;
166 fOutputContainer->AddAt(fhLSHighpt, 9) ;
167 fOutputContainer->AddAt(fhChi2, 10) ;
168 fOutputContainer->AddAt( fhChi2match, 11) ;
169}
170
171//______________________________________________________________________________
172void AliMUONQATask::Exec(Option_t *)
173{
174 // Processing of one event
175
176 fnevents++ ;
177
178 Long64_t entry = fChain->GetReadEntry() ;
179
180 if (!fESD) {
181 AliError("fESD is not connected to the input!") ;
182 return ;
183 }
184
185 if ( !((entry-1)%100) )
186 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
187
188 // ************************ MUON *************************************
189
190 const AliESDVertex* vertex = dynamic_cast<const AliESDVertex*>(fESD->GetVertex()) ;
191
192 Double_t zVertex = 0. ;
193 if (vertex)
194 zVertex = vertex->GetZv() ;
195
196 Int_t nTracks = fESD->GetNumberOfMuonTracks() ;
197
198 ULong64_t trigWord = fESD->GetTriggerMask() ;
199
200 if (trigWord & 0x80) {
201 fSLowpt++;
202 }
203 if (trigWord & 0x100){
204 fLSLowpt++;
205 }
206 if (trigWord & 0x200){
207 fLSHighpt++;
208 }
209 if (trigWord & 0x400){
210 fUSLowpt++;
211 }
212 if (trigWord & 0x800){
213 fUSHighpt++;
214 }
215
216 Int_t tracktrig = 0 ;
217 Int_t iTrack1 ;
218
219 for (iTrack1 = 0 ; iTrack1 < nTracks ; iTrack1++) { //1st loop
220 AliESDMuonTrack* muonTrack = fESD->GetMuonTrack(iTrack1) ;
221 ftracktot++ ;
222 fthetaX = muonTrack->GetThetaX();
223 fthetaY = muonTrack->GetThetaY();
224 fpYZ = 1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
225 fPzRec1 = - fpYZ / TMath::Sqrt(1.0 + TMath::Tan(fthetaY)*TMath::Tan(fthetaY));
226 fPxRec1 = fPzRec1 * TMath::Tan(fthetaX);
227 fPyRec1 = fPzRec1 * TMath::Tan(fthetaY);
228 fZ1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
229 fE1 = TMath::Sqrt(fmuonMass * fmuonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
230 fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
231
232 // -----------> transverse momentum
233 Float_t pt1 = fV1.Pt();
234 // ----------->Rapidity
235 Float_t y1 = fV1.Rapidity();
236
237 if(muonTrack->GetMatchTrigger()) {
238 fnTrackTrig++ ;
239 tracktrig++ ;
240 Float_t Chi2match = muonTrack->GetChi2MatchTrigger();
241 fhChi2match->Fill(Chi2match);
242 }
243
244 Float_t fitfmin = muonTrack->GetChi2();
245 Int_t ntrackhits = muonTrack->GetNHit();
246 Float_t Chi2= fitfmin / (2.0 * ntrackhits - 5);
247
248 fhChi2->Fill(Chi2);
249 fhPt->Fill(pt1);
250 fhY->Fill(y1);
251 }
252
253 fhMUONVertex->Fill(zVertex) ;
254 fhMUONMult->Fill(Float_t(nTracks)) ;
255
256 PostData(0, fOutputContainer);
257}
258
259//______________________________________________________________________________
260void AliMUONQATask::Terminate(Option_t *)
261{
262 // Processing when the event loop is ended
263 Int_t fSLowPt = fSLowpt;
264 if(fnevents){
265 fSLowPt = 100 * fSLowpt / fnevents ;
266 fhSLowpt->Fill(fSLowPt); }
267 Int_t fUSLowPt = fUSLowpt;
268 if(fnevents){
269 fUSLowPt = 100 * fUSLowpt / fnevents ;
270 fhUSLowpt->Fill(fUSLowPt); }
271 Int_t fUSHighPt = fUSHighpt;
272 if(fnevents){
273 fUSHighPt = 100 * fUSHighpt / fnevents ;
274 fhUSHighpt->Fill(fUSHighPt); }
275 Int_t fLSLowPt = fLSLowpt;
276 if(fnevents){
277 fLSLowPt = 100 * fLSLowpt / fnevents ;
278 fhLSLowpt->Fill(fLSLowPt); }
279 Int_t fLSHighPt = fLSHighpt;
280 if(fnevents){
281 fLSHighPt = 100 * fLSHighpt / fnevents ;
282 fhLSHighpt->Fill(fLSHighPt); }
283
284 Int_t effMatch = -1 ;
285 if (ftracktot){
286 effMatch = 100 * fnTrackTrig / ftracktot ;
287 fheffMatchT->Fill(effMatch);}
288
289 Bool_t problem = kFALSE ;
290 AliInfo(Form(" *** %s Report:", GetName())) ;
291
292 fOutputContainer = (TObjArray*)GetOutputData(0);
293 fhMUONVertex = (TH1F*)fOutputContainer->At(0);
294 fhMUONMult = (TH1F*)fOutputContainer->At(1);
295 fhPt = (TH1F*)fOutputContainer->At(2);
296 fhY = (TH1F*)fOutputContainer->At(3);
297 fheffMatchT=(TH1F*)fOutputContainer->At(4);
298 fhSLowpt=(TH1F*)fOutputContainer->At(5);
299 fhUSLowpt=(TH1F*)fOutputContainer->At(6);
300 fhUSHighpt=(TH1F*)fOutputContainer->At(7);
301 fhLSLowpt=(TH1F*)fOutputContainer->At(8);
302 fhLSHighpt=(TH1F*)fOutputContainer->At(9);
303
304 printf(" Total number of processed events %d \n", fnevents) ;
305 printf(" \n") ;
306 printf(" \n") ;
307 printf(" Table 1: \n") ;
308 printf(" ===================================================\n") ;
309 printf(" Global Trigger output Low pt High pt \n") ;
310 printf(" number of Single :\t");
311 printf(" %i\t", fSLowpt) ;
312 printf("\n");
313 printf(" number of UnlikeSign pair :\t");
314 printf(" %i\t%i\t", fUSLowpt, fUSHighpt) ;
315 printf("\n");
316 printf(" number of LikeSign pair :\t");
317 printf(" %i\t%i\t", fLSLowpt, fLSHighpt) ;
318 printf("\n");
319 printf(" matching efficiency with the trigger for single tracks = %2d %% \n", effMatch);
320 printf("\n") ;
321
322 TCanvas * cMUON1 = new TCanvas("cMUON1", "MUON ESD Vert & Mult", 400, 10, 600, 700) ;
323 cMUON1->Divide(1,2) ;
324 cMUON1->cd(1) ;
325 fhMUONVertex->SetXTitle("Vz (cm)");
326 fhMUONVertex->Draw() ;
327 cMUON1->cd(2) ;
328 fhMUONMult->SetXTitle(" Track Multiplicity");
329 fhMUONMult->Draw() ;
330 cMUON1->Print("MUON1.eps") ;
331
332 TCanvas * cMUON2 = new TCanvas("cMUON2", "MUON ESD Pt & Y", 400, 10, 600, 700) ;
333 cMUON2->Divide(1,2) ;
334 cMUON2->cd(1) ;
335 fhPt->SetXTitle("Pt (GeV)");
336 fhPt->Draw() ;
337 cMUON2->cd(2) ;
338 fhY->SetXTitle("Y");
339 fhY->Draw() ;
340 cMUON2->Print("MUON2.eps") ;
341
342 TCanvas * cMUON3 = new TCanvas("cMUON3", "Track Chi2 by dof and Chi2 of trigger/track matching ", 400, 10, 600, 700) ;
343 cMUON3->Divide(1,2) ;
344 cMUON3->cd(1) ;
345 fhChi2->SetXTitle("Chi2 by d.o.f.");
346 fhChi2->Draw();
347 cMUON3->cd(2) ;
348 fhChi2match->SetXTitle("Chi2 of trig/track matching");
349 fhChi2match->Draw();
350 cMUON3->Print("MUON3.eps") ;
351
352 TCanvas * cMUON4 = new TCanvas("cMUON4", "Trigger Matching and Trigger Response (%)", 400, 10, 600, 700) ;
353 cMUON4->Divide(2,3) ;
354 cMUON4->cd(1) ;
355 fheffMatchT->SetXTitle("%");
356 fheffMatchT->Draw() ;
357 cMUON4->cd(2) ;
358 fhSLowpt->SetXTitle("%");
359 fhSLowpt->Draw() ;
360 cMUON4->cd(3) ;
361 fhUSLowpt->SetXTitle("%");
362 fhUSLowpt->Draw() ;
363 cMUON4->cd(4) ;
364 fhUSHighpt->SetXTitle("%");
365 fhUSHighpt->Draw() ;
366 cMUON4->cd(5) ;
367 fhLSLowpt->SetXTitle("%");
368 fhLSLowpt->Draw() ;
369 cMUON4->cd(6) ;
370 fhLSHighpt->SetXTitle("%");
371 fhLSHighpt->Draw() ;
372 cMUON4->Print("MUON4.eps") ;
373
374 char line[1024] ;
375 sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ;
376 gROOT->ProcessLine(line);
377 sprintf(line, ".!rm -fR *.eps");
378 gROOT->ProcessLine(line);
379
380 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
381
382 TString report ;
383 if(problem)
384 report="Problems found, please check!!!";
385 else
386 report="OK";
387
388 AliInfo(Form("*** %s Summary Report: %s \n",GetName(), report.Data())) ;
389}