]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONCheck.cxx
Coding conventions
[u/mrichter/AliRoot.git] / MUON / AliMUONCheck.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 //-----------------------------------------------------------------------------
19 /// \class AliMUONCheck
20 ///
21 /// This class check the ESD tree, providing the matching with the trigger
22 /// response and designing useful plots (Pt, Y, ITS vertex, multiplicity).
23 ///  Note that there is a special flag to turn on for pdc06 production. 
24 /// It also checks the TR tree giving hit densities on the two first and 
25 /// last planes of the spectrometer as well as the time of flight on these planes.
26 /// MUONkine() provides event stack and check if the event are generated with 
27 /// at least one muon and two muons (for PDC06).
28 /// DumpDigit() as a replacement of the function from MUONCheck.C macro.
29 ///
30 /// \author Frederic Yermia, INFN Torino
31 //-----------------------------------------------------------------------------
32
33 #include "AliMUONCheck.h"
34 #include "AliMUONConstants.h"
35 #include "AliMUONMCDataInterface.h"
36 #include "AliMUONDataInterface.h"
37 #include "AliMpCDB.h"
38 #include "AliMpSegmentation.h"
39 #include "AliMpVSegmentation.h"
40 #include "AliMpDEManager.h"
41 #include "AliMpCDB.h"
42 #include "AliMUONVDigitStore.h"
43
44 #include "AliRunLoader.h"
45 #include "AliLoader.h"
46 #include "AliStack.h"
47 #include "AliTrackReference.h"
48 #include "AliTracker.h"
49 #include "AliESDEvent.h"
50 #include "AliESDMuonTrack.h"
51 #include "AliESDVertex.h"
52 #include "AliMagF.h"
53 #include "AliLog.h"
54
55 #include <TSystem.h>
56 #include <TCanvas.h>
57 #include <TLorentzVector.h>
58 #include <TFile.h>
59 #include <TH1.h>
60 #include <TParticle.h>
61
62 /// \cond CLASSIMP
63 ClassImp(AliMUONCheck)
64 /// \endcond
65
66 AliMUONCheck::AliMUONCheck(const char* galiceFile, const char* esdFile,Int_t firstEvent, Int_t lastEvent,const char* outDir) 
67 : TObject(),
68 fFileName(galiceFile),
69 fFileNameSim(),
70 fesdFileName(esdFile),
71 fkOutDir(outDir),
72 fFirstEvent(firstEvent),
73 fLastEvent(lastEvent)
74 {
75   /// ctor  
76 }
77
78 //_____________________________________________________________________________
79 AliMUONCheck::AliMUONCheck(const char* galiceFile, const char* galiceFileSim,
80                            const char* esdFile,Int_t firstEvent, Int_t lastEvent,
81                            const char* outDir) 
82 : TObject(),
83 fFileName(galiceFile),
84 fFileNameSim(galiceFileSim),
85 fesdFileName(esdFile),
86 fkOutDir(outDir),
87 fFirstEvent(firstEvent),
88 fLastEvent(lastEvent)
89 {
90   /// ctor
91 }
92
93 //_____________________________________________________________________________
94 AliMUONCheck::~AliMUONCheck()
95 {
96   /// Destructor
97 }
98
99 //_____________________________________________________________________________
100 void
101 AliMUONCheck::CheckESD(Bool_t pdc06TriggerResponse) 
102 {
103   /// Check ESD files
104   
105   // Histograms
106   TH1F * fhMUONVertex ; //! 
107   TH1F * fhMUONMult   ; //!
108   
109   // create histograms 
110   fhMUONVertex = new TH1F("hMUONVertex","ITS Vertex"                ,100, -25., 25.);
111   fhMUONMult   = new TH1F("hMUONMult"  ,"Multiplicity of ESD tracks",10,  -0.5, 9.5);
112   
113   TH1F *hY = new TH1F("hY","Rapidity",100,-5.,-1.);
114   TH1F *hPt = new TH1F("hPt","Pt",100, 0.,20.);
115   
116   // ------------->open the ESD file
117   TFile* esdFile = TFile::Open(fesdFileName.Data());
118   
119   if (!esdFile || !esdFile->IsOpen()) 
120   {
121     AliError(Form("Error opening %s file \n",fesdFileName.Data()));
122     return;
123   }
124   
125   Int_t fSPLowpt=0     ; //!
126   Int_t fSPHighpt=0    ; //!
127   Int_t fSPAllpt=0     ; //!
128   Int_t fSMLowpt=0     ; //!
129   Int_t fSMHighpt =0   ; //!
130   Int_t fSMAllpt=0     ; //!
131   Int_t fSULowpt=0     ; //!
132   Int_t fSUHighpt=0    ; //!
133   Int_t fSUAllpt=0     ; //!
134   Int_t fUSLowpt=0     ; //!
135   Int_t fUSHighpt=0    ; //!
136   Int_t fUSAllpt=0     ; //! 
137   Int_t fLSLowpt=0     ; //!
138   Int_t fLSHighpt=0    ; //! 
139   Int_t fLSAllpt=0     ; //!
140   
141   Int_t fSLowpt=0      ; //!
142   Int_t fSHighpt=0     ; //!
143   
144   Int_t fnTrackTrig=0  ; //!
145   Int_t ftracktot=0    ; //!
146   Int_t effMatch=0     ; //!
147   
148   TLorentzVector fV1;
149   Float_t muonMass = 0.105658389;
150   Double_t thetaX, thetaY, pYZ;
151   Double_t fPxRec1, fPyRec1, fPzRec1, fE1;
152   Int_t fZ1;
153   
154   AliESDEvent* fESD = new AliESDEvent();
155   TTree* tree = (TTree*) esdFile->Get("esdTree");
156   if (!tree) 
157   {
158     Error("CheckESD", "no ESD tree found");
159     AliError(Form("CheckESD", "no ESD tree found"));
160     return ;
161   }
162   fESD->ReadFromTree(tree);
163   
164   Int_t fnevents = tree->GetEntries();
165   Int_t endOfLoop = fLastEvent+1;
166   
167   if ( fLastEvent == -1 ) endOfLoop = fnevents;
168   Int_t ievent=0;
169   Int_t nev=0;
170   
171   for (ievent = fFirstEvent; ievent < endOfLoop; ++ievent ) 
172   {
173     nev++;    
174     tree->GetEvent(ievent);
175     if (!fESD) 
176     {
177       Error("CheckESD", "no ESD object found for event %d", ievent);
178       return ;
179     }
180     AliESDVertex* vertex = (AliESDVertex*) fESD->GetVertex();
181     
182     Double_t zVertex = 0. ;
183     if (vertex) zVertex = vertex->GetZv();
184     
185     Int_t nTracks = (Int_t)fESD->GetNumberOfMuonTracks() ;
186     ULong64_t trigword=fESD->GetTriggerMask();
187     
188     if(pdc06TriggerResponse)
189     {
190       if (trigword & 0x01) {
191         fSPLowpt++;
192       }
193       
194       if (trigword & 0x02){
195         fSPHighpt++;
196       }
197       if (trigword & 0x04){
198         fSPAllpt++;
199       } 
200       if (trigword & 0x08){
201         fSMLowpt++;
202       }  
203       if (trigword & 0x010){
204         fSMHighpt++;
205       }
206       if (trigword & 0x020){
207         fSMAllpt++;
208       } 
209       if (trigword & 0x040){
210         fSULowpt++;
211       }  
212       if (trigword & 0x080){
213         fSUHighpt++;
214       }   
215       if (trigword & 0x100){
216         fSUAllpt++;
217       }  
218       if (trigword & 0x200){
219         fUSLowpt++;
220       }
221       
222       if (trigword & 0x400){
223         fUSHighpt++;
224       }
225       if (trigword & 0x800){
226         fUSAllpt++;
227       }
228       if (trigword & 0x1000){
229         fLSLowpt++;
230       }
231       
232       if (trigword & 0x2000){
233         fLSHighpt++;
234       }
235       
236       if (trigword & 0x4000){
237         fLSAllpt++;
238       }
239     }// if pdc06TriggerResponse
240     else {
241       if (trigword & 0x01) {
242         fSLowpt++;
243       }
244       
245       if (trigword & 0x02){
246         fSHighpt++;
247       }
248       if (trigword & 0x04){
249         fLSLowpt++;
250       } 
251       if (trigword & 0x08){
252         fLSHighpt++;
253       }  
254       if (trigword & 0x010){
255         fUSLowpt++;
256       }
257       if (trigword & 0x020){
258         fUSHighpt++;
259       }
260     }
261     
262     Int_t tracktrig=0;
263     
264     for ( Int_t iTrack1 = 0; iTrack1<nTracks; ++iTrack1 ) 
265     { //1st loop
266       AliESDMuonTrack* muonTrack = fESD->GetMuonTrack(iTrack1);
267       
268       // skip fake tracks (ghosts)
269       if (!muonTrack->ContainTrackerData()) continue;
270       
271       ftracktot++;
272       
273       thetaX = muonTrack->GetThetaX();
274       thetaY = muonTrack->GetThetaY();
275       pYZ     =  1./TMath::Abs(muonTrack->GetInverseBendingMomentum());
276       
277       fPzRec1  = - pYZ / TMath::Sqrt(1.0 + TMath::Tan(thetaY)*TMath::Tan(thetaY));
278       fPxRec1  = fPzRec1 * TMath::Tan(thetaX);
279       fPyRec1  = fPzRec1 * TMath::Tan(thetaY);
280       fZ1 = Int_t(TMath::Sign(1.,muonTrack->GetInverseBendingMomentum()));
281       fE1 = TMath::Sqrt(muonMass * muonMass + fPxRec1 * fPxRec1 + fPyRec1 * fPyRec1 + fPzRec1 * fPzRec1);
282       fV1.SetPxPyPzE(fPxRec1, fPyRec1, fPzRec1, fE1);
283       // -----------> transverse momentum
284       Float_t pt1 = fV1.Pt();
285       // ----------->Rapidity
286       Float_t y1 = fV1.Rapidity();
287       
288       if(muonTrack->GetMatchTrigger()) 
289       {
290         fnTrackTrig++;
291         tracktrig++;
292       }
293       hY->Fill(y1);
294       hPt->Fill(pt1);
295     } // loop on track
296     
297     fhMUONVertex->Fill(zVertex) ;
298     fhMUONMult->Fill(Float_t(nTracks)) ;
299     
300   } // loop over events
301   
302   AliInfo(Form("Terminate %s:", GetName())) ;
303   
304   effMatch=100*fnTrackTrig/ftracktot;
305   
306   if(pdc06TriggerResponse)
307   {
308     printf("=================================================================\n") ;
309     printf("================  %s ESD SUMMARY    ==============\n", GetName()) ;
310     printf("                                                   \n") ;
311     printf("         Total number of processed events  %d      \n", nev) ;
312     printf("\n")  ;
313     printf("\n")  ;
314     printf("Table 1:                                         \n") ;
315     printf(" Global Trigger output       Low pt  High pt   All\n") ;
316     printf(" number of Single Plus      :\t");
317     printf("%i\t%i\t%i\t", fSPLowpt, fSPHighpt, fSPAllpt) ;
318     printf("\n");
319     printf(" number of Single Minus     :\t");
320     printf("%i\t%i\t%i\t", fSMLowpt, fSMHighpt, fSMAllpt) ;
321     printf("\n");
322     printf(" number of Single Undefined :\t"); 
323     printf("%i\t%i\t%i\t", fSULowpt, fSUHighpt, fSUAllpt) ;
324     printf("\n");
325     printf(" number of UnlikeSign pair  :\t"); 
326     printf("%i\t%i\t%i\t", fUSLowpt, fUSHighpt, fUSAllpt) ;
327     printf("\n");
328     printf(" number of LikeSign pair    :\t");  
329     printf("%i\t%i\t%i\t", fLSLowpt, fLSHighpt, fLSAllpt) ;
330     printf("\n");
331     printf("===================================================\n") ;
332     printf("\n") ;
333     printf("matching efficiency with the trigger for single tracks = %2d %% \n", effMatch);
334     printf("================================================================\n") ;  printf("\n") ;
335     
336   }//if(pdc06TriggerResponse)
337   
338   gSystem->cd(fkOutDir);
339   
340   FILE *outtxt=fopen("output.txt","a");
341   freopen("output.txt","a",outtxt);
342   
343   if(pdc06TriggerResponse){
344     fprintf(outtxt,"                                                   \n");
345     fprintf(outtxt,"===================================================\n");
346     fprintf(outtxt,"================      ESD SUMMARY    ==============\n");
347     fprintf(outtxt,"                                                   \n");
348     fprintf(outtxt,"    Total number of processed events  %d      \n", nev); 
349     fprintf(outtxt,"\n");
350     fprintf(outtxt,"\n");
351     fprintf(outtxt,"Table 1:                                         \n");
352     fprintf(outtxt," Global Trigger output       Low pt  High pt   All\n");
353     fprintf(outtxt," number of Single Plus      :\t");
354     fprintf(outtxt,"%i\t%i\t%i\t",fSPLowpt,fSPHighpt,fSPAllpt);
355     fprintf(outtxt,"\n");
356     fprintf(outtxt," number of Single Minus     :\t");
357     fprintf(outtxt,"%i\t%i\t%i\t",fSMLowpt,fSMHighpt,fSMAllpt);
358     fprintf(outtxt,"\n");
359     fprintf(outtxt," number of Single Undefined :\t"); 
360     fprintf(outtxt,"%i\t%i\t%i\t",fSULowpt,fSUHighpt,fSUAllpt);
361     fprintf(outtxt,"\n");
362     fprintf(outtxt," number of UnlikeSign pair  :\t"); 
363     fprintf(outtxt,"%i\t%i\t%i\t",fUSLowpt,fUSHighpt,fUSAllpt);
364     fprintf(outtxt,"\n");
365     fprintf(outtxt," number of LikeSign pair    :\t");  
366     fprintf(outtxt,"%i\t%i\t%i\t",fLSLowpt,fLSHighpt, fLSAllpt);
367     fprintf(outtxt,"\n");
368     fprintf(outtxt,"===================================================\n");
369     fprintf(outtxt,"\n");
370     fprintf(outtxt,"matching efficiency with the trigger for single tracks = %2d %% \n", effMatch);
371   }//if(pdc06TriggerResponse)
372   
373   else {
374     
375     fprintf(outtxt,"                                                   \n");
376     fprintf(outtxt,"===================================================\n");
377     fprintf(outtxt,"================      ESD SUMMARY    ==============\n");
378     fprintf(outtxt,"                                                   \n");
379     fprintf(outtxt,"    Total number of processed events  %d      \n", nev); 
380     fprintf(outtxt,"\n");
381     fprintf(outtxt,"\n");
382     fprintf(outtxt,"Table 1:                                         \n");
383     fprintf(outtxt," Global Trigger output       Low pt  High pt     \n");
384     fprintf(outtxt," number of Single       :\t");
385     fprintf(outtxt,"%i\t%i\t",fSLowpt,fSHighpt);
386     fprintf(outtxt,"\n");
387     fprintf(outtxt," number of UnlikeSign pair :\t"); 
388     fprintf(outtxt,"%i\t%i\t",fUSLowpt,fUSHighpt);
389     fprintf(outtxt,"\n");
390     fprintf(outtxt," number of LikeSign pair    :\t");  
391     fprintf(outtxt,"%i\t%i\t",fLSLowpt,fLSHighpt);
392     fprintf(outtxt,"\n");
393     fprintf(outtxt,"===================================================\n");
394     fprintf(outtxt,"\n");
395     fprintf(outtxt,"matching efficiency with the trigger for single tracks = %2d %% \n", effMatch);
396   }//else
397   fclose(outtxt);
398   
399   TCanvas * c1 = new TCanvas("c1", "ESD", 400, 10, 600, 700) ;
400   c1->Divide(1,2) ;
401   c1->cd(1) ;
402   fhMUONVertex->Draw() ;
403   c1->cd(2) ;
404   fhMUONMult->Draw() ;  
405   c1->Print("VertexAndMul.eps") ; 
406   TCanvas *c2 = new TCanvas("c2","ESD",400,10,600,700);
407   c2->Divide(1,2);
408   c2->cd(1);
409   hY->Draw();
410   c2->cd(2);
411   hPt->Draw();
412   c2->Print("YandPt.eps") ; 
413 }
414
415 //_____________________________________________________________________________
416 void
417 AliMUONCheck::CheckKine() 
418 {
419   /// Check Stack 
420   
421   AliMUONMCDataInterface diSim(fFileNameSim.Data());
422   if (!diSim.IsValid()) return;
423   
424   Int_t fnevents = diSim.NumberOfEvents();
425   
426   Int_t endOfLoop = fLastEvent+1;
427   
428   if ( fLastEvent == -1 ) endOfLoop = fnevents;
429   
430   Int_t nev=0;
431   Int_t nmu=0;
432   Int_t nonemu=0;
433   Int_t ndimu=0;
434   
435   for ( Int_t ievent = fFirstEvent; ievent < endOfLoop; ++ievent ) 
436   {
437     Int_t nmu2=0;
438     ++nev;  
439     
440     AliStack* stack = diSim.Stack(ievent);
441     Int_t npa = stack->GetNprimary();
442     Int_t npb = stack->GetNtrack(); 
443     printf("Primary particles  %i   \n",npa); 
444     printf("Sec particles  %i   \n",npb); 
445     printf("=================================================================\n") ;
446     printf("Primary particles listing:  \n"); 
447     printf("=================================================================\n") ;
448     for (Int_t i=0; i<npa; ++i) 
449     {
450       TParticle *p  = stack->Particle(i);
451       p->Print("");
452       Int_t pdg=p->GetPdgCode(); 
453       
454       if (abs(pdg) == 13) 
455       {
456         ++nmu2;
457       }
458     }
459     printf("=================================================================\n") ;
460     printf("=================================================================\n") ;
461     
462     printf("Secondaries particles listing:  \n"); 
463     printf("=================================================================\n") ;
464     for (Int_t i=npa; i<npb; ++i) 
465     {
466       stack->Particle(i)->Print("");
467     }
468     
469     printf("=================================================================\n") ; 
470     printf(">>> Event %d, Number of primary particles is %d \n",ievent, npa); 
471     printf(">>> Event %d, Number of secondary articles is %d \n",ievent, npb-npa); 
472     printf("=================================================================\n");
473     if(nmu2>0)
474     {
475       printf(">>> Okay!!! Event %d with at least one muon on primary stack! \n",ievent); 
476       ++nonemu;
477     }
478     
479     if(nmu2==0)
480     {
481       printf(">>> Warning!!! Event %d without muon on primary stack! \n",ievent);     
482       ++nmu;
483     }
484     
485     if(nmu2>1)
486     {
487       printf(">>> Okay!!! Event %d with at least two muons on primary stack! \n",ievent); 
488       ++ndimu; 
489     }
490     printf("=================================================================\n");  
491     printf("                                                                  \n");
492     printf("                                                                  \n") ;
493   }//ievent
494   
495   printf("=================================================================\n") ;
496   printf("               Total number of processed events  %d               \n", nev) ;
497   printf("                                                                 \n") ;
498   
499   if(nmu>0)
500   {
501     printf("--->                       WARNING!!!                       <---\n"); 
502     printf(" %i events without muon on primary stack \n",nmu); 
503   }
504   
505   if(nmu==0)
506   {
507     printf("--->                          OKAY!!!                        <---\n"); 
508     printf("  %i events generated with at least one muon on primary stack \n",nonemu);
509   }
510   
511   if(ndimu>0)
512   {
513     printf("--->                          OKAY!!!                        <---\n"); 
514     printf("  %i events generated with at least two muons on primary stack \n",ndimu); 
515   }
516   
517   printf("                                                                 \n") ;
518   printf("***                       Leaving MuonKine()                 *** \n");
519   printf("**************************************************************** \n");
520   
521   gSystem->cd(fkOutDir);
522   FILE *outtxt=fopen("output.txt","a");
523   freopen("output.txt","a",outtxt);
524   fprintf(outtxt,"                                                   \n");
525   fprintf(outtxt,"=================================================================\n");
526   fprintf(outtxt,"================         MUONkine SUMMARY        ================\n");
527   fprintf(outtxt,"\n");
528   fprintf(outtxt,"=================================================================\n");
529   fprintf(outtxt,"               Total number of processed events  %d              \n", nev) ;
530   fprintf(outtxt,"                                                                 \n");
531   
532   if(nmu>0)
533   {
534     fprintf(outtxt,"                        ---> WARNING!!! <---                     \n"); 
535     fprintf(outtxt,"  %i events without muon on primary stack \n",nmu); 
536   }
537   
538   if(nmu==0)
539   {
540     fprintf(outtxt,"                         ---> OKAY!!! <---                       \n"); 
541     fprintf(outtxt,"  %i events generated with at least one muon on primary stack \n",nonemu); 
542   }
543   
544   if(ndimu>0)
545   {
546     fprintf(outtxt,"                         ---> OKAY!!! <---                       \n"); 
547     fprintf(outtxt,"  %i events generated with at least two muons on primary stack \n",ndimu); 
548   }
549   
550   fprintf(outtxt,"                                                                 \n") ;
551   fprintf(outtxt,"***                       Leaving MuonKine()                 *** \n");
552   fprintf(outtxt,"**************************************************************** \n");
553   fclose(outtxt);
554 }
555
556 //_____________________________________________________________________________
557 void
558 AliMUONCheck::CheckTrackRef() 
559 {
560   /// Check TrackRef files
561   
562   AliMUONMCDataInterface diSim(fFileNameSim.Data());
563   if ( !diSim.IsValid() ) return;
564   
565   Int_t flag11=0,flag12=0,flag13=0,flag14=0;
566   
567   TH1F *tof01= new TH1F("tof01","TOF for first tracking plane",100,0.,100);
568   tof01->SetXTitle("tof (ns)");
569   TH1F *tof14= new TH1F("tof14","TOF for MT22",100,0.,100);
570   tof14->SetXTitle("tof (ns)");
571   
572   TH1F   *hitDensity[4];
573   hitDensity[0] =  new TH1F("TR_dhits01","",30,0,300);
574   hitDensity[0]->SetFillColor(3);
575   hitDensity[0]->SetXTitle("R (cm)");
576   hitDensity[1] =  new TH1F("TR_dhits10","",30,0,300);
577   hitDensity[1]->SetFillColor(3);
578   hitDensity[1]->SetXTitle("R (cm)");
579   hitDensity[2] =  new TH1F("TR_dhits11","",30,0,300);
580   hitDensity[2]->SetFillColor(3);
581   hitDensity[2]->SetXTitle("R (cm)");
582   hitDensity[3] =  new TH1F("TR_dhits14","",30,0,300);
583   hitDensity[3]->SetFillColor(3);
584   hitDensity[3]->SetXTitle("R (cm)");
585   
586   Int_t fnevents = diSim.NumberOfEvents();
587   
588   Int_t endOfLoop = fLastEvent+1;
589   
590   if ( fLastEvent == -1 ) endOfLoop = fnevents;
591   
592   Int_t nev=0;
593   Int_t ntot=fLastEvent+1-fFirstEvent;
594   
595   for ( Int_t ievent = fFirstEvent; ievent < endOfLoop; ++ievent ) 
596   {
597     Int_t  save=-99;
598     ++nev;  
599     
600     Int_t nentries = diSim.NumberOfTrackRefs(ievent);
601     
602     for ( Int_t l=0; l<nentries; ++l )
603     {
604       TClonesArray* trackRefs = diSim.TrackRefs(ievent,l);
605       if (!trackRefs) continue;
606       
607       Int_t nnn = trackRefs->GetEntriesFast();
608       
609       for ( Int_t k=0; k<nnn; ++k ) 
610       {
611         AliTrackReference *tref = static_cast<AliTrackReference*>(trackRefs->UncheckedAt(k));
612         Int_t label = tref->GetTrack();
613         Float_t x     =    tref->X();        // x-pos of hit
614         Float_t y     =    tref->Y();        // y-pos
615         Float_t z     = tref->Z();
616         
617         Float_t r=TMath::Sqrt(x*x+y*y);
618         Float_t time =    tref->GetTime();  
619         
620         Float_t wgt=1/(2*10*TMath::Pi()*r)/(ntot);
621         
622         if (save!=label){
623           save=label;
624           flag11=0;
625           flag12=0;
626           flag13=0;
627           flag14=0;
628         }
629         
630         if (save==label){
631           
632           //Ch 1, z=-526.16
633           if (z<=-521&& z>=-531&&flag11==0){
634             flag11=1;
635             hitDensity[0]->Fill(r,wgt);
636             tof01->Fill(1000000000*time,1);
637           };
638           
639           //Ch 10, z=-1437.6
640           if (z<=-1432&&z>=-1442&&flag12==0){
641             flag12=1;
642             hitDensity[1]->Fill(r,wgt);
643           }
644           
645           //Ch 11, z=-1603.5
646           if (z<=-1598&& z>=-1608&&flag13==0){
647             flag13=1;
648             hitDensity[2]->Fill(r,wgt);
649           };
650           
651           //ch 14 z=-1720.5    
652           if(z<=-1715&&z>=-1725&&flag14==0){
653             flag14=1;
654             hitDensity[3]->Fill(r,wgt);
655             tof14->Fill(1000000000*time,1);
656           }; 
657           
658         }//if save==label
659         
660       }//hits de tTR
661       
662     }//entree de tTR 
663     
664   }//evt loop
665   
666   gSystem->cd(fkOutDir);
667   TCanvas *c6 = new TCanvas("c6","TOF",400,10,600,700);
668   c6->Divide(1,2);
669   c6->cd(1);
670   
671   tof01->Draw();
672   c6->cd(2);
673   tof14->Draw();
674   c6->Print("tof_on_trigger.ps");
675   
676   TCanvas *c5 = new TCanvas("c5","TRef:Hits Density",400,10,600,700);
677   c5->Divide(2,2);
678   c5->cd(1);
679   hitDensity[0]->Draw();
680   c5->cd(2);
681   hitDensity[1]->Draw();
682   c5->cd(3);
683   hitDensity[2]->Draw();
684   c5->cd(4);
685   hitDensity[3]->Draw();
686   c5->Print("TR_Hit_densities.ps");
687   printf("=================================================================\n") ;
688   printf("================  %s Tref SUMMARY    ==============\n", GetName()) ;
689   printf("                                                   \n") ;
690   printf("         Total number of processed events  %d      \n", nev) ;
691   printf("***                Leaving TRef()               *** \n");
692   printf("*************************************************** \n");
693 }
694
695 //_____________________________________________________________________________
696 void 
697 AliMUONCheck::CheckOccupancy(Bool_t perDetEle) const
698 {
699   /// Check occupancy for the first event selected
700   
701   Int_t dEoccupancyBending[14][26];
702   Int_t dEoccupancyNonBending[14][26];
703   Int_t cHoccupancyBending[14];
704   Int_t cHoccupancyNonBending[14];
705   Int_t totaloccupancyBending =0;
706   Int_t totaloccupancyNonBending =0;
707   
708   Int_t dEchannelsBending[14][26];
709   Int_t dEchannelsNonBending[14][26];
710   Int_t cHchannelsBending[14];
711   Int_t cHchannelsNonBending[14];
712   Int_t totalchannelsBending =0;
713   Int_t totalchannelsNonBending =0;
714   
715   Int_t nchambers = AliMUONConstants::NCh();
716   
717   AliMUONDataInterface di(fFileNameSim);
718   
719   AliMUONVDigitStore* digitStore = di.DigitStore(fFirstEvent);
720   
721   // Compute values
722   for (Int_t ichamber=0; ichamber<nchambers; ++ichamber) 
723   {
724     cHchannelsBending[ichamber]=0;
725     cHchannelsNonBending[ichamber]=0;
726     cHoccupancyBending[ichamber]=0;
727     cHoccupancyNonBending[ichamber]=0;
728     
729     for (Int_t idetele=0; idetele<26; idetele++) 
730     {
731       Int_t detele = 100*(ichamber +1)+idetele;
732       
733       if ( AliMpDEManager::IsValidDetElemId(detele) ) 
734       {
735         Int_t cathode(0);
736         
737         const AliMpVSegmentation* segbend = AliMpSegmentation::Instance()
738         ->GetMpSegmentation(detele, AliMp::kCath0);
739         const AliMpVSegmentation* segnonbend = AliMpSegmentation::Instance()
740           ->GetMpSegmentation(detele, AliMp::kCath1);
741         
742         if (AliMpDEManager::GetPlaneType(detele, AliMp::kCath0) != AliMp::kBendingPlane ) 
743         {
744           const AliMpVSegmentation* tmp = segbend;
745           segbend    =  segnonbend;
746           segnonbend =  tmp;
747           cathode = 1;
748         }  
749         
750         Int_t nchannels = segbend->NofPads();
751         Int_t ndigits = digitStore->GetSize(detele,cathode);
752               dEchannelsBending[ichamber][idetele] = nchannels;
753         dEoccupancyBending[ichamber][idetele] = ndigits;
754               cHchannelsBending[ichamber] += nchannels;
755         cHoccupancyBending[ichamber] += ndigits;
756               totalchannelsBending += nchannels;
757         totaloccupancyBending += ndigits;
758         
759         nchannels = segnonbend->NofPads();
760         ndigits = digitStore->GetSize(detele,1-cathode);
761         
762               dEchannelsNonBending[ichamber][idetele] = nchannels;
763         dEoccupancyBending[ichamber][idetele] = ndigits;
764               cHchannelsNonBending[ichamber] += nchannels;
765               cHoccupancyNonBending[ichamber] += ndigits;
766               totalchannelsNonBending += nchannels;
767               totaloccupancyNonBending += ndigits;
768             }
769       if (perDetEle) 
770       {
771         printf(">>> Detection element %4d has %5d channels in bending and %5d channels in nonbending \n",
772                detele, dEchannelsBending[ichamber][idetele], dEchannelsNonBending[ichamber][idetele] ); 
773       }
774     }
775     printf(">>> Chamber %2d has %6d channels in bending and %6d channels in nonbending \n",
776            ichamber+1,  cHchannelsBending[ichamber], cHchannelsNonBending[ichamber]);
777   }
778   printf(">>Spectrometer has  %7d channels in bending and %7d channels in nonbending \n",
779          totalchannelsBending, totalchannelsNonBending);
780   
781   // Output values
782   
783   for ( Int_t ichamber = 0; ichamber < 14; ++ichamber ) 
784   {
785     printf(">>> Chamber %2d  nChannels Bending %5d  nChannels NonBending %5d \n", 
786          ichamber+1, 
787          cHoccupancyBending[ichamber],
788          cHoccupancyNonBending[ichamber]);           
789     printf(">>> Chamber %2d  Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n", 
790          ichamber+1, 
791          100.*((Float_t) cHoccupancyBending[ichamber])/((Float_t) cHchannelsBending[ichamber]),
792          100.*((Float_t) cHoccupancyNonBending[ichamber])/((Float_t) cHchannelsBending[ichamber])            );
793   
794     if ( perDetEle )
795     {
796       for(Int_t idetele=0; idetele<26; idetele++) 
797       {
798         Int_t detele = idetele + 100*(ichamber+1);
799         if ( AliMpDEManager::IsValidDetElemId(detele) ) 
800         {
801           printf(">>> DetEle %4d nChannels Bending %5d  nChannels NonBending %5d \n", 
802                idetele+100*(ichamber+1), 
803                dEoccupancyBending[ichamber][idetele],
804                dEoccupancyNonBending[ichamber][idetele]);  
805           printf(">>> DetEle %4d Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n", 
806                idetele+100*(ichamber+1), 
807                100.*((Float_t) dEoccupancyBending[ichamber][idetele])/((Float_t) dEchannelsBending[ichamber][idetele]),
808                100.*((Float_t) dEoccupancyNonBending[ichamber][idetele])/((Float_t) dEchannelsBending[ichamber][idetele]));  
809         }
810       }
811     }
812   }
813
814   printf(">>> Muon Spectrometer  Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n",  
815        100.*((Float_t) totaloccupancyBending)/((Float_t) totalchannelsBending),
816        100.*((Float_t) totaloccupancyNonBending)/((Float_t) totalchannelsNonBending)            );
817
818 }
819
820 //_____________________________________________________________________________
821 void AliMUONCheck::SetEventsToCheck(Int_t firstEvent, Int_t lastEvent)
822 {
823   /// Set first and last event number to check
824   
825   fFirstEvent = firstEvent;
826   fLastEvent = lastEvent;
827 }