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