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