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