]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONCheck.C
MUONrectrigger improved (Christophe)
[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 void MUONdigits(Int_t event2Check=0, char * filename="galice.root")
142 {
143   // Creating Run Loader and openning file containing Hits
144   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
145   if (RunLoader ==0x0) {
146     printf(">>> Error : Error Opening %s file \n",filename);
147     return;
148   }
149   // Loading MUON subsystem
150   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
151   MUONLoader->LoadDigits("READ");
152   // Creating MUON data container
153   AliMUONData muondata(MUONLoader,"MUON","MUON");
154
155   Int_t ievent, nevents;
156   nevents = RunLoader->GetNumberOfEvents();
157   AliMUONDigit * mDigit;
158   
159   for(ievent=0; ievent<nevents; ievent++) {
160     if (event2Check!=0) ievent=event2Check;
161     printf(">>> Event %d \n",ievent);
162     RunLoader->GetEvent(ievent);
163   
164     // Addressing
165     Int_t ichamber, nchambers;
166     nchambers = AliMUONConstants::NCh(); ;
167     muondata.SetTreeAddress("D,GLT");
168     //    char branchname[30];    
169  
170     muondata.GetDigits();
171     // Loop on chambers
172     for( ichamber=0; ichamber<nchambers; ichamber++) {
173       Int_t idigit, ndigits;
174       ndigits = (Int_t) muondata.Digits(ichamber)->GetEntriesFast();
175       for(idigit=0; idigit<ndigits; idigit++) {
176         mDigit = static_cast<AliMUONDigit*>(muondata.Digits(ichamber)->At(idigit));
177         Int_t PadX   = mDigit->PadX();     // Pad X number
178         Int_t PadY   = mDigit->PadY();     // Pad Y number
179         Int_t Signal = mDigit->Signal();   // Physics Signal
180         Int_t Physics= mDigit->Physics();  // Physics contribution to signal
181         //        Int_t Hit    = mDigit->Hit();      // iHit
182         Int_t Cathode= mDigit->Cathode();  // Cathode
183         Int_t Track0 = mDigit->Track(0);
184         Int_t Track1 = mDigit->Track(1); 
185         //Int_t Track2 = mDigit->Track(2);
186         Int_t TCharges0 = mDigit->TrackCharge(0);  //charge per track making this digit (up to 10)
187         Int_t TCharges1 = mDigit->TrackCharge(1);
188         //Int_t TCharges2 = mDigit->TrackCharge(2);
189         Int_t idDE = mDigit->DetElemId();
190         //        printf(">>> Cathode %d\n",Cathode);
191         
192         printf(">>> DetEle %4d Digit%4d Cath %1d (Ix,Iy)=(%3d,%3d) Signal=%4d Physics=%4d Track0=%4d Charge0=%4d Track1=%4d Charge1=%4d \n",
193                idDE, idigit, Cathode, PadX, PadY, Signal, Physics, Track0, 
194                TCharges0, Track1, TCharges1);
195       } // end digit loop
196     } // end chamber loop
197     muondata.ResetDigits();
198     //    } // end cathode loop
199     if (event2Check!=0) ievent=nevents;
200   }  // end event loop
201   MUONLoader->UnloadDigits();
202 }
203
204 void MUONoccupancy(Int_t event2Check=0,  Bool_t perDetEle =kFALSE, char * filename="galice.root") {
205   // Creating Run Loader and openning file containing Hits
206   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
207   if (RunLoader ==0x0) {
208     printf(">>> Error : Error Opening %s file \n",filename);
209     return;
210   }
211   // Loading MUON subsystem
212   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
213   MUONLoader->LoadDigits("READ");
214   // Creating MUON data container
215   AliMUONData muondata(MUONLoader,"MUON","MUON");
216   
217   Int_t ievent, nevents;
218   nevents = RunLoader->GetNumberOfEvents();
219   AliMUONDigit * mDigit =0x0;
220   AliMpVSegmentation * segbend = 0x0;
221   AliMpVSegmentation * segnonbend = 0x0;
222   AliMpIntPair pad(0,0);
223
224   Int_t dEoccupancy_bending[14][26];
225   Int_t dEoccupancy_nonbending[14][26];
226   Int_t cHoccupancy_bending[14];
227   Int_t cHoccupancy_nonbending[14];
228   Int_t totaloccupancy_bending =0;
229   Int_t totaloccupancy_nonbending =0;
230
231   Int_t dEchannels_bending[14][26];
232   Int_t dEchannels_nonbending[14][26];
233   Int_t cHchannels_bending[14];
234   Int_t cHchannels_nonbending[14];
235   Int_t totalchannels_bending =0;
236   Int_t totalchannels_nonbending =0;
237
238   Int_t ichamber, nchambers,idetele, detele, ix, iy;
239   nchambers = AliMUONConstants::NCh(); ;
240
241   AliMpSegFactory factory;
242
243   for (ichamber=0; ichamber<nchambers; ichamber++) {
244     cHchannels_bending[ichamber]=0;
245     cHchannels_nonbending[ichamber]=0;
246     for (idetele=0; idetele<26; idetele++) {
247       detele= 100*(ichamber +1)+idetele;
248       dEchannels_bending[ichamber][idetele]=0;
249       dEchannels_nonbending[ichamber][idetele]=0;
250       dEoccupancy_bending[ichamber][idetele]=0;
251       dEoccupancy_nonbending[ichamber][idetele]=0;
252       if ( AliMpDEManager::IsValidDetElemId(detele) ) {
253         
254         segbend    =  factory.CreateMpSegmentation(detele, 0);
255         segnonbend =  factory.CreateMpSegmentation(detele, 1);
256         if (AliMpDEManager::GetPlaneType(detele, 0) != kBendingPlane ) {
257           AliMpVSegmentation* tmp = segbend;
258           segbend    =  segnonbend;
259           segnonbend =  tmp;
260         }  
261           
262         for(ix=0; ix<=segbend->MaxPadIndexX(); ix++) {
263           for(iy=0; iy<=segbend->MaxPadIndexY(); iy++) {
264             pad.SetFirst(ix);
265             pad.SetSecond(iy);
266             if( segbend->HasPad(pad) )   {  
267               dEchannels_bending[ichamber][idetele]++;
268               cHchannels_bending[ichamber]++;
269               totalchannels_bending++;
270             }
271           }
272         }
273         for(ix=0; ix<=segnonbend->MaxPadIndexX(); ix++) {
274           for(iy=0; iy<=segnonbend->MaxPadIndexY(); iy++) {
275             pad.SetFirst(ix);
276             pad.SetSecond(iy);
277             if(segnonbend->HasPad(pad))  {
278               dEchannels_nonbending[ichamber][idetele]++;  
279               cHchannels_nonbending[ichamber]++;
280               totalchannels_nonbending++;
281             }
282           }
283         }
284         if (perDetEle) printf(">>> Detection element %4d has %5d channels in bending and %5d channels in nonbending \n",
285              detele, dEchannels_bending[ichamber][idetele], dEchannels_nonbending[ichamber][idetele] ); 
286       }
287     }
288     printf(">>> Chamber %2d has %6d channels in bending and %6d channels in nonbending \n",
289            ichamber+1,  cHchannels_bending[ichamber], cHchannels_nonbending[ichamber]);
290   }
291   printf(">>Spectrometer has  %7d channels in bending and %7d channels in nonbending \n",
292          totalchannels_bending, totalchannels_nonbending);
293
294   factory.DeleteSegmentations();
295
296   ievent=event2Check;
297   printf(">>> Event %d \n",ievent);
298   RunLoader->GetEvent(ievent);
299     
300   // Addressing
301   muondata.SetTreeAddress("D"); 
302   muondata.GetDigits();
303   // Loop on chambers
304   for( ichamber=0; ichamber<nchambers; ichamber++) {
305     cHoccupancy_bending[ichamber]   = 0;
306     cHoccupancy_nonbending[ichamber]= 0;
307     Int_t idigit, ndigits;
308     ndigits = (Int_t) muondata.Digits(ichamber)->GetEntriesFast();
309     for(idigit=0; idigit<ndigits; idigit++) {
310       mDigit = static_cast<AliMUONDigit*>(muondata.Digits(ichamber)->At(idigit));
311       Int_t detele = mDigit->DetElemId();
312       Int_t idetele = detele-(ichamber+1)*100;
313       if ( mDigit->Cathode() == 0 ) {
314
315         cHoccupancy_bending[ichamber]++;
316         dEoccupancy_bending[ichamber][idetele]++;
317         totaloccupancy_bending++;
318       }
319       else {
320         cHoccupancy_nonbending[ichamber]++;
321         dEoccupancy_nonbending[ichamber][idetele]++;
322         totaloccupancy_nonbending++;
323       }
324     } // end digit loop    
325
326     printf(">>> Chamber %2d  nChannels Bending %5d  nChannels NonBending %5d \n", 
327            ichamber+1, 
328            cHoccupancy_bending[ichamber],
329            cHoccupancy_nonbending[ichamber]);           
330     printf(">>> Chamber %2d  Occupancy Bending %5.2f \%  Occupancy NonBending %5.2f \% \n", 
331            ichamber+1, 
332            100.*((Float_t) cHoccupancy_bending[ichamber])/((Float_t) cHchannels_bending[ichamber]),
333            100.*((Float_t) cHoccupancy_nonbending[ichamber])/((Float_t) cHchannels_bending[ichamber])            );
334
335
336     for(Int_t idetele=0; idetele<26; idetele++) {
337       Int_t detele = idetele + 100*(ichamber+1);
338       if ( AliMpDEManager::IsValidDetElemId(detele) ) {
339         if (perDetEle) {
340           printf(">>> DetEle %4d nChannels Bending %5d  nChannels NonBending %5d \n", 
341                  idetele+100*(ichamber+1), 
342                  dEoccupancy_bending[ichamber][idetele],
343                  dEoccupancy_nonbending[ichamber][idetele]);  
344           printf(">>> DetEle %4d Occupancy Bending %5.2f \%  Occupancy NonBending %5.2f \% \n", 
345                  idetele+100*(ichamber+1), 
346                  100.*((Float_t) dEoccupancy_bending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]),
347                  100.*((Float_t) dEoccupancy_nonbending[ichamber][idetele])/((Float_t) dEchannels_bending[ichamber][idetele]));  
348         }
349       }
350     }
351   } // end chamber loop
352   printf(">>> Muon Spectrometer  Occupancy Bending %5.2f\%  Occupancy NonBending %5.2f\% \n",  
353            100.*((Float_t) totaloccupancy_bending)/((Float_t) totalchannels_bending),
354          100.*((Float_t) totaloccupancy_nonbending)/((Float_t) totalchannels_nonbending)            );
355   muondata.ResetDigits();
356   //    } // end cathode loop
357   MUONLoader->UnloadDigits();
358 }
359
360 void MUONrecpoints(Int_t event2Check=0, char * filename="galice.root") {
361
362   // Creating Run Loader and openning file containing Hits
363   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
364   if (RunLoader ==0x0) {
365     printf(">>> Error : Error Opening %s file \n",filename);
366     return;
367   }
368   // Getting MUONloader
369   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
370   MUONLoader->LoadRecPoints("READ");
371   // Creating MUON data container
372   AliMUONData muondata(MUONLoader,"MUON","MUON");
373
374   Int_t ievent, nevents;
375   nevents = RunLoader->GetNumberOfEvents();
376   AliMUONRawCluster * mRecPoint = 0;
377   
378   for(ievent=0; ievent<nevents; ievent++) {
379     if (event2Check!=0) ievent=event2Check;
380     printf(">>> Event %d \n",ievent);
381     RunLoader->GetEvent(ievent);
382     // Addressing
383     Int_t ichamber, nchambers;
384     nchambers = AliMUONConstants::NTrackingCh();
385     muondata.SetTreeAddress("RC,TC"); 
386     char branchname[30];    
387     muondata.GetRawClusters();
388     // Loop on chambers
389     for( ichamber=0; ichamber<nchambers; ichamber++) {
390       sprintf(branchname,"MUONRawClusters%d",ichamber+1);
391       //printf(">>>  branchname %s\n",branchname);
392       Int_t irecpoint, nrecpoints;
393       nrecpoints = (Int_t) muondata.RawClusters(ichamber)->GetEntriesFast();
394       // printf(">>> Chamber %2d, Number of recpoints = %6d \n",ichamber+1, nrecpoints);
395       for(irecpoint=0; irecpoint<nrecpoints; irecpoint++) {
396         mRecPoint = static_cast<AliMUONRawCluster*>(muondata.RawClusters(ichamber)->At(irecpoint));
397         Int_t Track0 = mRecPoint->GetTrack(0);
398         Int_t Track1 = mRecPoint->GetTrack(1); 
399         Int_t Track2 = mRecPoint->GetTrack(2);
400         Int_t Q0 = mRecPoint->GetCharge(0);
401         Int_t Q1 = mRecPoint->GetCharge(1);
402         Float_t x0 = mRecPoint->GetX(0);
403         Float_t x1 = mRecPoint->GetX(1);
404         Float_t y0 = mRecPoint->GetY(0);
405         Float_t y1 = mRecPoint->GetY(1);
406         Float_t z0 = mRecPoint->GetZ(0);
407         Float_t z1 = mRecPoint->GetZ(1);
408         Float_t chi2_0 =  mRecPoint->GetChi2(0);
409         //Float_t chi2_1 =  mRecPoint->GetChi2(1);
410         Int_t de = mRecPoint->GetDetElemId();
411         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",
412                irecpoint,de,x0,y0,z0,Q0,Q1,Track0, Track1, Track2, chi2_0);
413         if( (x0!=x1) || (y0!=y1) || (z0!=z1) )
414           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); 
415       } // end recpoint loop
416     } // end chamber loop
417     muondata.ResetRawClusters();
418     if (event2Check!=0) ievent=nevents;
419   }  // end event loop
420   MUONLoader->UnloadRecPoints();
421 }
422
423
424
425 void MUONrectrigger (Int_t event2Check=0, char * filename="galice.root"){
426   // reads and dumps trigger objects from MUON.RecPoints.root
427   TClonesArray * globalTrigger;
428   TClonesArray * localTrigger;
429   
430   // Do NOT print out all the info if the loop runs over all events 
431   Int_t PRINTOUT = (event2Check == 0 ) ? 0 : 1 ;  
432   
433   // Book a ntuple for more detailled studies
434   TNtuple *Tgtuple = new TNtuple("Tgtuple","Trigger Ntuple","ev:global:spapt:smapt:undefapt:uplpt:uphpt:upapt:suapt");
435   Int_t WRITE = 0;
436
437   // counters
438   Int_t SPLowpt=0,SPHighpt=0,SPAllpt=0;
439   Int_t SMLowpt=0,SMHighpt=0,SMAllpt=0;
440   Int_t SULowpt=0,SUHighpt=0,SUAllpt=0;
441   Int_t USLowpt=0,USHighpt=0,USAllpt=0;
442   Int_t LSLowpt=0,LSHighpt=0,LSAllpt=0;
443
444   // Creating Run Loader and openning file containing Hits
445   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
446   if (RunLoader ==0x0) {
447     printf(">>> Error : Error Opening %s file \n",filename);
448     return;
449   }
450   
451   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
452   MUONLoader->LoadDigits("READ");
453   // Creating MUON data container
454   AliMUONData muondata(MUONLoader,"MUON","MUON");
455   
456   
457   Int_t ievent, nevents;
458   nevents = RunLoader->GetNumberOfEvents();
459   
460   AliMUONGlobalTrigger *gloTrg;
461   AliMUONLocalTrigger *locTrg;
462   
463   for (ievent=0; ievent<nevents; ievent++) {
464     if (event2Check!=0) ievent=event2Check;
465     if (ievent%100==0 || event2Check) cout << "Processing event " << ievent << endl;
466     RunLoader->GetEvent(ievent);
467     
468     muondata.SetTreeAddress("GLT"); 
469     muondata.GetTriggerD();
470     
471     globalTrigger = muondata.GlobalTrigger();
472     localTrigger = muondata.LocalTrigger();
473     
474     Int_t nglobals = (Int_t) globalTrigger->GetEntriesFast(); // should be 1
475     Int_t nlocals  = (Int_t) localTrigger->GetEntriesFast(); // up to 234
476     if (PRINTOUT) printf("###################################################\n");
477     if (PRINTOUT) {cout << " event " << ievent 
478                         << " nglobals nlocals: " << nglobals << " " << nlocals << "\n"; }
479     
480     for (Int_t iglobal=0; iglobal<nglobals; iglobal++) { // Global Trigger
481       gloTrg = static_cast<AliMUONGlobalTrigger*>(globalTrigger->At(iglobal));
482       
483       SPLowpt+=gloTrg->SinglePlusLpt() ;
484       SPHighpt+=gloTrg->SinglePlusHpt() ;
485       SPAllpt+=gloTrg->SinglePlusApt() ;
486       SMLowpt+=gloTrg->SingleMinusLpt();
487       SMHighpt+=gloTrg->SingleMinusHpt();
488       SMAllpt+=gloTrg->SingleMinusApt();
489       SULowpt+=gloTrg->SingleUndefLpt();
490       SUHighpt+=gloTrg->SingleUndefHpt();
491       SUAllpt+=gloTrg->SingleUndefApt();
492       USLowpt+=gloTrg->PairUnlikeLpt(); 
493       USHighpt+=gloTrg->PairUnlikeHpt();
494       USAllpt+=gloTrg->PairUnlikeApt();
495       LSLowpt+=gloTrg->PairLikeLpt(); 
496       LSHighpt+=gloTrg->PairLikeHpt();
497       LSAllpt+=gloTrg->PairLikeApt();
498
499       if (PRINTOUT) {
500         printf("===================================================\n");
501         printf(" Global Trigger output       Low pt  High pt   All\n");
502         printf(" number of Single Plus      :\t");
503         printf("%i\t%i\t%i\t",gloTrg->SinglePlusLpt(),
504                gloTrg->SinglePlusHpt(),gloTrg->SinglePlusApt());
505         printf("\n");
506
507         printf(" number of Single Minus     :\t");
508         printf("%i\t%i\t%i\t",gloTrg->SingleMinusLpt(),
509                gloTrg->SingleMinusHpt(),gloTrg->SingleMinusApt());
510         printf("\n");
511
512         printf(" number of Single Undefined :\t"); 
513         printf("%i\t%i\t%i\t",gloTrg->SingleUndefLpt(),
514                gloTrg->SingleUndefHpt(),gloTrg->SingleUndefApt());
515         printf("\n");
516         
517         printf(" number of UnlikeSign pair  :\t"); 
518         printf("%i\t%i\t%i\t",gloTrg->PairUnlikeLpt(),
519                gloTrg->PairUnlikeHpt(),gloTrg->PairUnlikeApt());
520         printf("\n");
521         
522         printf(" number of LikeSign pair    :\t");  
523         printf("%i\t%i\t%i\t",gloTrg->PairLikeLpt(),
524                gloTrg->PairLikeHpt(),gloTrg->PairLikeApt());
525         printf("\n");
526         
527         printf("===================================================\n");
528       }
529       
530     } // end of loop on Global Trigger
531     
532     for (Int_t ilocal=0; ilocal<nlocals; ilocal++) { // Local Trigger
533       if (PRINTOUT) cout << " >>> Output for Local Trigger " << ilocal << "\n";
534       
535       locTrg = static_cast<AliMUONLocalTrigger*>(localTrigger->At(ilocal));
536       
537       if (PRINTOUT){
538         cout << "Circuit StripX Dev StripY: " 
539              << locTrg->LoCircuit() << " " 
540              << locTrg->LoStripX() << " " 
541              << locTrg->LoDev() << " " 
542              << locTrg->LoStripY() 
543              << "\n";
544         cout << "Lpt Hpt Apt: "     
545              << locTrg->LoLpt() << " "   
546              << locTrg->LoHpt() << " "  
547              << locTrg->LoApt() << "\n";
548       }
549     } // end of loop on Local Trigger
550
551
552     // fill ntuple
553     //TNtuple *Tgtuple = new TNtuple("Tgtuple","Trigger Ntuple","ev:global:spapt:smapt:undefapt:uplpt:uphpt:upapt:suapt");
554        Tgtuple->Fill(ievent,nglobals,gloTrg->SinglePlusApt(),gloTrg->SingleMinusApt(),gloTrg->SingleUndefApt(),gloTrg->PairUnlikeLpt(),gloTrg->PairUnlikeHpt(),gloTrg->PairUnlikeApt(),gloTrg->SingleUndefApt());
555
556
557     muondata.ResetTrigger();
558     if (event2Check!=0) ievent=nevents;
559   } // end loop on event  
560   
561   // Print out summary if loop ran over all event
562   if (!event2Check){
563     printf("\n");
564     printf("===================================================\n");
565     printf("===================  SUMMARY  =====================\n");
566     printf("\n");
567     printf("Total number of events processed %d \n", (event2Check==0) ? nevents : 1);
568     printf("\n");
569     printf(" Global Trigger output       Low pt  High pt   All\n");
570     printf(" number of Single Plus      :\t");
571     printf("%i\t%i\t%i\t",SPLowpt,SPHighpt,SPAllpt);
572     printf("\n");
573     printf(" number of Single Minus     :\t");
574     printf("%i\t%i\t%i\t",SMLowpt,SMHighpt,SMAllpt);
575     printf("\n");
576     printf(" number of Single Undefined :\t"); 
577     printf("%i\t%i\t%i\t",SULowpt,SUHighpt,SUAllpt);
578     printf("\n");
579     printf(" number of UnlikeSign pair  :\t"); 
580     printf("%i\t%i\t%i\t",USLowpt,USHighpt,USAllpt);
581     printf("\n");
582     printf(" number of LikeSign pair    :\t");  
583     printf("%i\t%i\t%i\t",LSLowpt,LSHighpt, LSAllpt);
584     printf("\n");
585     printf("===================================================\n");
586   }
587   
588   if (WRITE){
589     TFile *myFile = new TFile("TriggerCheck.root", "RECREATE");
590     Tgtuple->Write();
591     myFile->Close();
592   }
593
594
595   MUONLoader->UnloadRecPoints();
596 }
597
598 void MUONrectracks (Int_t event2Check=0, char * filename="galice.root"){
599 // reads and dumps trigger objects from MUON.RecPoints.root
600   TClonesArray * RecTracks;
601   
602   // Creating Run Loader and openning file containing Hits
603   AliRunLoader * RunLoader = AliRunLoader::Open(filename,"MUONFolder","READ");
604   if (RunLoader ==0x0) {
605     printf(">>> Error : Error Opening %s file \n",filename);
606     return;
607   }
608   
609   AliLoader * MUONLoader = RunLoader->GetLoader("MUONLoader");
610   MUONLoader->LoadTracks("READ");
611   // Creating MUON data container
612   AliMUONData muondata(MUONLoader,"MUON","MUON");
613   
614     Int_t ievent, nevents;
615   nevents = RunLoader->GetNumberOfEvents();
616   
617   //  AliMUONTrack * rectrack;
618   
619   for (ievent=0; ievent<nevents; ievent++) {
620     if (event2Check!=0) ievent=event2Check;
621     RunLoader->GetEvent(ievent);
622     
623     muondata.SetTreeAddress("RT");
624     muondata.GetRecTracks();
625     RecTracks = muondata.RecTracks();
626     
627     
628     Int_t nrectracks = (Int_t) RecTracks->GetEntriesFast(); //
629
630     printf(">>> Event %d, Number of Recconstructed tracks %d \n",ievent, nrectracks);
631     // loop over tracks
632  
633  
634     Int_t nTrackHits;// nPrimary;
635     Double_t fitFmin;
636     Double_t bendingSlope, nonBendingSlope, inverseBendingMomentum;
637     Double_t xRec, yRec, zRec, chi2MatchTrigger;
638     Bool_t matchTrigger;
639
640   // setting pointer for tracks, triggertracks & trackparam at vertex
641     AliMUONTrack* recTrack = 0;
642     AliMUONTrackParam* trackParam = 0;
643
644     for (Int_t iRecTracks = 0; iRecTracks <  nrectracks;  iRecTracks++) {
645     // reading info from tracks
646       recTrack = (AliMUONTrack*) RecTracks->At(iRecTracks);
647       trackParam = (AliMUONTrackParam*) (recTrack->GetTrackParamAtHit())->First();
648       trackParam->ExtrapToZ(0.0);
649       bendingSlope            = trackParam->GetBendingSlope();
650       nonBendingSlope         = trackParam->GetNonBendingSlope();
651       inverseBendingMomentum = trackParam->GetInverseBendingMomentum();
652       xRec  = trackParam->GetNonBendingCoor();
653       yRec  = trackParam->GetBendingCoor();
654       zRec  = trackParam->GetZ();
655
656       nTrackHits       = recTrack->GetNTrackHits();
657       fitFmin          = recTrack->GetFitFMin();
658       matchTrigger     = recTrack->GetMatchTrigger();
659       chi2MatchTrigger = recTrack->GetChi2MatchTrigger();
660       
661       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);
662     } // end loop tracks
663
664     muondata.ResetRecTracks();
665     if (event2Check!=0) ievent=nevents;
666   } // end loop on event  
667   MUONLoader->UnloadTracks();
668 }
669
670
671
672
673
674
675
676
677
678