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