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