]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Updated comments for Doxygen
[u/mrichter/AliRoot.git] / MUON / MUONRecoCheck.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 // ROOT includes
17 #include "TClonesArray.h"
18 #include "TH1.h"
19 #include "TParticle.h"
20 #include "TFile.h"
21
22 // STEER includes
23 #include "AliRun.h"
24 #include "AliHeader.h"
25 #include "AliMC.h"
26 #include "AliStack.h"
27 #include "AliRunLoader.h"
28
29 // MUON includes
30 #include "AliMUON.h"
31 #include "AliMUONConstants.h"
32 #include "AliMUONTrack.h"
33 #include "AliMUONRecoCheck.h"
34 #include "AliMUONTrackParam.h"
35
36 Int_t TrackCheck( Bool_t *compTrack);
37
38 void MUONRecoCheck (Int_t nEvent = 1, char * filename="galice.root"){
39
40   // Utility macro to check the muon reconstruction. Reconstructed tracks are compared
41   // to reference tracks. The reference tracks are built from AliTrackReference for the
42   // hit in chamber (0..9) and from kinematics (TreeK) for the vertex parameters.     
43   
44   Int_t nTrackReco, nTrackRef;
45   AliMUONTrack *trackReco, *trackRef;
46   Bool_t *compTrack;
47   Bool_t compTrackOK[10];
48   Int_t nHitOK = 0;
49   Int_t testTrack = 0;  
50   Int_t iTrack = 0;
51   Int_t indexOK = 0;
52   Int_t trackID = 0;
53   Double_t sigma2Cut = 16;  // 4 sigmas cut, sigma2Cut = 4*4
54   AliMUONTrackParam *trackParam;
55   TClonesArray *trackParamAtHit;
56   Double_t x1,y1,z1,pX1,pY1,pZ1,p1;
57   Double_t x2,y2,z2,pX2,pY2,pZ2,p2;
58   TParticle* particle = new TParticle();
59
60   TClonesArray *trackRecoArray = NULL;
61   TClonesArray *trackRefArray = NULL;
62
63
64   // File for histograms and histogram booking
65   TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
66
67   TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5);
68   TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
69   TH1F *hNHitComp = new TH1F("hNHitComp"," Nb of compatible hits / track ",15,-0.5,14.5);
70   TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5);
71   TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
72   
73   TH1F *hResMomVertex = new TH1F("hMomVertex"," delta P vertex (GeV/c)",100,-10.,10);
74   TH1F *hResMomFirstHit = new TH1F("hMomFirstHit"," delta P first hit (GeV/c)",100,-10.,10);
75
76
77   
78   AliMUONRecoCheck rc(filename);
79   AliRunLoader *runLoader = rc.GetRunLoader();
80
81     
82   Int_t nevents = runLoader->GetNumberOfEvents();
83   
84   if (nevents < nEvent) nEvent = nevents;
85   
86   Int_t ievent;
87   Int_t nReconstructibleTracks = 0;
88   Int_t nReconstructedTracks = 0;
89   Int_t nReconstructibleTracksCheck = 0;
90
91   for (ievent=0; ievent<nEvent; ievent++) {
92     if (!(ievent%10)) printf(" **** event # %d  \n",ievent);
93     runLoader->GetEvent(ievent);
94     rc.ResetTracks();
95     rc.MakeTrackRef(); // make reconstructible tracks
96 //     rc.PrintEvent();
97
98     
99     trackRecoArray = rc.GetTrackReco();
100     trackRefArray = rc.GetMuonTrackRef();
101     
102     nTrackRef = trackRefArray->GetEntriesFast();
103     nTrackReco = trackRecoArray->GetEntriesFast();
104
105     hReconstructible->Fill(rc.GetNumberOfReconstuctibleTracks());
106     hReco->Fill(rc.GetNumberOfRecoTracks());
107
108     nReconstructibleTracks += rc.GetNumberOfReconstuctibleTracks();
109     nReconstructedTracks += rc.GetNumberOfRecoTracks();
110
111     for (Int_t index1 = 0; index1 < nTrackRef; index1++) {  
112       trackRef = (AliMUONTrack *)trackRefArray->At(index1);
113
114       testTrack = 0;
115       indexOK = 0;
116       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
117         compTrackOK[ch] = kFALSE;
118       for (Int_t index2 = 0; index2 < nTrackReco; index2++) {
119         trackReco = (AliMUONTrack *)trackRecoArray->At(index2);
120
121         // check if trackRef is compatible with trackReco
122         compTrack = trackRef->CompatibleTrack(trackReco,sigma2Cut);
123
124         iTrack = TrackCheck(compTrack);
125         
126         if (iTrack > testTrack) {
127           nHitOK = 0;
128           for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
129             if (compTrack[ch]) nHitOK++;
130             compTrackOK[ch] = compTrack[ch];
131           }
132           testTrack = iTrack;
133           indexOK = index2;
134         }
135       }
136
137       hTestTrack->Fill(testTrack);
138       trackID = trackRef->GetTrackID();
139       hTrackRefID->Fill(trackID);
140
141       if (testTrack == 4) {     // tracking requirements verified, track is found
142         nReconstructibleTracksCheck++;
143         hNHitComp->Fill(nHitOK);
144         particle = runLoader->GetHeader()->Stack()->Particle(trackID);
145 //      printf(" trackID: %d , PDG code: %d \n",trackID,particle->GetPdgCode());
146         trackParam = trackRef->GetTrackParamAtVertex();
147         x1 = trackParam->GetNonBendingCoor();
148         y1 = trackParam->GetBendingCoor();
149         z1 = trackParam->GetZ();
150         pX1 = trackParam->Px();
151         pY1 = trackParam->Py();
152         pZ1 = trackParam->Pz();
153         p1  = trackParam->P();
154         
155 //      printf(" Ref. track at vertex: x,y,z: %f %f %f px,py,pz,p: %f %f %f %f \n",x1,y1,z1,pX1,pY1,pZ1,p1);
156         
157         trackParam = ((AliMUONTrack *)trackRecoArray->At(indexOK))->GetTrackParamAtVertex();
158         x2 = trackParam->GetNonBendingCoor();
159         y2 = trackParam->GetBendingCoor();
160         z2 = trackParam->GetZ();
161         pX2 = trackParam->Px();
162         pY2 = trackParam->Py();
163         pZ2 = trackParam->Pz();
164         p2  = trackParam->P();
165 //      printf(" Reconst. track at vertex: x,y,z: %f %f %f px,py,pz: %f %f %f %f \n",x2,y2,z2,pX2,pY2,pZ2,p2);
166         
167         hResMomVertex->Fill(p2-p1);
168
169         trackParamAtHit =  trackRef->GetTrackParamAtHit();
170         trackParam = (AliMUONTrackParam*) trackParamAtHit->First();
171         x1 = trackParam->GetNonBendingCoor();
172         y1 = trackParam->GetBendingCoor();
173         z1 = trackParam->GetZ();
174         pX1 = trackParam->Px();
175         pY1 = trackParam->Py();
176         pZ1 = trackParam->Pz();
177         p1  = trackParam->P();
178 //      printf(" Ref. track at 1st hit: x,y,z: %f %f %f px,py,pz: %f %f %f \n",x1,y1,z1,pX1,pY1,pZ1);
179         trackParamAtHit =  ((AliMUONTrack *) trackRecoArray->At(indexOK))->GetTrackParamAtHit();
180         trackParam = (AliMUONTrackParam*) trackParamAtHit->First();
181         x2 = trackParam->GetNonBendingCoor();
182         y2 = trackParam->GetBendingCoor();
183         z2 = trackParam->GetZ();
184         pX2 = trackParam->Px();
185         pY2 = trackParam->Py();
186         pZ2 = trackParam->Pz();
187         p2  = trackParam->P();
188 //      printf(" Reconst. track at 1st hit: x,y,z: %f %f %f px,py,pz: %f %f %f \n",x2,y2,z2,pX2,pY2,pZ2);
189
190         hResMomFirstHit->Fill(p2-p1);
191                
192       }
193     } // end loop track ref.
194
195   } // end loop on event  
196
197   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
198   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
199   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
200   
201   histoFile->Write();
202   histoFile->Close();
203 }
204
205
206 Int_t TrackCheck( Bool_t *compTrack)
207 {
208   // Apply reconstruction requirements
209   // Return number of validated conditions 
210   // If all the tests are verified then TrackCheck = 4 (good track)
211   Int_t iTrack = 0;
212   Int_t hitsInLastStations = 0;
213   
214   // apply reconstruction requirements
215   if (compTrack[0] || compTrack[1]) iTrack++; // at least one hit in st. 0
216   if (compTrack[2] || compTrack[3]) iTrack++; // at least one hit in st. 1
217   if (compTrack[4] || compTrack[5]) iTrack++; // at least one hit in st. 2
218   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
219     if (compTrack[ch]) hitsInLastStations++; 
220   }
221   if (hitsInLastStations > 2) iTrack++; // at least 3 hits in st. 3 & 4
222   
223   return iTrack;
224
225 }
226
227
228
229
230
231
232