]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Main changes:
[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
37 #include "AliMUONVTrackStore.h"
38
39 Int_t TrackCheck( Bool_t *compTrack);
40
41 void MUONRecoCheck (Int_t nEvent = 1, char* geoFilename = "geometry.root", 
42                     char * pathSim="./generated/", char * filename="galice.root"){
43   
44   // Utility macro to check the muon reconstruction. Reconstructed tracks are compared
45   // to reference tracks. The reference tracks are built from AliTrackReference for the
46   // hit in chamber (0..9) and from kinematics (TreeK) for the vertex parameters.     
47   
48   Bool_t *compTrack;
49   Bool_t compTrackOK[10];
50   Int_t nClusterOk = 0;
51   Int_t testTrack = 0;  
52   Int_t iTrack = 0;
53   AliMUONTrack* trackOK(0x0);
54   Int_t trackID = 0;
55   Double_t sigma2Cut = 16;  // 4 sigmas cut, sigma2Cut = 4*4
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(filename,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       testTrack = 0;
119       trackOK = 0x0;
120       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
121       {
122         compTrackOK[ch] = kFALSE;
123       }
124       
125       TIter next2(trackStore->CreateIterator());
126       AliMUONTrack* trackReco;
127       
128       while ( ( trackReco = static_cast<AliMUONTrack*>(next2()) ) )
129       {
130         // check if trackRef is compatible with trackReco
131         compTrack = trackRef->CompatibleTrack(trackReco,sigma2Cut);
132         
133         iTrack = TrackCheck(compTrack);
134         
135         if (iTrack > testTrack) 
136         {
137           nClusterOk = 0;
138           for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
139           {
140             if (compTrack[ch]) nClusterOk++;
141             compTrackOK[ch] = compTrack[ch];
142           }
143           testTrack = iTrack;
144           trackOK = trackReco;
145         }
146       }
147       
148       hTestTrack->Fill(testTrack);
149       trackID = trackRef->GetTrackID();
150       hTrackRefID->Fill(trackID);
151       
152       if (testTrack == 4) {     // tracking requirements verified, track is found
153         nReconstructibleTracksCheck++;
154         hNClusterComp->Fill(nClusterOk);
155         trackParam = trackRef->GetTrackParamAtVertex();
156         x1 = trackParam->GetNonBendingCoor();
157         y1 = trackParam->GetBendingCoor();
158         z1 = trackParam->GetZ();
159         pX1 = trackParam->Px();
160         pY1 = trackParam->Py();
161         pZ1 = trackParam->Pz();
162         p1  = trackParam->P();
163         
164         //      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);
165         trackReco = trackOK;
166         trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtCluster()->First())));
167         AliMUONTrackExtrap::ExtrapToVertex(trackParam,x1,y1,z1);
168         x2 = trackParam->GetNonBendingCoor();
169         y2 = trackParam->GetBendingCoor();
170         z2 = trackParam->GetZ();
171         pX2 = trackParam->Px();
172         pY2 = trackParam->Py();
173         pZ2 = trackParam->Pz();
174         p2  = trackParam->P();
175         delete trackParam;
176         //      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);
177         
178         hResMomVertex->Fill(p2-p1);
179         
180         trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
181         x1 = trackParam->GetNonBendingCoor();
182         y1 = trackParam->GetBendingCoor();
183         z1 = trackParam->GetZ();
184         pX1 = trackParam->Px();
185         pY1 = trackParam->Py();
186         pZ1 = trackParam->Pz();
187         p1  = trackParam->P();
188         //      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);
189         trackParam = (AliMUONTrackParam*) trackOK->GetTrackParamAtCluster()->First();
190         x2 = trackParam->GetNonBendingCoor();
191         y2 = trackParam->GetBendingCoor();
192         z2 = trackParam->GetZ();
193         pX2 = trackParam->Px();
194         pY2 = trackParam->Py();
195         pZ2 = trackParam->Pz();
196         p2  = trackParam->P();
197         //      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);
198         
199         hResMomFirstCluster->Fill(p2-p1);
200                
201       }
202     } // end loop track ref.
203
204   } // end loop on event  
205   
206   delete field;
207   
208   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
209   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
210   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
211   
212   histoFile->Write();
213   histoFile->Close();
214 }
215
216
217 Int_t TrackCheck( Bool_t *compTrack)
218 {
219   // Apply reconstruction requirements
220   // Return number of validated conditions 
221   // If all the tests are verified then TrackCheck = 4 (good track)
222   Int_t iTrack = 0;
223   Int_t nCompClustersInLastStations = 0;
224   
225   // apply reconstruction requirements
226   if (compTrack[0] || compTrack[1]) iTrack++; // at least one compatible cluster in st. 0
227   if (compTrack[2] || compTrack[3]) iTrack++; // at least one compatible cluster in st. 1
228   if (compTrack[4] || compTrack[5]) iTrack++; // at least one compatible cluster in st. 2
229   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
230     if (compTrack[ch]) nCompClustersInLastStations++; 
231   }
232   if (nCompClustersInLastStations > 2) iTrack++; // at least 3 compatible clusters in st. 3 & 4
233   
234   return iTrack;
235   
236 }
237