]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Fill ESD MUON clusters with additional informations (including pads)
[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 "TFile.h"
20 #include <TGeoManager.h>
21
22 // STEER includes
23 #include "AliRun.h"
24 #include "AliHeader.h"
25 #include "AliMC.h"
26 #include "AliStack.h"
27 #include "AliMagFMaps.h"
28 #include "AliTracker.h"
29
30 // MUON includes
31 #include "AliMUONConstants.h"
32 #include "AliMUONTrack.h"
33 #include "AliMUONRecoCheck.h"
34 #include "AliMUONTrackParam.h"
35 #include "AliMUONTrackExtrap.h"
36 #include "AliMUONVTrackStore.h"
37
38 Int_t TrackCheck( Bool_t *compTrack);
39
40 void MUONRecoCheck (Int_t nEvent = 1, char* geoFilename = "geometry.root", 
41                     char * pathSim="./generated/", char * esdFileName="AliESDs.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   Bool_t *compTrack;
48   Bool_t compTrackOK[10];
49   Int_t nClusterOk = 0;
50   Int_t testTrack = 0;  
51   Int_t iTrack = 0;
52   AliMUONTrack* trackOK(0x0);
53   Int_t trackID = 0;
54   Double_t sigmaCut = 4.;  // 4 sigmas cut
55   Double_t maxChi2 = 999.;
56   AliMUONTrackParam *trackParam;
57   Double_t x1,y1,z1,pX1,pY1,pZ1,p1;
58   Double_t x2,y2,z2,pX2,pY2,pZ2,p2;
59   
60   // File for histograms and histogram booking
61   TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
62   
63   TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5);
64   TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
65   TH1F *hNClusterComp = new TH1F("hNClusterComp"," Nb of compatible clusters / track ",15,-0.5,14.5);
66   TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5);
67   TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
68   
69   TH1F *hResMomVertex = new TH1F("hResMomVertex"," delta P vertex (GeV/c)",100,-10.,10);
70   TH1F *hResMomFirstCluster = new TH1F("hResMomFirstCluster"," delta P first cluster (GeV/c)",100,-10.,10);
71   
72   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
73   if (!gGeoManager) {
74     TGeoManager::Import(geoFilename);
75     if (!gGeoManager) {
76       Error("MUONmass_ESD", "getting geometry from file %s failed", geoFilename);
77       return;
78     }
79   }
80   
81   // set  mag field 
82   // waiting for mag field in CDB 
83   printf("Loading field map...\n");
84   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
85   AliTracker::SetFieldMap(field, kFALSE);
86   // set the magnetic field for track extrapolations
87   AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
88     
89   AliMUONRecoCheck rc(esdFileName, pathSim);
90   
91   Int_t nevents = rc.NumberOfEvents();
92   
93   if (nevents < nEvent) nEvent = nevents;
94   
95   Int_t ievent;
96   Int_t nReconstructibleTracks = 0;
97   Int_t nReconstructedTracks = 0;
98   Int_t nReconstructibleTracksCheck = 0;
99   
100   for (ievent=0; ievent<nEvent; ievent++)
101   {
102     if (!(ievent%10)) printf(" **** event # %d  \n",ievent);
103     
104     AliMUONVTrackStore* trackStore = rc.ReconstructedTracks(ievent);
105     AliMUONVTrackStore* trackRefStore = rc.ReconstructibleTracks(ievent);
106     
107     hReconstructible->Fill(trackRefStore->GetSize());
108     hReco->Fill(trackStore->GetSize());
109     
110     nReconstructibleTracks += trackRefStore->GetSize();
111     nReconstructedTracks += trackStore->GetSize();
112     
113     TIter next(trackRefStore->CreateIterator());
114     AliMUONTrack* trackRef;
115     
116     while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) )
117     {
118       maxChi2 = 999.;
119       testTrack = 0;
120       trackOK = 0x0;
121       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
122       {
123         compTrackOK[ch] = kFALSE;
124       }
125       
126       TIter next2(trackStore->CreateIterator());
127       AliMUONTrack* trackReco;
128       
129       while ( ( trackReco = static_cast<AliMUONTrack*>(next2()) ) )
130       {
131         // check if trackRef is compatible with trackReco
132         if (trackReco->GetNClusters() > 1) {
133           
134           // check cluster by cluster if trackReco contain info at each cluster
135           compTrack = trackRef->CompatibleTrack(trackReco,sigmaCut);
136           
137           iTrack = TrackCheck(compTrack);
138           
139           if (iTrack > testTrack) 
140           {
141             nClusterOk = 0;
142             for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
143             {
144               if (compTrack[ch]) nClusterOk++;
145               compTrackOK[ch] = compTrack[ch];
146             }
147             testTrack = iTrack;
148             trackOK = trackReco;
149           }
150           
151         } else {
152           
153           // check only parameters at the z position of the first trackRef
154           AliMUONTrackParam *refParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
155           AliMUONTrackParam recoParam(*((AliMUONTrackParam*) trackReco->GetTrackParamAtCluster()->First()));
156           AliMUONTrackExtrap::ExtrapToZCov(&recoParam, refParam->GetZ());
157           Double_t chi2;
158           if (refParam->CompatibleTrackParam(recoParam, sigmaCut, chi2)) {
159             
160             if (chi2 < maxChi2) {
161               maxChi2 = chi2;
162               trackOK = trackReco;
163             }
164             
165           }
166           
167         }
168         
169       }
170       
171       hTestTrack->Fill(testTrack);
172       trackID = trackRef->GetTrackID();
173       hTrackRefID->Fill(trackID);
174       
175       if (testTrack == 4 || maxChi2 < 5.*sigmaCut*sigmaCut) {     // tracking requirements verified, track is found
176         nReconstructibleTracksCheck++;
177         hNClusterComp->Fill(nClusterOk);
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 = trackOK;
189         trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtCluster()->First())));
190         AliMUONTrackExtrap::ExtrapToVertex(trackParam,x1,y1,z1,0.,0.);
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         trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
204         x1 = trackParam->GetNonBendingCoor();
205         y1 = trackParam->GetBendingCoor();
206         z1 = trackParam->GetZ();
207         pX1 = trackParam->Px();
208         pY1 = trackParam->Py();
209         pZ1 = trackParam->Pz();
210         p1  = trackParam->P();
211         //      printf(" Ref. track at 1st cluster: x,y,z: %f %f %f px,py,pz: %f %f %f \n",x1,y1,z1,pX1,pY1,pZ1);
212         trackParam = (AliMUONTrackParam*) trackOK->GetTrackParamAtCluster()->First();
213         x2 = trackParam->GetNonBendingCoor();
214         y2 = trackParam->GetBendingCoor();
215         z2 = trackParam->GetZ();
216         pX2 = trackParam->Px();
217         pY2 = trackParam->Py();
218         pZ2 = trackParam->Pz();
219         p2  = trackParam->P();
220         //      printf(" Reconst. track at 1st cluster: x,y,z: %f %f %f px,py,pz: %f %f %f \n",x2,y2,z2,pX2,pY2,pZ2);
221         
222         hResMomFirstCluster->Fill(p2-p1);
223                
224       }
225     } // end loop track ref.
226
227   } // end loop on event  
228   
229   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
230   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
231   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
232   
233   histoFile->Write();
234   histoFile->Close();
235 }
236
237
238 Int_t TrackCheck( Bool_t *compTrack)
239 {
240   // Apply reconstruction requirements
241   // Return number of validated conditions 
242   // If all the tests are verified then TrackCheck = 4 (good track)
243   Int_t iTrack = 0;
244   Int_t nCompClustersInLastStations = 0;
245   
246   // apply reconstruction requirements
247   if (compTrack[0] || compTrack[1]) iTrack++; // at least one compatible cluster in st. 0
248   if (compTrack[2] || compTrack[3]) iTrack++; // at least one compatible cluster in st. 1
249   if (compTrack[4] || compTrack[5]) iTrack++; // at least one compatible cluster in st. 2
250   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
251     if (compTrack[ch]) nCompClustersInLastStations++; 
252   }
253   if (nCompClustersInLastStations > 2) iTrack++; // at least 3 compatible clusters in st. 3 & 4
254   
255   return iTrack;
256   
257 }
258