]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
Adding CreateIterator(void) and GetNeighbours() pure virtual methods,
[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 "TParticle.h"
20 #include "TFile.h"
21
22 // STEER includes
23 #include "AliRun.h"
24 #include "AliHeader.h"
25 #include "AliMC.h"
26 #include "AliStack.h"
27 #include "AliRunLoader.h"
28 #include "AliMagFMaps.h"
29
30 // MUON includes
31 #include "AliMUON.h"
32 #include "AliMUONConstants.h"
33 #include "AliMUONTrack.h"
34 #include "AliMUONRecoCheck.h"
35 #include "AliMUONTrackParam.h"
36 #include "AliTracker.h"
37
38 Int_t TrackCheck( Bool_t *compTrack);
39
40 void MUONRecoCheck (Int_t nEvent = 1, char * filename="galice.root"){
41
42   // Utility macro to check the muon reconstruction. Reconstructed tracks are compared
43   // to reference tracks. The reference tracks are built from AliTrackReference for the
44   // hit in chamber (0..9) and from kinematics (TreeK) for the vertex parameters.     
45   
46   Int_t nTrackReco, nTrackRef;
47   AliMUONTrack *trackReco, *trackRef;
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   Int_t indexOK = 0;
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   TParticle* particle = new TParticle();
61
62   TClonesArray *trackRecoArray = NULL;
63   TClonesArray *trackRefArray = NULL;
64
65
66   // File for histograms and histogram booking
67   TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
68
69   TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks ",15,-0.5,14.5);
70   TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
71   TH1F *hNHitComp = new TH1F("hNHitComp"," Nb of compatible hits / track ",15,-0.5,14.5);
72   TH1F *hTestTrack = new TH1F("hTestTrack"," Reconstruction requirement / track",15,-0.5,14.5);
73   TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
74   
75   TH1F *hResMomVertex = new TH1F("hMomVertex"," delta P vertex (GeV/c)",100,-10.,10);
76   TH1F *hResMomFirstHit = new TH1F("hMomFirstHit"," delta P first hit (GeV/c)",100,-10.,10);
77
78   // set  mag field 
79   // waiting for mag field in CDB 
80   printf("Loading field map...\n");
81   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 1, 1., 10., AliMagFMaps::k5kG);
82   AliTracker::SetFieldMap(field, kFALSE);
83
84   AliRunLoader* runLoader = AliRunLoader::Open(filename,"read");
85   AliLoader * MUONLoader = runLoader->GetLoader("MUONLoader");
86   AliMUONData * MUONData = new AliMUONData(MUONLoader,"MUON","MUON");
87
88   runLoader->LoadKinematics("READ");
89   runLoader->LoadTrackRefs("READ");
90   MUONLoader->LoadTracks("READ");
91   
92   AliMUONRecoCheck rc(runLoader,MUONData);
93     
94   Int_t nevents = runLoader->GetNumberOfEvents();
95   
96   if (nevents < nEvent) nEvent = nevents;
97   
98   Int_t ievent;
99   Int_t nReconstructibleTracks = 0;
100   Int_t nReconstructedTracks = 0;
101   Int_t nReconstructibleTracksCheck = 0;
102
103   for (ievent=0; ievent<nEvent; ievent++) {
104     if (!(ievent%10)) printf(" **** event # %d  \n",ievent);
105     runLoader->GetEvent(ievent);
106     rc.ResetTracks();
107     rc.MakeTrackRef(); // make reconstructible tracks
108 //     rc.PrintEvent();
109
110     
111     trackRecoArray = rc.GetTrackReco();
112     trackRefArray = rc.GetMuonTrackRef();
113     
114     nTrackRef = trackRefArray->GetEntriesFast();
115     nTrackReco = trackRecoArray->GetEntriesFast();
116
117     hReconstructible->Fill(rc.GetNumberOfReconstuctibleTracks());
118     hReco->Fill(rc.GetNumberOfRecoTracks());
119
120     nReconstructibleTracks += rc.GetNumberOfReconstuctibleTracks();
121     nReconstructedTracks += rc.GetNumberOfRecoTracks();
122
123     for (Int_t index1 = 0; index1 < nTrackRef; index1++) {  
124       trackRef = (AliMUONTrack *)trackRefArray->At(index1);
125
126       testTrack = 0;
127       indexOK = 0;
128       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
129         compTrackOK[ch] = kFALSE;
130       for (Int_t index2 = 0; index2 < nTrackReco; index2++) {
131         trackReco = (AliMUONTrack *)trackRecoArray->At(index2);
132
133         // check if trackRef is compatible with trackReco
134         compTrack = trackRef->CompatibleTrack(trackReco,sigma2Cut);
135
136         iTrack = TrackCheck(compTrack);
137         
138         if (iTrack > testTrack) {
139           nHitOK = 0;
140           for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) {
141             if (compTrack[ch]) nHitOK++;
142             compTrackOK[ch] = compTrack[ch];
143           }
144           testTrack = iTrack;
145           indexOK = index2;
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         particle = runLoader->GetHeader()->Stack()->Particle(trackID);
157 //      printf(" trackID: %d , PDG code: %d \n",trackID,particle->GetPdgCode());
158         trackParam = trackRef->GetTrackParamAtVertex();
159         x1 = trackParam->GetNonBendingCoor();
160         y1 = trackParam->GetBendingCoor();
161         z1 = trackParam->GetZ();
162         pX1 = trackParam->Px();
163         pY1 = trackParam->Py();
164         pZ1 = trackParam->Pz();
165         p1  = trackParam->P();
166         
167 //      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);
168         
169         trackParam = ((AliMUONTrack *)trackRecoArray->At(indexOK))->GetTrackParamAtVertex();
170         x2 = trackParam->GetNonBendingCoor();
171         y2 = trackParam->GetBendingCoor();
172         z2 = trackParam->GetZ();
173         pX2 = trackParam->Px();
174         pY2 = trackParam->Py();
175         pZ2 = trackParam->Pz();
176         p2  = trackParam->P();
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 =  ((AliMUONTrack *) trackRecoArray->At(indexOK))->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   MUONLoader->UnloadTracks();
210   runLoader->UnloadKinematics();
211   runLoader->UnloadTrackRefs();
212   runLoader->Delete();
213   field->Delete();
214
215   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
216   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
217   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
218   
219   histoFile->Write();
220   histoFile->Close();
221 }
222
223
224 Int_t TrackCheck( Bool_t *compTrack)
225 {
226   // Apply reconstruction requirements
227   // Return number of validated conditions 
228   // If all the tests are verified then TrackCheck = 4 (good track)
229   Int_t iTrack = 0;
230   Int_t hitsInLastStations = 0;
231   
232   // apply reconstruction requirements
233   if (compTrack[0] || compTrack[1]) iTrack++; // at least one hit in st. 0
234   if (compTrack[2] || compTrack[3]) iTrack++; // at least one hit in st. 1
235   if (compTrack[4] || compTrack[5]) iTrack++; // at least one hit in st. 2
236   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
237     if (compTrack[ch]) hitsInLastStations++; 
238   }
239   if (hitsInLastStations > 2) iTrack++; // at least 3 hits in st. 3 & 4
240   
241   return iTrack;
242
243 }
244
245
246
247
248
249
250