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