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