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