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