]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Oups. We probably do not want the memory leak checks part of the production version
[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
33 // STEER includes
34 #include "AliCDBManager.h"
35
36 // MUON includes
37 #include "AliMUONCDB.h"
38 #include "AliMUONConstants.h"
39 #include "AliMUONTrack.h"
40 #include "AliMUONRecoCheck.h"
41 #include "AliMUONTrackParam.h"
42 #include "AliMUONTrackExtrap.h"
43 #include "AliMUONRecoParam.h"
44 #include "AliMUONVTrackStore.h"
45
46 Int_t TrackCheck( Bool_t *compTrack);
47
48 void MUONRecoCheck (Int_t nEvent = -1, char * pathSim="./generated/", char * esdFileName="AliESDs.root",
49                     const char* ocdbPath = "local://$ALICE_ROOT/OCDB")
50 {
51   
52   Bool_t compTrack[10];
53   Bool_t compTrackOK[10];
54   Int_t nClusterOk = 0;
55   Int_t testTrack = 0;  
56   Int_t iTrack = 0;
57   AliMUONTrack* trackOK(0x0);
58   Int_t trackID = 0;
59   Double_t maxChi2 = 999.;
60   AliMUONTrackParam *trackParam;
61   Double_t x1,y1,z1,pX1,pY1,pZ1,p1;
62   Double_t x2,y2,z2,pX2,pY2,pZ2,p2;
63   
64   // File for histograms and histogram booking
65   TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
66   
67   TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5);
68   TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
69   TH1F *hNClusterComp = new TH1F("hNClusterComp"," Nb of compatible clusters / track ",15,-0.5,14.5);
70   TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5);
71   TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
72   
73   TH1F *hResMomVertex = new TH1F("hResMomVertex"," delta P vertex (GeV/c)",100,-10.,10);
74   TH1F *hResMomFirstCluster = new TH1F("hResMomFirstCluster"," delta P first cluster (GeV/c)",100,-10.,10);
75   
76   AliMUONRecoCheck rc(esdFileName, pathSim);
77   
78   // load necessary data from OCDB
79   AliCDBManager::Instance()->SetDefaultStorage(ocdbPath);
80   AliCDBManager::Instance()->SetRun(rc.GetRunNumber());
81   AliMUONCDB::LoadField();
82   AliMUONRecoParam* recoParam = AliMUONCDB::LoadRecoParam();
83   if (!recoParam) return;
84   
85   // set the magnetic field for track extrapolations
86   AliMUONTrackExtrap::SetField();
87
88   // get sigma cut from recoParam to associate clusters with TrackRefs
89   Double_t sigmaCut = (recoParam->ImproveTracks()) ? recoParam->GetSigmaCutForImprovement() : recoParam->GetSigmaCutForTracking();
90   
91   Int_t nevents = rc.NumberOfEvents();
92   
93   if (nevents < nEvent || nEvent < 0) 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, kFALSE);
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           trackRef->CompatibleTrack(trackReco,sigmaCut,compTrack);
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           Double_t chi2;
157           if (AliMUONTrackExtrap::ExtrapToZCov(&recoParam, refParam->GetZ()) &&
158               refParam->CompatibleTrackParam(recoParam, sigmaCut, chi2) && chi2 < maxChi2)
159           {
160               maxChi2 = chi2;
161               trackOK = trackReco;
162           }
163           
164         }
165         
166       }
167       
168       hTestTrack->Fill(testTrack);
169       trackID = trackRef->GetMCLabel();
170       hTrackRefID->Fill(trackID);
171       
172       if (testTrack == 4 || maxChi2 < 5.*sigmaCut*sigmaCut) {     // tracking requirements verified, track is found
173         nReconstructibleTracksCheck++;
174         hNClusterComp->Fill(nClusterOk);
175         trackParam = trackRef->GetTrackParamAtVertex();
176         x1 = trackParam->GetNonBendingCoor();
177         y1 = trackParam->GetBendingCoor();
178         z1 = trackParam->GetZ();
179         pX1 = trackParam->Px();
180         pY1 = trackParam->Py();
181         pZ1 = trackParam->Pz();
182         p1  = trackParam->P();
183         
184         //      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);
185         trackReco = trackOK;
186         trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtVertex())));
187         x2 = trackParam->GetNonBendingCoor();
188         y2 = trackParam->GetBendingCoor();
189         z2 = trackParam->GetZ();
190         pX2 = trackParam->Px();
191         pY2 = trackParam->Py();
192         pZ2 = trackParam->Pz();
193         p2  = trackParam->P();
194         delete trackParam;
195         //      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);
196         
197         hResMomVertex->Fill(p2-p1);
198         
199         trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
200         x1 = trackParam->GetNonBendingCoor();
201         y1 = trackParam->GetBendingCoor();
202         z1 = trackParam->GetZ();
203         pX1 = trackParam->Px();
204         pY1 = trackParam->Py();
205         pZ1 = trackParam->Pz();
206         p1  = trackParam->P();
207         //      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);
208         trackParam = (AliMUONTrackParam*) trackOK->GetTrackParamAtCluster()->First();
209         x2 = trackParam->GetNonBendingCoor();
210         y2 = trackParam->GetBendingCoor();
211         z2 = trackParam->GetZ();
212         pX2 = trackParam->Px();
213         pY2 = trackParam->Py();
214         pZ2 = trackParam->Pz();
215         p2  = trackParam->P();
216         //      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);
217         
218         hResMomFirstCluster->Fill(p2-p1);
219                
220       }
221     } // end loop track ref.
222
223   } // end loop on event  
224   
225   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
226   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
227   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
228   
229   histoFile->Write();
230   histoFile->Close();
231 }
232
233
234 Int_t TrackCheck( Bool_t *compTrack)
235 {
236   // Apply reconstruction requirements
237   // Return number of validated conditions 
238   // If all the tests are verified then TrackCheck = 4 (good track)
239   Int_t iTrack = 0;
240   Int_t nCompClustersInLastStations = 0;
241   
242   // apply reconstruction requirements
243   if (compTrack[0] || compTrack[1]) iTrack++; // at least one compatible cluster in st. 0
244   if (compTrack[2] || compTrack[3]) iTrack++; // at least one compatible cluster in st. 1
245   if (compTrack[4] || compTrack[5]) iTrack++; // at least one compatible cluster in st. 2
246   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
247     if (compTrack[ch]) nCompClustersInLastStations++; 
248   }
249   if (nCompClustersInLastStations > 2) iTrack++; // at least 3 compatible clusters in st. 3 & 4
250   
251   return iTrack;
252   
253 }
254