]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONCheck.C
Additional protection in case of negative indexes. More investigation is needed
[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 #include "TNtuple.h"
32
33 // STEER includes
34 #include "AliRun.h"
35 #include "AliRunLoader.h"
36 #include "AliHeader.h"
37 #include "AliLoader.h"
38 #include "AliStack.h"
39
40 // MUON includes
41 #include "AliMUON.h"
42 #include "AliMUONData.h"
43 #include "AliMUONHit.h"
44 #include "AliMUONConstants.h"
45 #include "AliMUONDigit.h"
46 #include "AliMUONRawCluster.h"
47 #include "AliMUONGlobalTrigger.h"
48 #include "AliMUONLocalTrigger.h"
49 #include "AliMUONTrack.h"
50 #include "AliMUONTrackParam.h"
51
52 #include "AliMpVSegmentation.h"
53 #include "AliMpIntPair.h"
54 #include "AliMpDEManager.h"
55 #include "AliMpSegFactory.h"
56 #endif
57
58
59 void MUONkine(Int_t event2Check=0, char * filename="galice.root")
60 {
61   // Stack of particle for each event
62   AliStack* stack;
63   // Creating Run Loader and openning file containing Hits
64   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
65   if (RunLoader ==0x0) {
66     printf(">>> Error : Error Opening %s file \n",filename);
67     return;
68   }
69
70   RunLoader->LoadKinematics("READ");
71   Int_t ievent, nevents;
72   nevents = RunLoader->GetNumberOfEvents();
73
74   for(ievent=0; ievent<nevents; ievent++) {  // Event loop
75     if (event2Check!=0) ievent=event2Check;
76     Int_t iparticle, nparticles;
77     // Getting event ievent
78     RunLoader->GetEvent(ievent); 
79     stack = RunLoader->Stack();
80     nparticles = (Int_t) stack->GetNtrack();
81     printf(">>> Event %d, Number of particles is %d \n",ievent, nparticles);
82     for(iparticle=0; iparticle<nparticles; iparticle++) {
83       stack->Particle(iparticle)->Print("");  
84     }
85     if (event2Check!=0) ievent=nevents;
86   }
87   RunLoader->UnloadKinematics();
88 }
89
90
91 void MUONhits(Int_t event2Check=0, char * filename="galice.root")
92 {
93   // Creating Run Loader and openning file containing Hits
94   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
95   if (RunLoader ==0x0) {
96     printf(">>> Error : Error Opening %s file \n",filename);
97     return;
98   }
99   // Loading MUON subsystem
100   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
101   MUONLoader->LoadHits("READ");  // Loading Tree of hits for MUON
102   AliMUONData muondata(MUONLoader,"MUON","MUON");  // Creating MUON data container
103   Int_t ievent, nevents;
104   nevents = RunLoader->GetNumberOfEvents();
105
106   for(ievent=0; ievent<nevents; ievent++) {  // Event loop
107     if (event2Check!=0) ievent=event2Check;
108     printf(">>> Event %d \n",ievent);
109     // Getting event ievent
110     RunLoader->GetEvent(ievent); 
111     muondata.SetTreeAddress("H");
112     Int_t itrack, ntracks;
113     ntracks = (Int_t) muondata.GetNtracks();
114     for (itrack=0; itrack<ntracks; itrack++) { // Track loop
115       //Getting List of Hits of Track itrack
116       muondata.GetTrack(itrack);
117       Int_t ihit, nhits;
118       nhits = (Int_t) muondata.Hits()->GetEntriesFast();
119       printf(">>> Track %d, Number of hits %d \n",itrack,nhits);
120       AliMUONHit* mHit;
121       for(ihit=0; ihit<nhits; ihit++) {
122         mHit = static_cast<AliMUONHit*>(muondata.Hits()->At(ihit));
123         Int_t detele   = mHit-> DetElemId(); // Detection element if defined
124         Int_t hittrack = mHit->Track();
125         Float_t x      = mHit->X();
126         Float_t y      = mHit->Y();
127         Float_t z      = mHit->Z();
128         Float_t elos   = mHit->Eloss();
129         Float_t momentum = mHit->Momentum();
130         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",
131                ihit, detele, hittrack,x,y,z,elos,momentum);
132       }
133       muondata.ResetHits();
134     } // end track loop
135     if (event2Check!=0) ievent=nevents;
136   }  // end event loop
137   MUONLoader->UnloadHits();
138 }
139
140
141
142 void MUONdigits(Int_t event2Check=0, char * filename="galice.root")
143 {
144   // Creating Run Loader and openning file containing Hits
145   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
146   if (RunLoader ==0x0) {
147     printf(">>> Error : Error Opening %s file \n",filename);
148     return;
149   }
150   // Loading MUON subsystem
151   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
152   MUONLoader->LoadDigits("READ");
153   // Creating MUON data container
154   AliMUONData muondata(MUONLoader,"MUON","MUON");
155   
156   Int_t ievent, nevents;
157   nevents = RunLoader->GetNumberOfEvents();
158   AliMUONDigit * mDigit;
159   
160   for(ievent=0; ievent<nevents; ievent++) {
161     if (event2Check!=0) ievent=event2Check;
162     printf(">>> Event %d \n",ievent);
163     RunLoader->GetEvent(ievent);
164     
165     // Addressing
166     Int_t ichamber, nchambers;
167     nchambers = AliMUONConstants::NCh(); ;
168     muondata.SetTreeAddress("D,GLT");
169     
170     muondata.GetDigits();
171     // Loop on chambers
172     for( ichamber=0; ichamber<nchambers; ichamber++) {
173       Int_t idigit, ndigits;
174       TClonesArray* digits = muondata.Digits(ichamber);
175       digits->Sort();
176       ndigits = (Int_t)digits->GetEntriesFast();
177       for(idigit=0; idigit<ndigits; idigit++) {
178         mDigit = static_cast<AliMUONDigit*>(digits->At(idigit));
179         mDigit->Print("tracks");
180       } // end digit loop
181     } // end chamber loop
182     muondata.ResetDigits();
183     if (event2Check!=0) ievent=nevents;
184   }  // end event loop
185   MUONLoader->UnloadDigits();
186 }
187
188 void MUONsdigits(Int_t event2Check=0, char * filename="galice.root")
189 {
190   // Creating Run Loader and openning file containing Hits
191   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
192   if (RunLoader ==0x0) {
193     printf(">>> Error : Error Opening %s file \n",filename);
194     return;
195   }
196   // Loading MUON subsystem
197   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
198   MUONLoader->LoadSDigits("READ");
199   // Creating MUON data container
200   AliMUONData muondata(MUONLoader,"MUON","MUON");
201   
202   Int_t ievent, nevents;
203   nevents = RunLoader->GetNumberOfEvents();
204   AliMUONDigit * mDigit;
205   
206   for(ievent=0; ievent<nevents; ievent++) {
207     if (event2Check!=0) ievent=event2Check;
208     printf(">>> Event %d \n",ievent);
209     RunLoader->GetEvent(ievent);
210     
211     // Addressing
212     Int_t ichamber, nchambers;
213     nchambers = AliMUONConstants::NCh(); ;
214     muondata.SetTreeAddress("S");
215     
216     muondata.GetSDigits();
217     // Loop on chambers
218     for( ichamber=0; ichamber<nchambers; ichamber++) {
219       Int_t idigit, ndigits;
220       TClonesArray* digits = muondata.SDigits(ichamber);
221       ndigits = (Int_t)digits->GetEntriesFast();
222       for(idigit=0; idigit<ndigits; idigit++) {
223         mDigit = static_cast<AliMUONDigit*>(digits->At(idigit));
224         mDigit->Print("tracks");
225       } // end digit loop
226     } // end chamber loop
227     muondata.ResetSDigits();
228     if (event2Check!=0) ievent=nevents;
229   }  // end event loop
230   MUONLoader->UnloadSDigits();
231 }
232
233 void MUONoccupancy(Int_t event2Check=0,  Bool_t perDetEle =kFALSE, char * filename="galice.root") {
234   // Creating Run Loader and openning file containing Hits
235   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
236   if (RunLoader ==0x0) {
237     printf(">>> Error : Error Opening %s file \n",filename);
238     return;
239   }
240   // Loading MUON subsystem
241   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
242   MUONLoader->LoadDigits("READ");
243   // Creating MUON data container
244   AliMUONData muondata(MUONLoader,"MUON","MUON");
245   
246   Int_t ievent, nevents;
247   nevents = RunLoader->GetNumberOfEvents();
248   AliMUONDigit * mDigit =0x0;
249   AliMpVSegmentation * segbend = 0x0;
250   AliMpVSegmentation * segnonbend = 0x0;
251   AliMpIntPair pad(0,0);
252
253   Int_t dEoccupancy_bending[14][26];
254   Int_t dEoccupancy_nonbending[14][26];
255   Int_t cHoccupancy_bending[14];
256   Int_t cHoccupancy_nonbending[14];
257   Int_t totaloccupancy_bending =0;
258   Int_t totaloccupancy_nonbending =0;
259
260   Int_t dEchannels_bending[14][26];
261   Int_t dEchannels_nonbending[14][26];
262   Int_t cHchannels_bending[14];
263   Int_t cHchannels_nonbending[14];
264   Int_t totalchannels_bending =0;
265   Int_t totalchannels_nonbending =0;
266
267   Int_t ichamber, nchambers,idetele, detele, ix, iy;
268   nchambers = AliMUONConstants::NCh(); ;
269
270   AliMpSegFactory factory;
271
272   for (ichamber=0; ichamber<nchambers; ichamber++) {
273     cHchannels_bending[ichamber]=0;
274     cHchannels_nonbending[ichamber]=0;
275     for (idetele=0; idetele<26; idetele++) {
276       detele= 100*(ichamber +1)+idetele;
277       dEchannels_bending[ichamber][idetele]=0;
278       dEchannels_nonbending[ichamber][idetele]=0;
279       dEoccupancy_bending[ichamber][idetele]=0;
280       dEoccupancy_nonbending[ichamber][idetele]=0;
281       if ( AliMpDEManager::IsValidDetElemId(detele) ) {
282         
283         segbend    =  factory.CreateMpSegmentation(detele, 0);
284         segnonbend =  factory.CreateMpSegmentation(detele, 1);
285         if (AliMpDEManager::GetPlaneType(detele, 0) != kBendingPlane ) {
286           AliMpVSegmentation* tmp = segbend;
287           segbend    =  segnonbend;
288           segnonbend =  tmp;
289         }  
290           
291         for(ix=0; ix<=segbend->MaxPadIndexX(); ix++) {
292           for(iy=0; iy<=segbend->MaxPadIndexY(); iy++) {
293             pad.SetFirst(ix);
294             pad.SetSecond(iy);
295             if( segbend->HasPad(pad) )   {  
296               dEchannels_bending[ichamber][idetele]++;
297               cHchannels_bending[ichamber]++;
298               totalchannels_bending++;
299             }
300           }
301         }
302         for(ix=0; ix<=segnonbend->MaxPadIndexX(); ix++) {
303           for(iy=0; iy<=segnonbend->MaxPadIndexY(); iy++) {
304             pad.SetFirst(ix);
305             pad.SetSecond(iy);
306             if(segnonbend->HasPad(pad))  {
307               dEchannels_nonbending[ichamber][idetele]++;  
308               cHchannels_nonbending[ichamber]++;
309               totalchannels_nonbending++;
310             }
311           }
312         }
313         if (perDetEle) printf(">>> Detection element %4d has %5d channels in bending and %5d channels in nonbending \n",
314              detele, dEchannels_bending[ichamber][idetele], dEchannels_nonbending[ichamber][idetele] ); 
315       }
316     }
317     printf(">>> Chamber %2d has %6d channels in bending and %6d channels in nonbending \n",
318            ichamber+1,  cHchannels_bending[ichamber], cHchannels_nonbending[ichamber]);
319   }
320   printf(">>Spectrometer has  %7d channels in bending and %7d channels in nonbending \n",
321          totalchannels_bending, totalchannels_nonbending);
322
323   factory.DeleteSegmentations();
324
325   ievent=event2Check;
326   printf(">>> Event %d \n",ievent);
327   RunLoader->GetEvent(ievent);
328     
329   // Addressing
330   muondata.SetTreeAddress("D"); 
331   muondata.GetDigits();
332   // Loop on chambers
333   for( ichamber=0; ichamber<nchambers; ichamber++) {
334     cHoccupancy_bending[ichamber]   = 0;
335     cHoccupancy_nonbending[ichamber]= 0;
336     Int_t idigit, ndigits;
337     ndigits = (Int_t) muondata.Digits(ichamber)->GetEntriesFast();
338     for(idigit=0; idigit<ndigits; idigit++) {
339       mDigit = static_cast<AliMUONDigit*>(muondata.Digits(ichamber)->At(idigit));
340       Int_t detele = mDigit->DetElemId();
341       Int_t idetele = detele-(ichamber+1)*100;
342       if ( mDigit->Cathode() == 0 ) {
343
344         cHoccupancy_bending[ichamber]++;
345         dEoccupancy_bending[ichamber][idetele]++;
346         totaloccupancy_bending++;
347       }
348       else {
349         cHoccupancy_nonbending[ichamber]++;
350         dEoccupancy_nonbending[ichamber][idetele]++;
351         totaloccupancy_nonbending++;
352       }
353     } // end digit loop    
354
355     printf(">>> Chamber %2d  nChannels Bending %5d  nChannels NonBending %5d \n", 
356            ichamber+1, 
357            cHoccupancy_bending[ichamber],
358            cHoccupancy_nonbending[ichamber]);           
359     printf(">>> Chamber %2d  Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n", 
360            ichamber+1, 
361            100.*((Float_t) cHoccupancy_bending[ichamber])/((Float_t) cHchannels_bending[ichamber]),
362            100.*((Float_t) cHoccupancy_nonbending[ichamber])/((Float_t) cHchannels_bending[ichamber])            );
363
364
365     for(Int_t idetele=0; idetele<26; idetele++) {
366       Int_t detele = idetele + 100*(ichamber+1);
367       if ( AliMpDEManager::IsValidDetElemId(detele) ) {
368         if (perDetEle) {
369           printf(">>> DetEle %4d nChannels Bending %5d  nChannels NonBending %5d \n", 
370                  idetele+100*(ichamber+1), 
371                  dEoccupancy_bending[ichamber][idetele],
372                  dEoccupancy_nonbending[ichamber][idetele]);  
373           printf(">>> DetEle %4d Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n", 
374                  idetele+100*(ichamber+1), 
375                  100.*((Float_t) dEoccupancy_bending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]),
376                  100.*((Float_t) dEoccupancy_nonbending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]));  
377         }
378       }
379     }
380   } // end chamber loop
381   printf(">>> Muon Spectrometer  Occupancy Bending %5.2f %%  Occupancy NonBending %5.2f %% \n",  
382            100.*((Float_t) totaloccupancy_bending)/((Float_t) totalchannels_bending),
383          100.*((Float_t) totaloccupancy_nonbending)/((Float_t) totalchannels_nonbending)            );
384   muondata.ResetDigits();
385   //    } // end cathode loop
386   MUONLoader->UnloadDigits();
387 }
388
389 void MUONrecpoints(Int_t event2Check=0, char * filename="galice.root") {
390
391   // Creating Run Loader and openning file containing Hits
392   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
393   if (RunLoader ==0x0) {
394     printf(">>> Error : Error Opening %s file \n",filename);
395     return;
396   }
397   // Getting MUONloader
398   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
399   MUONLoader->LoadRecPoints("READ");
400   // Creating MUON data container
401   AliMUONData muondata(MUONLoader,"MUON","MUON");
402
403   Int_t ievent, nevents;
404   nevents = RunLoader->GetNumberOfEvents();
405   AliMUONRawCluster * mRecPoint = 0;
406   
407   for(ievent=0; ievent<nevents; ievent++) {
408     if (event2Check!=0) ievent=event2Check;
409     printf(">>> Event %d \n",ievent);
410     RunLoader->GetEvent(ievent);
411     // Addressing
412     Int_t ichamber, nchambers;
413     nchambers = AliMUONConstants::NTrackingCh();
414     muondata.SetTreeAddress("RC,TC"); 
415     char branchname[30];    
416     muondata.GetRawClusters();
417     // Loop on chambers
418     for( ichamber=0; ichamber<nchambers; ichamber++) {
419       sprintf(branchname,"MUONRawClusters%d",ichamber+1);
420       //printf(">>>  branchname %s\n",branchname);
421       Int_t irecpoint, nrecpoints;
422       nrecpoints = (Int_t) muondata.RawClusters(ichamber)->GetEntriesFast();
423       // printf(">>> Chamber %2d, Number of recpoints = %6d \n",ichamber+1, nrecpoints);
424       for(irecpoint=0; irecpoint<nrecpoints; irecpoint++) {
425         mRecPoint = static_cast<AliMUONRawCluster*>(muondata.RawClusters(ichamber)->At(irecpoint));
426         Int_t Track0 = mRecPoint->GetTrack(0);
427         Int_t Track1 = mRecPoint->GetTrack(1); 
428         Int_t Track2 = mRecPoint->GetTrack(2);
429         Int_t Q0 = mRecPoint->GetCharge(0);
430         Int_t Q1 = mRecPoint->GetCharge(1);
431         Float_t x0 = mRecPoint->GetX(0);
432         Float_t x1 = mRecPoint->GetX(1);
433         Float_t y0 = mRecPoint->GetY(0);
434         Float_t y1 = mRecPoint->GetY(1);
435         Float_t z0 = mRecPoint->GetZ(0);
436         Float_t z1 = mRecPoint->GetZ(1);
437         Float_t chi2_0 =  mRecPoint->GetChi2(0);
438         //Float_t chi2_1 =  mRecPoint->GetChi2(1);
439         Int_t de = mRecPoint->GetDetElemId();
440         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",
441                irecpoint,de,x0,y0,z0,Q0,Q1,Track0, Track1, Track2, chi2_0);
442         if( (x0!=x1) || (y0!=y1) || (z0!=z1) )
443           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); 
444       } // end recpoint loop
445     } // end chamber loop
446     muondata.ResetRawClusters();
447     if (event2Check!=0) ievent=nevents;
448   }  // end event loop
449   MUONLoader->UnloadRecPoints();
450 }
451
452 void MUONrectrigger (Int_t event2Check=0, char * filename="galice.root", Int_t WRITE = 0)
453 {
454   // reads and dumps trigger objects from MUON.RecPoints.root
455   TClonesArray * globalTrigger;
456   TClonesArray * localTrigger;
457   
458   // Do NOT print out all the info if the loop runs over all events 
459   Int_t PRINTOUT = (event2Check == 0 ) ? 0 : 1 ;  
460
461   // Book a ntuple for more detailled studies
462   TNtuple *TgtupleGlo = new TNtuple("TgtupleGlo","Global Trigger Ntuple","ev:global:spapt:smapt:undefapt:uplpt:uphpt:upapt:lpapt");
463   TNtuple *TgtupleLoc = new TNtuple("TgtupleLoc","Local Trigger Ntuple","LoCircuit:LoStripX:LoDev:StripY:LoLpt:LoHpt:LoApt");
464
465   // counters
466   Int_t SPLowpt=0,SPHighpt=0,SPAllpt=0;
467   Int_t SMLowpt=0,SMHighpt=0,SMAllpt=0;
468   Int_t SULowpt=0,SUHighpt=0,SUAllpt=0;
469   Int_t USLowpt=0,USHighpt=0,USAllpt=0;
470   Int_t LSLowpt=0,LSHighpt=0,LSAllpt=0;
471
472   // Creating Run Loader and openning file containing Hits
473   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
474   if (RunLoader ==0x0) {
475     printf(">>> Error : Error Opening %s file \n",filename);
476     return;
477   }
478   
479   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
480   MUONLoader->LoadDigits("READ");
481   // Creating MUON data container
482   AliMUONData muondata(MUONLoader,"MUON","MUON");
483   
484   
485   Int_t ievent, nevents;
486   nevents = RunLoader->GetNumberOfEvents();
487   
488   AliMUONGlobalTrigger *gloTrg(0x0);
489   AliMUONLocalTrigger *locTrg(0x0);
490   
491   for (ievent=0; ievent<nevents; ievent++) {
492     if (event2Check!=0) ievent=event2Check;
493     if (ievent%100==0 || event2Check) cout << "Processing event " << ievent << endl;
494     RunLoader->GetEvent(ievent);
495     
496     muondata.SetTreeAddress("D,GLT"); 
497     muondata.GetTriggerD();
498     
499     globalTrigger = muondata.GlobalTrigger();
500     localTrigger = muondata.LocalTrigger();
501     
502     Int_t nglobals = (Int_t) globalTrigger->GetEntriesFast(); // should be 1
503     Int_t nlocals  = (Int_t) localTrigger->GetEntriesFast(); // up to 234
504     if (PRINTOUT) printf("###################################################\n");
505     if (PRINTOUT) {cout << " event " << ievent 
506                         << " nglobals nlocals: " << nglobals << " " << nlocals << "\n"; }
507     
508     for (Int_t iglobal=0; iglobal<nglobals; iglobal++) { // Global Trigger
509       gloTrg = static_cast<AliMUONGlobalTrigger*>(globalTrigger->At(iglobal));
510       
511       SPLowpt+=gloTrg->SinglePlusLpt() ;
512       SPHighpt+=gloTrg->SinglePlusHpt() ;
513       SPAllpt+=gloTrg->SinglePlusApt() ;
514       SMLowpt+=gloTrg->SingleMinusLpt();
515       SMHighpt+=gloTrg->SingleMinusHpt();
516       SMAllpt+=gloTrg->SingleMinusApt();
517       SULowpt+=gloTrg->SingleUndefLpt();
518       SUHighpt+=gloTrg->SingleUndefHpt();
519       SUAllpt+=gloTrg->SingleUndefApt();
520       USLowpt+=gloTrg->PairUnlikeLpt(); 
521       USHighpt+=gloTrg->PairUnlikeHpt();
522       USAllpt+=gloTrg->PairUnlikeApt();
523       LSLowpt+=gloTrg->PairLikeLpt(); 
524       LSHighpt+=gloTrg->PairLikeHpt();
525       LSAllpt+=gloTrg->PairLikeApt();
526
527       if (PRINTOUT) {
528         printf("===================================================\n");
529         printf(" Global Trigger output       Low pt  High pt   All\n");
530         printf(" number of Single Plus      :\t");
531         printf("%i\t%i\t%i\t",gloTrg->SinglePlusLpt(),
532                gloTrg->SinglePlusHpt(),gloTrg->SinglePlusApt());
533         printf("\n");
534
535         printf(" number of Single Minus     :\t");
536         printf("%i\t%i\t%i\t",gloTrg->SingleMinusLpt(),
537                gloTrg->SingleMinusHpt(),gloTrg->SingleMinusApt());
538         printf("\n");
539
540         printf(" number of Single Undefined :\t"); 
541         printf("%i\t%i\t%i\t",gloTrg->SingleUndefLpt(),
542                gloTrg->SingleUndefHpt(),gloTrg->SingleUndefApt());
543         printf("\n");
544         
545         printf(" number of UnlikeSign pair  :\t"); 
546         printf("%i\t%i\t%i\t",gloTrg->PairUnlikeLpt(),
547                gloTrg->PairUnlikeHpt(),gloTrg->PairUnlikeApt());
548         printf("\n");
549         
550         printf(" number of LikeSign pair    :\t");  
551         printf("%i\t%i\t%i\t",gloTrg->PairLikeLpt(),
552                gloTrg->PairLikeHpt(),gloTrg->PairLikeApt());
553         printf("\n");
554         
555         printf("===================================================\n");
556       }
557       
558     } // end of loop on Global Trigger
559     
560     for (Int_t ilocal=0; ilocal<nlocals; ilocal++) { // Local Trigger
561       if (PRINTOUT) cout << " >>> Output for Local Trigger " << ilocal << "\n";
562       
563       locTrg = static_cast<AliMUONLocalTrigger*>(localTrigger->At(ilocal));
564       
565       if (PRINTOUT){
566         cout << "Circuit StripX Dev StripY: " 
567              << locTrg->LoCircuit() << " " 
568              << locTrg->LoStripX() << " " 
569              << locTrg->LoDev() << " " 
570              << locTrg->LoStripY() 
571              << "\n";
572         cout << "Lpt Hpt Apt: "     
573              << locTrg->LoLpt() << " "   
574              << locTrg->LoHpt() << " "  
575              << locTrg->LoApt() << "\n";
576
577       }
578         TgtupleLoc->Fill(locTrg->LoCircuit(),locTrg->LoStripX(),locTrg->LoDev(),locTrg->LoStripY(),locTrg->LoLpt(),locTrg->LoHpt(),locTrg->LoApt());
579     } // end of loop on Local Trigger
580
581
582     // fill ntuple
583     //TNtuple *Tgtuple = new TNtuple("Tgtuple","Trigger Ntuple","ev:global:spapt:smapt:undefapt:uplpt:uphpt:upapt:lpapt");
584        TgtupleGlo->Fill(ievent,nglobals,gloTrg->SinglePlusApt(),gloTrg->SingleMinusApt(),gloTrg->SingleUndefApt(),gloTrg->PairUnlikeLpt(),gloTrg->PairUnlikeHpt(),gloTrg->PairUnlikeApt(),gloTrg->PairLikeApt());
585
586     muondata.ResetTrigger();
587     if (event2Check!=0) ievent=nevents;
588   } // end loop on event  
589   
590   // Print out summary if loop ran over all event
591   if (!event2Check){
592     printf("\n");
593     printf("===================================================\n");
594     printf("===================  SUMMARY  =====================\n");
595     printf("\n");
596     printf("Total number of events processed %d \n", (event2Check==0) ? nevents : 1);
597     printf("\n");
598     printf(" Global Trigger output       Low pt  High pt   All\n");
599     printf(" number of Single Plus      :\t");
600     printf("%i\t%i\t%i\t",SPLowpt,SPHighpt,SPAllpt);
601     printf("\n");
602     printf(" number of Single Minus     :\t");
603     printf("%i\t%i\t%i\t",SMLowpt,SMHighpt,SMAllpt);
604     printf("\n");
605     printf(" number of Single Undefined :\t"); 
606     printf("%i\t%i\t%i\t",SULowpt,SUHighpt,SUAllpt);
607     printf("\n");
608     printf(" number of UnlikeSign pair  :\t"); 
609     printf("%i\t%i\t%i\t",USLowpt,USHighpt,USAllpt);
610     printf("\n");
611     printf(" number of LikeSign pair    :\t");  
612     printf("%i\t%i\t%i\t",LSLowpt,LSHighpt, LSAllpt);
613     printf("\n");
614     printf("===================================================\n");
615   }
616   
617   if (WRITE){
618     TFile *myFile = new TFile("TriggerCheck.root", "RECREATE");  
619     TgtupleGlo->Write();
620     TgtupleLoc->Write();
621     myFile->Close();
622   }
623
624
625   MUONLoader->UnloadRecPoints();
626 }
627
628
629 void MUONrectracks (Int_t event2Check=0, char * filename="galice.root"){
630 // reads and dumps trigger objects from MUON.RecPoints.root
631   TClonesArray * RecTracks;
632   
633   // Creating Run Loader and openning file containing Hits
634   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
635   if (RunLoader ==0x0) {
636     printf(">>> Error : Error Opening %s file \n",filename);
637     return;
638   }
639   
640   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
641   MUONLoader->LoadTracks("READ");
642   // Creating MUON data container
643   AliMUONData muondata(MUONLoader,"MUON","MUON");
644   
645     Int_t ievent, nevents;
646   nevents = RunLoader->GetNumberOfEvents();
647   
648   //  AliMUONTrack * rectrack;
649   
650   for (ievent=0; ievent<nevents; ievent++) {
651     if (event2Check!=0) ievent=event2Check;
652     RunLoader->GetEvent(ievent);
653     
654     muondata.SetTreeAddress("RT");
655     muondata.GetRecTracks();
656     RecTracks = muondata.RecTracks();
657     
658     
659     Int_t nrectracks = (Int_t) RecTracks->GetEntriesFast(); //
660
661     printf(">>> Event %d, Number of Recconstructed tracks %d \n",ievent, nrectracks);
662     // loop over tracks
663  
664  
665     Int_t nTrackHits;// nPrimary;
666     Double_t fitFmin;
667     Double_t bendingSlope, nonBendingSlope, inverseBendingMomentum;
668     Double_t xRec, yRec, zRec, chi2MatchTrigger;
669     Bool_t matchTrigger;
670     Double_t Pz,Px,Py,Pt,Ptot,Eta ;
671
672   // setting pointer for tracks, triggertracks & trackparam at vertex
673     AliMUONTrack* recTrack = 0;
674     AliMUONTrackParam* trackParam = 0;
675
676     for (Int_t iRecTracks = 0; iRecTracks <  nrectracks;  iRecTracks++) {
677     // reading info from tracks
678       recTrack = (AliMUONTrack*) RecTracks->At(iRecTracks);
679       trackParam = (AliMUONTrackParam*) (recTrack->GetTrackParamAtHit())->First();
680       trackParam->ExtrapToZ(0.0);
681       bendingSlope            = trackParam->GetBendingSlope();
682       nonBendingSlope         = trackParam->GetNonBendingSlope();
683       inverseBendingMomentum = trackParam->GetInverseBendingMomentum();
684       xRec  = trackParam->GetNonBendingCoor();
685       yRec  = trackParam->GetBendingCoor();
686       zRec  = trackParam->GetZ();
687
688       nTrackHits       = recTrack->GetNTrackHits();
689       fitFmin          = recTrack->GetFitFMin();
690       matchTrigger     = recTrack->GetMatchTrigger();
691       chi2MatchTrigger = recTrack->GetChi2MatchTrigger();
692       
693       Px = trackParam->Px();
694       Py = trackParam->Py(); 
695       Pz = trackParam->Pz(); 
696       Pt = TMath::Sqrt(Px*Px + Py*Py );
697       Ptot = TMath::Sqrt(Px*Px + Py*Py + Pz*Pz);
698       Eta =  (Pt!=0) ? 0.5*log( (Ptot+Pz)/(Ptot-Pz) ) : 999999999.999 ;
699        
700       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);
701       printf("    Px=%f  Py =%f  Pz =%f   Pt=%f  Ptot=%f   PseudoRap=%f  \n",Px,Py,Pz,Pt,Ptot,Eta);
702     } // end loop tracks
703
704     muondata.ResetRecTracks();
705     if (event2Check!=0) ievent=nevents;
706   } // end loop on event  
707   MUONLoader->UnloadTracks();
708 }
709
710
711
712
713
714
715
716
717
718