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