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