]>
Commit | Line | Data |
---|---|---|
b8dc484b | 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" | |
a0bfad0e | 21 | #include <TGeoManager.h> |
b8dc484b | 22 | |
23 | // STEER includes | |
24 | #include "AliRun.h" | |
25 | #include "AliHeader.h" | |
26 | #include "AliMC.h" | |
27 | #include "AliStack.h" | |
28 | #include "AliRunLoader.h" | |
45a05f07 | 29 | #include "AliMagFMaps.h" |
b8dc484b | 30 | |
31 | // MUON includes | |
32 | #include "AliMUON.h" | |
33 | #include "AliMUONConstants.h" | |
34 | #include "AliMUONTrack.h" | |
35 | #include "AliMUONRecoCheck.h" | |
36 | #include "AliMUONTrackParam.h" | |
8cde4af5 | 37 | #include "AliMUONTrackExtrap.h" |
45a05f07 | 38 | #include "AliTracker.h" |
b8dc484b | 39 | |
40 | Int_t TrackCheck( Bool_t *compTrack); | |
41 | ||
a0bfad0e | 42 | void MUONRecoCheck (Int_t nEvent = 1, char* geoFilename = "geometry.root", char * filename="galice.root"){ |
b8dc484b | 43 | |
44 | // Utility macro to check the muon reconstruction. Reconstructed tracks are compared | |
45 | // to reference tracks. The reference tracks are built from AliTrackReference for the | |
46 | // hit in chamber (0..9) and from kinematics (TreeK) for the vertex parameters. | |
47 | ||
48 | Int_t nTrackReco, nTrackRef; | |
49 | AliMUONTrack *trackReco, *trackRef; | |
50 | Bool_t *compTrack; | |
51 | Bool_t compTrackOK[10]; | |
52 | Int_t nHitOK = 0; | |
53 | Int_t testTrack = 0; | |
54 | Int_t iTrack = 0; | |
55 | Int_t indexOK = 0; | |
56 | Int_t trackID = 0; | |
57 | Double_t sigma2Cut = 16; // 4 sigmas cut, sigma2Cut = 4*4 | |
58 | AliMUONTrackParam *trackParam; | |
59 | TClonesArray *trackParamAtHit; | |
60 | Double_t x1,y1,z1,pX1,pY1,pZ1,p1; | |
61 | Double_t x2,y2,z2,pX2,pY2,pZ2,p2; | |
62 | TParticle* particle = new TParticle(); | |
63 | ||
64 | TClonesArray *trackRecoArray = NULL; | |
65 | TClonesArray *trackRefArray = NULL; | |
66 | ||
67 | ||
68 | // File for histograms and histogram booking | |
69 | TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE"); | |
70 | ||
71 | TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5); | |
72 | TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5); | |
73 | TH1F *hNHitComp = new TH1F("hNHitComp"," Nb of compatible hits / track ",15,-0.5,14.5); | |
74 | TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5); | |
75 | TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5); | |
76 | ||
77 | TH1F *hResMomVertex = new TH1F("hMomVertex"," delta P vertex (GeV/c)",100,-10.,10); | |
78 | TH1F *hResMomFirstHit = new TH1F("hMomFirstHit"," delta P first hit (GeV/c)",100,-10.,10); | |
79 | ||
a0bfad0e | 80 | // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex) |
81 | if (!gGeoManager) { | |
82 | TGeoManager::Import(geoFilename); | |
83 | if (!gGeoManager) { | |
84 | Error("MUONmass_ESD", "getting geometry from file %s failed", geoFilename); | |
85 | return; | |
86 | } | |
87 | } | |
88 | ||
45a05f07 | 89 | // set mag field |
90 | // waiting for mag field in CDB | |
91 | printf("Loading field map...\n"); | |
b97b210c | 92 | AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG); |
45a05f07 | 93 | AliTracker::SetFieldMap(field, kFALSE); |
a0bfad0e | 94 | // set the magnetic field for track extrapolations |
95 | AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap()); | |
b8dc484b | 96 | |
45a05f07 | 97 | AliRunLoader* runLoader = AliRunLoader::Open(filename,"read"); |
98 | AliLoader * MUONLoader = runLoader->GetLoader("MUONLoader"); | |
99 | AliMUONData * MUONData = new AliMUONData(MUONLoader,"MUON","MUON"); | |
b8dc484b | 100 | |
45a05f07 | 101 | runLoader->LoadKinematics("READ"); |
102 | runLoader->LoadTrackRefs("READ"); | |
103 | MUONLoader->LoadTracks("READ"); | |
104 | ||
105 | AliMUONRecoCheck rc(runLoader,MUONData); | |
b8dc484b | 106 | |
107 | Int_t nevents = runLoader->GetNumberOfEvents(); | |
108 | ||
109 | if (nevents < nEvent) nEvent = nevents; | |
110 | ||
111 | Int_t ievent; | |
112 | Int_t nReconstructibleTracks = 0; | |
113 | Int_t nReconstructedTracks = 0; | |
114 | Int_t nReconstructibleTracksCheck = 0; | |
115 | ||
116 | for (ievent=0; ievent<nEvent; ievent++) { | |
117 | if (!(ievent%10)) printf(" **** event # %d \n",ievent); | |
118 | runLoader->GetEvent(ievent); | |
119 | rc.ResetTracks(); | |
120 | rc.MakeTrackRef(); // make reconstructible tracks | |
121 | // rc.PrintEvent(); | |
122 | ||
123 | ||
124 | trackRecoArray = rc.GetTrackReco(); | |
125 | trackRefArray = rc.GetMuonTrackRef(); | |
126 | ||
127 | nTrackRef = trackRefArray->GetEntriesFast(); | |
128 | nTrackReco = trackRecoArray->GetEntriesFast(); | |
129 | ||
130 | hReconstructible->Fill(rc.GetNumberOfReconstuctibleTracks()); | |
131 | hReco->Fill(rc.GetNumberOfRecoTracks()); | |
132 | ||
133 | nReconstructibleTracks += rc.GetNumberOfReconstuctibleTracks(); | |
134 | nReconstructedTracks += rc.GetNumberOfRecoTracks(); | |
135 | ||
136 | for (Int_t index1 = 0; index1 < nTrackRef; index1++) { | |
137 | trackRef = (AliMUONTrack *)trackRefArray->At(index1); | |
138 | ||
139 | testTrack = 0; | |
140 | indexOK = 0; | |
141 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) | |
142 | compTrackOK[ch] = kFALSE; | |
143 | for (Int_t index2 = 0; index2 < nTrackReco; index2++) { | |
144 | trackReco = (AliMUONTrack *)trackRecoArray->At(index2); | |
145 | ||
146 | // check if trackRef is compatible with trackReco | |
147 | compTrack = trackRef->CompatibleTrack(trackReco,sigma2Cut); | |
148 | ||
149 | iTrack = TrackCheck(compTrack); | |
150 | ||
151 | if (iTrack > testTrack) { | |
152 | nHitOK = 0; | |
153 | for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) { | |
154 | if (compTrack[ch]) nHitOK++; | |
155 | compTrackOK[ch] = compTrack[ch]; | |
156 | } | |
157 | testTrack = iTrack; | |
158 | indexOK = index2; | |
159 | } | |
160 | } | |
161 | ||
162 | hTestTrack->Fill(testTrack); | |
163 | trackID = trackRef->GetTrackID(); | |
164 | hTrackRefID->Fill(trackID); | |
165 | ||
166 | if (testTrack == 4) { // tracking requirements verified, track is found | |
167 | nReconstructibleTracksCheck++; | |
168 | hNHitComp->Fill(nHitOK); | |
169 | particle = runLoader->GetHeader()->Stack()->Particle(trackID); | |
170 | // printf(" trackID: %d , PDG code: %d \n",trackID,particle->GetPdgCode()); | |
171 | trackParam = trackRef->GetTrackParamAtVertex(); | |
172 | x1 = trackParam->GetNonBendingCoor(); | |
173 | y1 = trackParam->GetBendingCoor(); | |
174 | z1 = trackParam->GetZ(); | |
175 | pX1 = trackParam->Px(); | |
176 | pY1 = trackParam->Py(); | |
177 | pZ1 = trackParam->Pz(); | |
178 | p1 = trackParam->P(); | |
179 | ||
180 | // 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); | |
8cde4af5 | 181 | trackReco = (AliMUONTrack *)trackRecoArray->At(indexOK); |
182 | trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtHit()->First()))); | |
183 | AliMUONTrackExtrap::ExtrapToVertex(trackParam,x1,y1,z1); | |
b8dc484b | 184 | x2 = trackParam->GetNonBendingCoor(); |
185 | y2 = trackParam->GetBendingCoor(); | |
186 | z2 = trackParam->GetZ(); | |
187 | pX2 = trackParam->Px(); | |
188 | pY2 = trackParam->Py(); | |
189 | pZ2 = trackParam->Pz(); | |
190 | p2 = trackParam->P(); | |
8cde4af5 | 191 | delete trackParam; |
b8dc484b | 192 | // 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); |
193 | ||
194 | hResMomVertex->Fill(p2-p1); | |
195 | ||
196 | trackParamAtHit = trackRef->GetTrackParamAtHit(); | |
197 | trackParam = (AliMUONTrackParam*) trackParamAtHit->First(); | |
198 | x1 = trackParam->GetNonBendingCoor(); | |
199 | y1 = trackParam->GetBendingCoor(); | |
200 | z1 = trackParam->GetZ(); | |
201 | pX1 = trackParam->Px(); | |
202 | pY1 = trackParam->Py(); | |
203 | pZ1 = trackParam->Pz(); | |
204 | p1 = trackParam->P(); | |
205 | // 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); | |
206 | trackParamAtHit = ((AliMUONTrack *) trackRecoArray->At(indexOK))->GetTrackParamAtHit(); | |
207 | trackParam = (AliMUONTrackParam*) trackParamAtHit->First(); | |
208 | x2 = trackParam->GetNonBendingCoor(); | |
209 | y2 = trackParam->GetBendingCoor(); | |
210 | z2 = trackParam->GetZ(); | |
211 | pX2 = trackParam->Px(); | |
212 | pY2 = trackParam->Py(); | |
213 | pZ2 = trackParam->Pz(); | |
214 | p2 = trackParam->P(); | |
215 | // 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); | |
216 | ||
217 | hResMomFirstHit->Fill(p2-p1); | |
218 | ||
219 | } | |
220 | } // end loop track ref. | |
221 | ||
222 | } // end loop on event | |
223 | ||
45a05f07 | 224 | MUONLoader->UnloadTracks(); |
225 | runLoader->UnloadKinematics(); | |
226 | runLoader->UnloadTrackRefs(); | |
a0bfad0e | 227 | delete runLoader; |
228 | delete field; | |
229 | delete MUONData; | |
45a05f07 | 230 | |
b8dc484b | 231 | printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks); |
232 | printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks); | |
233 | printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck); | |
234 | ||
235 | histoFile->Write(); | |
236 | histoFile->Close(); | |
237 | } | |
238 | ||
239 | ||
240 | Int_t TrackCheck( Bool_t *compTrack) | |
241 | { | |
242 | // Apply reconstruction requirements | |
243 | // Return number of validated conditions | |
244 | // If all the tests are verified then TrackCheck = 4 (good track) | |
245 | Int_t iTrack = 0; | |
246 | Int_t hitsInLastStations = 0; | |
247 | ||
248 | // apply reconstruction requirements | |
249 | if (compTrack[0] || compTrack[1]) iTrack++; // at least one hit in st. 0 | |
250 | if (compTrack[2] || compTrack[3]) iTrack++; // at least one hit in st. 1 | |
251 | if (compTrack[4] || compTrack[5]) iTrack++; // at least one hit in st. 2 | |
252 | for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) { | |
253 | if (compTrack[ch]) hitsInLastStations++; | |
254 | } | |
255 | if (hitsInLastStations > 2) iTrack++; // at least 3 hits in st. 3 & 4 | |
256 | ||
257 | return iTrack; | |
258 | ||
259 | } | |
260 | ||
261 | ||
262 | ||
263 | ||
264 | ||
265 | ||
266 |