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