]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Major rewriting using VStores (Laurent)
[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 * filenameSim="galice_sim.root", 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 nHitOK = 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   TClonesArray *trackParamAtHit;
58   Double_t x1,y1,z1,pX1,pY1,pZ1,p1;
59   Double_t x2,y2,z2,pX2,pY2,pZ2,p2;
60   
61   // File for histograms and histogram booking
62   TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
63   
64   TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5);
65   TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
66   TH1F *hNHitComp = new TH1F("hNHitComp"," Nb of compatible hits / track ",15,-0.5,14.5);
67   TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5);
68   TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
69   
70   TH1F *hResMomVertex = new TH1F("hMomVertex"," delta P vertex (GeV/c)",100,-10.,10);
71   TH1F *hResMomFirstHit = new TH1F("hMomFirstHit"," delta P first hit (GeV/c)",100,-10.,10);
72   
73   // Import TGeo geometry (needed by AliMUONTrackExtrap::ExtrapToVertex)
74   if (!gGeoManager) {
75     TGeoManager::Import(geoFilename);
76     if (!gGeoManager) {
77       Error("MUONmass_ESD", "getting geometry from file %s failed", geoFilename);
78       return;
79     }
80   }
81   
82   // set  mag field 
83   // waiting for mag field in CDB 
84   printf("Loading field map...\n");
85   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
86   AliTracker::SetFieldMap(field, kFALSE);
87   // set the magnetic field for track extrapolations
88   AliMUONTrackExtrap::SetField(AliTracker::GetFieldMap());
89     
90   AliMUONRecoCheck rc(filename,filenameSim);
91   
92   Int_t nevents = rc.NumberOfEvents();
93   
94   if (nevents < nEvent) nEvent = nevents;
95   
96   Int_t ievent;
97   Int_t nReconstructibleTracks = 0;
98   Int_t nReconstructedTracks = 0;
99   Int_t nReconstructibleTracksCheck = 0;
100   
101   for (ievent=0; ievent<nEvent; ievent++)
102   {
103     if (!(ievent%10)) printf(" **** event # %d  \n",ievent);
104     
105     AliMUONVTrackStore* trackStore = rc.ReconstructedTracks(ievent);
106     AliMUONVTrackStore* trackRefStore = rc.ReconstructibleTracks(ievent);
107     
108     hReconstructible->Fill(trackRefStore->GetSize());
109     hReco->Fill(trackStore->GetSize());
110     
111     nReconstructibleTracks += trackRefStore->GetSize();
112     nReconstructedTracks += trackStore->GetSize();
113     
114     TIter next(trackRefStore->CreateIterator());
115     AliMUONTrack* trackRef;
116     
117     while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) )
118     {
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         compTrack = trackRef->CompatibleTrack(trackReco,sigma2Cut);
133         
134         iTrack = TrackCheck(compTrack);
135         
136         if (iTrack > testTrack) 
137         {
138           nHitOK = 0;
139           for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
140           {
141             if (compTrack[ch]) nHitOK++;
142             compTrackOK[ch] = compTrack[ch];
143           }
144           testTrack = iTrack;
145           trackOK = trackReco;
146         }
147       }
148       
149       hTestTrack->Fill(testTrack);
150       trackID = trackRef->GetTrackID();
151       hTrackRefID->Fill(trackID);
152       
153       if (testTrack == 4) {     // tracking requirements verified, track is found
154         nReconstructibleTracksCheck++;
155         hNHitComp->Fill(nHitOK);
156         trackParam = trackRef->GetTrackParamAtVertex();
157         x1 = trackParam->GetNonBendingCoor();
158         y1 = trackParam->GetBendingCoor();
159         z1 = trackParam->GetZ();
160         pX1 = trackParam->Px();
161         pY1 = trackParam->Py();
162         pZ1 = trackParam->Pz();
163         p1  = trackParam->P();
164         
165         //      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);
166         trackReco = trackOK;
167         trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtHit()->First())));
168         AliMUONTrackExtrap::ExtrapToVertex(trackParam,x1,y1,z1);
169         x2 = trackParam->GetNonBendingCoor();
170         y2 = trackParam->GetBendingCoor();
171         z2 = trackParam->GetZ();
172         pX2 = trackParam->Px();
173         pY2 = trackParam->Py();
174         pZ2 = trackParam->Pz();
175         p2  = trackParam->P();
176         delete trackParam;
177         //      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);
178         
179         hResMomVertex->Fill(p2-p1);
180         
181         trackParamAtHit =  trackRef->GetTrackParamAtHit();
182         trackParam = (AliMUONTrackParam*) trackParamAtHit->First();
183         x1 = trackParam->GetNonBendingCoor();
184         y1 = trackParam->GetBendingCoor();
185         z1 = trackParam->GetZ();
186         pX1 = trackParam->Px();
187         pY1 = trackParam->Py();
188         pZ1 = trackParam->Pz();
189         p1  = trackParam->P();
190         //      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);
191         trackParamAtHit =  trackOK->GetTrackParamAtHit();
192         trackParam = (AliMUONTrackParam*) trackParamAtHit->First();
193         x2 = trackParam->GetNonBendingCoor();
194         y2 = trackParam->GetBendingCoor();
195         z2 = trackParam->GetZ();
196         pX2 = trackParam->Px();
197         pY2 = trackParam->Py();
198         pZ2 = trackParam->Pz();
199         p2  = trackParam->P();
200         //      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);
201         
202         hResMomFirstHit->Fill(p2-p1);
203                
204       }
205     } // end loop track ref.
206
207   } // end loop on event  
208   
209   delete field;
210   
211   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
212   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
213   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
214   
215   histoFile->Write();
216   histoFile->Close();
217 }
218
219
220 Int_t TrackCheck( Bool_t *compTrack)
221 {
222   // Apply reconstruction requirements
223   // Return number of validated conditions 
224   // If all the tests are verified then TrackCheck = 4 (good track)
225   Int_t iTrack = 0;
226   Int_t hitsInLastStations = 0;
227   
228   // apply reconstruction requirements
229   if (compTrack[0] || compTrack[1]) iTrack++; // at least one hit in st. 0
230   if (compTrack[2] || compTrack[3]) iTrack++; // at least one hit in st. 1
231   if (compTrack[4] || compTrack[5]) iTrack++; // at least one hit in st. 2
232   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
233     if (compTrack[ch]) hitsInLastStations++; 
234   }
235   if (hitsInLastStations > 2) iTrack++; // at least 3 hits in st. 3 & 4
236   
237   return iTrack;
238   
239 }
240
241
242
243
244
245
246