]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONCheck.C
New macro for checking new misaligner class.
[u/mrichter/AliRoot.git] / MUON / MUONCheck.C
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 // Macro for checking aliroot output and associated files contents
20 // Gines Martinez, Subatech June 2003
21 //
22 #if !defined(__CINT__) || defined(__MAKECINT__)
23 // ROOT includes
24 #include "TBranch.h"
25 #include "TClonesArray.h"
26 #include "TFile.h"
27 #include "TH1.h"
28 #include "TMath.h"
29 #include "TParticle.h"
30 #include "TTree.h"
31
32 // STEER includes
33 #include "AliRun.h"
34 #include "AliRunLoader.h"
35 #include "AliHeader.h"
36 #include "AliLoader.h"
37 #include "AliStack.h"
38
39 // MUON includes
40 #include "AliMUON.h"
41 #include "AliMUONData.h"
42 #include "AliMUONHit.h"
43 #include "AliMUONConstants.h"
44 #include "AliMUONDigit.h"
45 #include "AliMUONRawCluster.h"
46 #include "AliMUONGlobalTrigger.h"
47 #include "AliMUONLocalTrigger.h"
48 #include "AliMUONSegmentationManager.h"
49 #include "AliMUONTrack.h"
50 #include "AliMUONTrackParam.h"
51
52 #include "AliMpVSegmentation.h"
53 #include "AliMpIntPair.h"
54 #endif
55
56
57 void MUONkine(Int_t event2Check=0, char * filename="galice.root")
58 {
59   // Stack of particle for each event
60   AliStack* stack;
61   // Creating Run Loader and openning file containing Hits
62   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
63   if (RunLoader ==0x0) {
64     printf(">>> Error : Error Opening %s file \n",filename);
65     return;
66   }
67
68   RunLoader->LoadKinematics("READ");
69   Int_t ievent, nevents;
70   nevents = RunLoader->GetNumberOfEvents();
71
72   for(ievent=0; ievent<nevents; ievent++) {  // Event loop
73     if (event2Check!=0) ievent=event2Check;
74     Int_t iparticle, nparticles;
75     // Getting event ievent
76     RunLoader->GetEvent(ievent); 
77     stack = RunLoader->Stack();
78     nparticles = (Int_t) stack->GetNtrack();
79     printf(">>> Event %d, Number of particles is %d \n",ievent, nparticles);
80     for(iparticle=0; iparticle<nparticles; iparticle++) {
81       stack->Particle(iparticle)->Print("");  
82     }
83     if (event2Check!=0) ievent=nevents;
84   }
85   RunLoader->UnloadKinematics();
86 }
87
88
89 void MUONhits(Int_t event2Check=0, char * filename="galice.root")
90 {
91   // Creating Run Loader and openning file containing Hits
92   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
93   if (RunLoader ==0x0) {
94     printf(">>> Error : Error Opening %s file \n",filename);
95     return;
96   }
97   // Loading MUON subsystem
98   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
99   MUONLoader->LoadHits("READ");  // Loading Tree of hits for MUON
100   AliMUONData muondata(MUONLoader,"MUON","MUON");  // Creating MUON data container
101   Int_t ievent, nevents;
102   nevents = RunLoader->GetNumberOfEvents();
103
104   for(ievent=0; ievent<nevents; ievent++) {  // Event loop
105     if (event2Check!=0) ievent=event2Check;
106     printf(">>> Event %d \n",ievent);
107     // Getting event ievent
108     RunLoader->GetEvent(ievent); 
109     muondata.SetTreeAddress("H");
110     Int_t itrack, ntracks;
111     ntracks = (Int_t) muondata.GetNtracks();
112     for (itrack=0; itrack<ntracks; itrack++) { // Track loop
113       //Getting List of Hits of Track itrack
114       muondata.GetTrack(itrack);
115       Int_t ihit, nhits;
116       nhits = (Int_t) muondata.Hits()->GetEntriesFast();
117       printf(">>> Track %d, Number of hits %d \n",itrack,nhits);
118       AliMUONHit* mHit;
119       for(ihit=0; ihit<nhits; ihit++) {
120         mHit = static_cast<AliMUONHit*>(muondata.Hits()->At(ihit));
121         Int_t detele   = mHit-> DetElemId(); // Detection element if defined
122         Int_t hittrack = mHit->Track();
123         Float_t x      = mHit->X();
124         Float_t y      = mHit->Y();
125         Float_t z      = mHit->Z();
126         Float_t elos   = mHit->Eloss();
127         Float_t momentum = mHit->Momentum();
128         printf(">>> >>>  Hit%4d DetEle %4d Track%4d (X,Y,Z)=(%7.2f,%7.2f,%8.2f)cm Elost=%7.2gGeV  P=%6.1fGeV/c\n",
129                ihit, detele, hittrack,x,y,z,elos,momentum);
130       }
131       muondata.ResetHits();
132     } // end track loop
133     if (event2Check!=0) ievent=nevents;
134   }  // end event loop
135   MUONLoader->UnloadHits();
136 }
137
138
139 void MUONdigits(Int_t event2Check=0, char * filename="galice.root")
140 {
141   // Creating Run Loader and openning file containing Hits
142   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
143   if (RunLoader ==0x0) {
144     printf(">>> Error : Error Opening %s file \n",filename);
145     return;
146   }
147   // Loading MUON subsystem
148   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
149   MUONLoader->LoadDigits("READ");
150   // Creating MUON data container
151   AliMUONData muondata(MUONLoader,"MUON","MUON");
152
153   Int_t ievent, nevents;
154   nevents = RunLoader->GetNumberOfEvents();
155   AliMUONDigit * mDigit;
156   
157   for(ievent=0; ievent<nevents; ievent++) {
158     if (event2Check!=0) ievent=event2Check;
159     printf(">>> Event %d \n",ievent);
160     RunLoader->GetEvent(ievent);
161   
162     // Addressing
163     Int_t ichamber, nchambers;
164     nchambers = AliMUONConstants::NCh(); ;
165     muondata.SetTreeAddress("D"); 
166     //    char branchname[30];    
167  
168     muondata.GetDigits();
169     // Loop on chambers
170     for( ichamber=0; ichamber<nchambers; ichamber++) {
171       Int_t idigit, ndigits;
172       ndigits = (Int_t) muondata.Digits(ichamber)->GetEntriesFast();
173       for(idigit=0; idigit<ndigits; idigit++) {
174         mDigit = static_cast<AliMUONDigit*>(muondata.Digits(ichamber)->At(idigit));
175         Int_t PadX   = mDigit->PadX();     // Pad X number
176         Int_t PadY   = mDigit->PadY();     // Pad Y number
177         Int_t Signal = mDigit->Signal();   // Physics Signal
178         Int_t Physics= mDigit->Physics();  // Physics contribution to signal
179         //        Int_t Hit    = mDigit->Hit();      // iHit
180         Int_t Cathode= mDigit->Cathode();  // Cathode
181         Int_t Track0 = mDigit->Track(0);
182         Int_t Track1 = mDigit->Track(1); 
183         //Int_t Track2 = mDigit->Track(2);
184         Int_t TCharges0 = mDigit->TrackCharge(0);  //charge per track making this digit (up to 10)
185         Int_t TCharges1 = mDigit->TrackCharge(1);
186         //Int_t TCharges2 = mDigit->TrackCharge(2);
187         Int_t idDE = mDigit->DetElemId();
188         //        printf(">>> Cathode %d\n",Cathode);
189         
190         printf(">>> DetEle %4d Digit%4d Cath %1d (Ix,Iy)=(%3d,%3d) Signal=%4d Physics=%4d Track0=%4d Charge0=%4d Track1=%4d Charge1=%4d \n",
191                idDE, idigit, Cathode, PadX, PadY, Signal, Physics, Track0, 
192                TCharges0, Track1, TCharges1);
193       } // end digit loop
194     } // end chamber loop
195     muondata.ResetDigits();
196     //    } // end cathode loop
197     if (event2Check!=0) ievent=nevents;
198   }  // end event loop
199   MUONLoader->UnloadDigits();
200 }
201
202 void MUONoccupancy(Int_t event2Check=0,  Bool_t perDetEle =kFALSE, char * filename="galice.root") {
203   // Creating Run Loader and openning file containing Hits
204   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
205   if (RunLoader ==0x0) {
206     printf(">>> Error : Error Opening %s file \n",filename);
207     return;
208   }
209   // Loading MUON subsystem
210   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
211   MUONLoader->LoadDigits("READ");
212   // Creating MUON data container
213   AliMUONData muondata(MUONLoader,"MUON","MUON");
214   
215   Int_t ievent, nevents;
216   nevents = RunLoader->GetNumberOfEvents();
217   AliMUONDigit * mDigit =0x0;
218   AliMpVSegmentation * segbend = 0x0;
219   AliMpVSegmentation * segnonbend = 0x0;
220   AliMpIntPair pad(0,0);
221
222   Int_t dEoccupancy_bending[14][26];
223   Int_t dEoccupancy_nonbending[14][26];
224   Int_t cHoccupancy_bending[14];
225   Int_t cHoccupancy_nonbending[14];
226   Int_t totaloccupancy_bending =0;
227   Int_t totaloccupancy_nonbending =0;
228
229   Int_t dEchannels_bending[14][26];
230   Int_t dEchannels_nonbending[14][26];
231   Int_t cHchannels_bending[14];
232   Int_t cHchannels_nonbending[14];
233   Int_t totalchannels_bending =0;
234   Int_t totalchannels_nonbending =0;
235
236   Int_t ichamber, nchambers,idetele, detele, ix, iy;
237   nchambers = AliMUONConstants::NCh(); ;
238
239   for (ichamber=0; ichamber<nchambers; ichamber++) {
240     cHchannels_bending[ichamber]=0;
241     cHchannels_nonbending[ichamber]=0;
242     for (idetele=0; idetele<26; idetele++) {
243       detele= 100*(ichamber +1)+idetele;
244       dEchannels_bending[ichamber][idetele]=0;
245       dEchannels_nonbending[ichamber][idetele]=0;
246       dEoccupancy_bending[ichamber][idetele]=0;
247       dEoccupancy_nonbending[ichamber][idetele]=0;
248       if ( AliMUONSegmentationManager::IsValidDetElemId(detele) ) {
249         segbend    =  AliMUONSegmentationManager::Segmentation(detele, kBendingPlane);
250         segnonbend =  AliMUONSegmentationManager::Segmentation(detele, kNonBendingPlane);
251         for(ix=0; ix<=segbend->MaxPadIndexX(); ix++) {
252           for(iy=0; iy<=segbend->MaxPadIndexY(); iy++) {
253             pad.SetFirst(ix);
254             pad.SetSecond(iy);
255             if( segbend->HasPad(pad) )   {  
256               dEchannels_bending[ichamber][idetele]++;
257               cHchannels_bending[ichamber]++;
258               totalchannels_bending++;
259             }
260           }
261         }
262         for(ix=0; ix<=segnonbend->MaxPadIndexX(); ix++) {
263           for(iy=0; iy<=segnonbend->MaxPadIndexY(); iy++) {
264             pad.SetFirst(ix);
265             pad.SetSecond(iy);
266             if(segnonbend->HasPad(pad))  {
267               dEchannels_nonbending[ichamber][idetele]++;  
268               cHchannels_nonbending[ichamber]++;
269               totalchannels_nonbending++;
270             }
271           }
272         }
273         if (perDetEle) printf(">>> Detection element %4d has %5d channels in bending and %5d channels in nonbending \n",
274              detele, dEchannels_bending[ichamber][idetele], dEchannels_nonbending[ichamber][idetele] ); 
275       }
276     }
277     printf(">>> Chamber %2d has %6d channels in bending and %6d channels in nonbending \n",
278            ichamber+1,  cHchannels_bending[ichamber], cHchannels_nonbending[ichamber]);
279   }
280   printf(">>Spectrometer has  %7d channels in bending and %7d channels in nonbending \n",
281          totalchannels_bending, totalchannels_nonbending);
282
283   ievent=event2Check;
284   printf(">>> Event %d \n",ievent);
285   RunLoader->GetEvent(ievent);
286     
287   // Addressing
288   muondata.SetTreeAddress("D"); 
289   muondata.GetDigits();
290   // Loop on chambers
291   for( ichamber=0; ichamber<nchambers; ichamber++) {
292     cHoccupancy_bending[ichamber]   = 0;
293     cHoccupancy_nonbending[ichamber]= 0;
294     Int_t idigit, ndigits;
295     ndigits = (Int_t) muondata.Digits(ichamber)->GetEntriesFast();
296     for(idigit=0; idigit<ndigits; idigit++) {
297       mDigit = static_cast<AliMUONDigit*>(muondata.Digits(ichamber)->At(idigit));
298       Int_t detele = mDigit->DetElemId();
299       Int_t idetele = detele-(ichamber+1)*100;
300       if ( mDigit->Cathode() == 0 ) {
301
302         cHoccupancy_bending[ichamber]++;
303         dEoccupancy_bending[ichamber][idetele]++;
304         totaloccupancy_bending++;
305       }
306       else {
307         cHoccupancy_nonbending[ichamber]++;
308         dEoccupancy_nonbending[ichamber][idetele]++;
309         totaloccupancy_nonbending++;
310       }
311     } // end digit loop    
312
313     printf(">>> Chamber %2d  nChannels Bending %5d  nChannels NonBending %5d \n", 
314            ichamber+1, 
315            cHoccupancy_bending[ichamber],
316            cHoccupancy_nonbending[ichamber]);           
317     printf(">>> Chamber %2d  Occupancy Bending %5.2f \%  Occupancy NonBending %5.2f \% \n", 
318            ichamber+1, 
319            100.*((Float_t) cHoccupancy_bending[ichamber])/((Float_t) cHchannels_bending[ichamber]),
320            100.*((Float_t) cHoccupancy_nonbending[ichamber])/((Float_t) cHchannels_bending[ichamber])            );
321
322
323     for(Int_t idetele=0; idetele<26; idetele++) {
324       Int_t detele = idetele + 100*(ichamber+1);
325       if ( AliMUONSegmentationManager::IsValidDetElemId(detele) ) {
326         if (perDetEle) {
327           printf(">>> DetEle %4d nChannels Bending %5d  nChannels NonBending %5d \n", 
328                  idetele+100*(ichamber+1), 
329                  dEoccupancy_bending[ichamber][idetele],
330                  dEoccupancy_nonbending[ichamber][idetele]);  
331           printf(">>> DetEle %4d Occupancy Bending %5.2f \%  Occupancy NonBending %5.2f \% \n", 
332                  idetele+100*(ichamber+1), 
333                  100.*((Float_t) dEoccupancy_bending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]),
334                  100.*((Float_t) dEoccupancy_nonbending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]));  
335         }
336       }
337     }
338   } // end chamber loop
339   printf(">>> Muon Spectrometer  Occupancy Bending %5.2f\%  Occupancy NonBending %5.2f\% \n",  
340            100.*((Float_t) totaloccupancy_bending)/((Float_t) totalchannels_bending),
341          100.*((Float_t) totaloccupancy_nonbending)/((Float_t) totalchannels_nonbending)            );
342   muondata.ResetDigits();
343   //    } // end cathode loop
344   MUONLoader->UnloadDigits();
345 }
346
347 void MUONrecpoints(Int_t event2Check=0, char * filename="galice.root") {
348
349   // Creating Run Loader and openning file containing Hits
350   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
351   if (RunLoader ==0x0) {
352     printf(">>> Error : Error Opening %s file \n",filename);
353     return;
354   }
355   // Getting MUONloader
356   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
357   MUONLoader->LoadRecPoints("READ");
358   // Creating MUON data container
359   AliMUONData muondata(MUONLoader,"MUON","MUON");
360
361   Int_t ievent, nevents;
362   nevents = RunLoader->GetNumberOfEvents();
363   AliMUONRawCluster * mRecPoint = 0;
364   
365   for(ievent=0; ievent<nevents; ievent++) {
366     if (event2Check!=0) ievent=event2Check;
367     printf(">>> Event %d \n",ievent);
368     RunLoader->GetEvent(ievent);
369     // Addressing
370     Int_t ichamber, nchambers;
371     nchambers = AliMUONConstants::NTrackingCh();
372     muondata.SetTreeAddress("RC"); 
373     char branchname[30];    
374     muondata.GetRawClusters();
375     // Loop on chambers
376     for( ichamber=0; ichamber<nchambers; ichamber++) {
377       sprintf(branchname,"MUONRawClusters%d",ichamber+1);
378       //printf(">>>  branchname %s\n",branchname);
379       Int_t irecpoint, nrecpoints;
380       nrecpoints = (Int_t) muondata.RawClusters(ichamber)->GetEntriesFast();
381       // printf(">>> Chamber %2d, Number of recpoints = %6d \n",ichamber+1, nrecpoints);
382       for(irecpoint=0; irecpoint<nrecpoints; irecpoint++) {
383         mRecPoint = static_cast<AliMUONRawCluster*>(muondata.RawClusters(ichamber)->At(irecpoint));
384         Int_t Track0 = mRecPoint->GetTrack(0);
385         Int_t Track1 = mRecPoint->GetTrack(1); 
386         Int_t Track2 = mRecPoint->GetTrack(2);
387         Int_t Q0 = mRecPoint->GetCharge(0);
388         Int_t Q1 = mRecPoint->GetCharge(1);
389         Float_t x0 = mRecPoint->GetX(0);
390         Float_t x1 = mRecPoint->GetX(1);
391         Float_t y0 = mRecPoint->GetY(0);
392         Float_t y1 = mRecPoint->GetY(1);
393         Float_t z0 = mRecPoint->GetZ(0);
394         Float_t z1 = mRecPoint->GetZ(1);
395         Float_t chi2_0 =  mRecPoint->GetChi2(0);
396         //Float_t chi2_1 =  mRecPoint->GetChi2(1);
397         Int_t de = mRecPoint->GetDetElemId();
398         printf(">>> >>> RecPoint %4d  DetEle %4d (X,Y,Z)=(%7.2f,%7.2f,%8.2f)cm  Q0=%4d  Q1=%4d Hit=%4d Track1=%4d Track2=%4d Chi2=%6.3f \n",
399                irecpoint,de,x0,y0,z0,Q0,Q1,Track0, Track1, Track2, chi2_0);
400         if( (x0!=x1) || (y0!=y1) || (z0!=z1) )
401           printf(">>> >>> Warning (X0,Y0,Z0)=(%7.2f, %7.2f, %8.2f)cm != (X1,Y1,Z1)=(%7.2f,%7.2f,%8.2f)cm \n",x0,y0,z0,x1,y1,z1); 
402       } // end recpoint loop
403     } // end chamber loop
404     muondata.ResetRawClusters();
405     if (event2Check!=0) ievent=nevents;
406   }  // end event loop
407   MUONLoader->UnloadRecPoints();
408 }
409
410 void MUONrectrigger (Int_t event2Check=0, char * filename="galice.root"){
411 // reads and dumps trigger objects from MUON.RecPoints.root
412   TClonesArray * globalTrigger;
413   TClonesArray * localTrigger;
414   
415   // Creating Run Loader and openning file containing Hits
416   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
417   if (RunLoader ==0x0) {
418     printf(">>> Error : Error Opening %s file \n",filename);
419     return;
420   }
421   
422   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
423   MUONLoader->LoadDigits("READ");
424   // Creating MUON data container
425   AliMUONData muondata(MUONLoader,"MUON","MUON");
426   
427   
428   Int_t ievent, nevents;
429   nevents = RunLoader->GetNumberOfEvents();
430   
431   AliMUONGlobalTrigger *gloTrg;
432   AliMUONLocalTrigger *locTrg;
433   
434   for (ievent=0; ievent<nevents; ievent++) {
435     if (event2Check!=0) ievent=event2Check;
436     RunLoader->GetEvent(ievent);
437     
438     muondata.SetTreeAddress("GLT"); 
439     muondata.GetTriggerD();
440     
441     globalTrigger = muondata.GlobalTrigger();
442     localTrigger = muondata.LocalTrigger();
443     
444     Int_t nglobals = (Int_t) globalTrigger->GetEntriesFast(); // should be 1
445     Int_t nlocals  = (Int_t) localTrigger->GetEntriesFast(); // up to 234
446     printf("###################################################\n");
447     cout << " event " << ievent 
448          << " nglobals nlocals: " << nglobals << " " << nlocals << "\n"; 
449     
450     for (Int_t iglobal=0; iglobal<nglobals; iglobal++) { // Global Trigger
451       gloTrg = static_cast<AliMUONGlobalTrigger*>(globalTrigger->At(iglobal));
452       
453       printf("===================================================\n");
454       printf(" Global Trigger output       Low pt  High pt   All\n");
455       printf(" number of Single Plus      :\t");
456       printf("%i\t%i\t%i\t",gloTrg->SinglePlusLpt(),
457              gloTrg->SinglePlusHpt(),gloTrg->SinglePlusApt());
458       printf("\n");
459       printf(" number of Single Minus     :\t");
460       printf("%i\t%i\t%i\t",gloTrg->SingleMinusLpt(),
461              gloTrg->SingleMinusHpt(),gloTrg->SingleMinusApt());
462       printf("\n");
463       printf(" number of Single Undefined :\t"); 
464       printf("%i\t%i\t%i\t",gloTrg->SingleUndefLpt(),
465              gloTrg->SingleUndefHpt(),gloTrg->SingleUndefApt());
466       printf("\n");
467       printf(" number of UnlikeSign pair  :\t"); 
468       printf("%i\t%i\t%i\t",gloTrg->PairUnlikeLpt(),
469              gloTrg->PairUnlikeHpt(),gloTrg->PairUnlikeApt());
470       printf("\n");
471       printf(" number of LikeSign pair    :\t");  
472       printf("%i\t%i\t%i\t",gloTrg->PairLikeLpt(),
473              gloTrg->PairLikeHpt(),gloTrg->PairLikeApt());
474       printf("\n");
475       printf("===================================================\n");
476       
477     } // end of loop on Global Trigger
478     
479     for (Int_t ilocal=0; ilocal<nlocals; ilocal++) { // Local Trigger
480       cout << " >>> Output for Local Trigger " << ilocal << "\n";
481       
482       locTrg = static_cast<AliMUONLocalTrigger*>(localTrigger->At(ilocal));
483       
484       cout << "Circuit StripX Dev StripY: " 
485            << locTrg->LoCircuit() << " " 
486            << locTrg->LoStripX() << " " 
487            << locTrg->LoDev() << " " 
488            << locTrg->LoStripY() 
489            << "\n";
490       cout << "Lpt Hpt Apt: "     
491            << locTrg->LoLpt() << " "   
492            << locTrg->LoHpt() << " "  
493            << locTrg->LoApt() << "\n";
494       
495     } // end of loop on Local Trigger
496     muondata.ResetTrigger();
497     if (event2Check!=0) ievent=nevents;
498   } // end loop on event  
499   MUONLoader->UnloadRecPoints();
500 }
501
502
503
504 void MUONrectracks (Int_t event2Check=0, char * filename="galice.root"){
505 // reads and dumps trigger objects from MUON.RecPoints.root
506   TClonesArray * RecTracks;
507   
508   // Creating Run Loader and openning file containing Hits
509   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
510   if (RunLoader ==0x0) {
511     printf(">>> Error : Error Opening %s file \n",filename);
512     return;
513   }
514   
515   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
516   MUONLoader->LoadTracks("READ");
517   // Creating MUON data container
518   AliMUONData muondata(MUONLoader,"MUON","MUON");
519   
520     Int_t ievent, nevents;
521   nevents = RunLoader->GetNumberOfEvents();
522   
523   //  AliMUONTrack * rectrack;
524   
525   for (ievent=0; ievent<nevents; ievent++) {
526     if (event2Check!=0) ievent=event2Check;
527     RunLoader->GetEvent(ievent);
528     
529     muondata.SetTreeAddress("RT");
530     muondata.GetRecTracks();
531     RecTracks = muondata.RecTracks();
532     
533     
534     Int_t nrectracks = (Int_t) RecTracks->GetEntriesFast(); //
535
536     printf(">>> Event %d, Number of Recconstructed tracks %d \n",ievent, nrectracks);
537     // loop over tracks
538  
539  
540     Int_t nTrackHits;// nPrimary;
541     Double_t fitFmin;
542     Double_t bendingSlope, nonBendingSlope, inverseBendingMomentum;
543     Double_t xRec, yRec, zRec, chi2MatchTrigger;
544     Bool_t matchTrigger;
545
546   // setting pointer for tracks, triggertracks & trackparam at vertex
547     AliMUONTrack* recTrack = 0;
548     AliMUONTrackParam* trackParam = 0;
549
550     for (Int_t iRecTracks = 0; iRecTracks <  nrectracks;  iRecTracks++) {
551     // reading info from tracks
552       recTrack = (AliMUONTrack*) RecTracks->At(iRecTracks);
553       trackParam = (AliMUONTrackParam*) (recTrack->GetTrackParamAtHit())->First();
554       trackParam->ExtrapToZ(0.0);
555       bendingSlope            = trackParam->GetBendingSlope();
556       nonBendingSlope         = trackParam->GetNonBendingSlope();
557       inverseBendingMomentum = trackParam->GetInverseBendingMomentum();
558       xRec  = trackParam->GetNonBendingCoor();
559       yRec  = trackParam->GetBendingCoor();
560       zRec  = trackParam->GetZ();
561
562       nTrackHits       = recTrack->GetNTrackHits();
563       fitFmin          = recTrack->GetFitFMin();
564       matchTrigger     = recTrack->GetMatchTrigger();
565       chi2MatchTrigger = recTrack->GetChi2MatchTrigger();
566       
567       printf(">>> RecTrack %4d  NofClusters=%2d BendMomentum=%7.2f NonBendSlope=%5.2f  BendSlope=%5.2f Match2Trig=%1d (vertex@z=0)=(%5.2f,%5.2f,%5.1f)cm \n", iRecTracks, nTrackHits, 1/inverseBendingMomentum , nonBendingSlope*180./TMath::Pi(), bendingSlope*180./TMath::Pi(),  matchTrigger, xRec,yRec,zRec);
568     } // end loop tracks
569
570     muondata.ResetRecTracks();
571     if (event2Check!=0) ievent=nevents;
572   } // end loop on event  
573   MUONLoader->UnloadTracks();
574 }
575
576
577
578
579
580
581
582
583
584