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