]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/MUONRecoCheck.C
QA ref defaut storage setter in sim and rec
[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 "AliMagF.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[10];
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   if (!TGeoGlobalMagField::Instance()->GetField()) {
93     printf("Loading field map...\n");
94     AliMagF* field = new AliMagF("Maps","Maps",2,1.,1., 10.,AliMagF::k5kG);
95     TGeoGlobalMagField::Instance()->SetField(field);
96   }
97   // set the magnetic field for track extrapolations
98   AliMUONTrackExtrap::SetField();
99
100   AliMUONRecoCheck rc(esdFileName, pathSim);
101   
102   Int_t nevents = rc.NumberOfEvents();
103   
104   if (nevents < nEvent) nEvent = nevents;
105   
106   Int_t ievent;
107   Int_t nReconstructibleTracks = 0;
108   Int_t nReconstructedTracks = 0;
109   Int_t nReconstructibleTracksCheck = 0;
110   
111   for (ievent=0; ievent<nEvent; ievent++)
112   {
113     if (!(ievent%10)) printf(" **** event # %d  \n",ievent);
114     
115     AliMUONVTrackStore* trackStore = rc.ReconstructedTracks(ievent);
116     AliMUONVTrackStore* trackRefStore = rc.ReconstructibleTracks(ievent);
117     
118     hReconstructible->Fill(trackRefStore->GetSize());
119     hReco->Fill(trackStore->GetSize());
120     
121     nReconstructibleTracks += trackRefStore->GetSize();
122     nReconstructedTracks += trackStore->GetSize();
123     
124     TIter next(trackRefStore->CreateIterator());
125     AliMUONTrack* trackRef;
126     
127     while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) )
128     {
129       maxChi2 = 999.;
130       testTrack = 0;
131       trackOK = 0x0;
132       for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
133       {
134         compTrackOK[ch] = kFALSE;
135       }
136       
137       TIter next2(trackStore->CreateIterator());
138       AliMUONTrack* trackReco;
139       
140       while ( ( trackReco = static_cast<AliMUONTrack*>(next2()) ) )
141       {
142         // check if trackRef is compatible with trackReco
143         if (trackReco->GetNClusters() > 1) {
144           
145           // check cluster by cluster if trackReco contain info at each cluster
146           trackRef->CompatibleTrack(trackReco,sigmaCut,compTrack);
147           
148           iTrack = TrackCheck(compTrack);
149           
150           if (iTrack > testTrack) 
151           {
152             nClusterOk = 0;
153             for (Int_t ch = 0; ch < AliMUONConstants::NTrackingCh(); ch++) 
154             {
155               if (compTrack[ch]) nClusterOk++;
156               compTrackOK[ch] = compTrack[ch];
157             }
158             testTrack = iTrack;
159             trackOK = trackReco;
160           }
161           
162         } else {
163           
164           // check only parameters at the z position of the first trackRef
165           AliMUONTrackParam *refParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
166           AliMUONTrackParam recoParam(*((AliMUONTrackParam*) trackReco->GetTrackParamAtCluster()->First()));
167           AliMUONTrackExtrap::ExtrapToZCov(&recoParam, refParam->GetZ());
168           Double_t chi2;
169           if (refParam->CompatibleTrackParam(recoParam, sigmaCut, chi2)) {
170             
171             if (chi2 < maxChi2) {
172               maxChi2 = chi2;
173               trackOK = trackReco;
174             }
175             
176           }
177           
178         }
179         
180       }
181       
182       hTestTrack->Fill(testTrack);
183       trackID = trackRef->GetTrackID();
184       hTrackRefID->Fill(trackID);
185       
186       if (testTrack == 4 || maxChi2 < 5.*sigmaCut*sigmaCut) {     // tracking requirements verified, track is found
187         nReconstructibleTracksCheck++;
188         hNClusterComp->Fill(nClusterOk);
189         trackParam = trackRef->GetTrackParamAtVertex();
190         x1 = trackParam->GetNonBendingCoor();
191         y1 = trackParam->GetBendingCoor();
192         z1 = trackParam->GetZ();
193         pX1 = trackParam->Px();
194         pY1 = trackParam->Py();
195         pZ1 = trackParam->Pz();
196         p1  = trackParam->P();
197         
198         //      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);
199         trackReco = trackOK;
200         trackParam = new AliMUONTrackParam(*((AliMUONTrackParam*)(trackReco->GetTrackParamAtCluster()->First())));
201         AliMUONTrackExtrap::ExtrapToVertex(trackParam,x1,y1,z1,0.,0.);
202         x2 = trackParam->GetNonBendingCoor();
203         y2 = trackParam->GetBendingCoor();
204         z2 = trackParam->GetZ();
205         pX2 = trackParam->Px();
206         pY2 = trackParam->Py();
207         pZ2 = trackParam->Pz();
208         p2  = trackParam->P();
209         delete trackParam;
210         //      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);
211         
212         hResMomVertex->Fill(p2-p1);
213         
214         trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
215         x1 = trackParam->GetNonBendingCoor();
216         y1 = trackParam->GetBendingCoor();
217         z1 = trackParam->GetZ();
218         pX1 = trackParam->Px();
219         pY1 = trackParam->Py();
220         pZ1 = trackParam->Pz();
221         p1  = trackParam->P();
222         //      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);
223         trackParam = (AliMUONTrackParam*) trackOK->GetTrackParamAtCluster()->First();
224         x2 = trackParam->GetNonBendingCoor();
225         y2 = trackParam->GetBendingCoor();
226         z2 = trackParam->GetZ();
227         pX2 = trackParam->Px();
228         pY2 = trackParam->Py();
229         pZ2 = trackParam->Pz();
230         p2  = trackParam->P();
231         //      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);
232         
233         hResMomFirstCluster->Fill(p2-p1);
234                
235       }
236     } // end loop track ref.
237
238   } // end loop on event  
239   
240   printf(" nb of reconstructible tracks: %d \n", nReconstructibleTracks);
241   printf(" nb of reconstructed tracks: %d \n", nReconstructedTracks);
242   printf(" nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
243   
244   histoFile->Write();
245   histoFile->Close();
246 }
247
248
249 Int_t TrackCheck( Bool_t *compTrack)
250 {
251   // Apply reconstruction requirements
252   // Return number of validated conditions 
253   // If all the tests are verified then TrackCheck = 4 (good track)
254   Int_t iTrack = 0;
255   Int_t nCompClustersInLastStations = 0;
256   
257   // apply reconstruction requirements
258   if (compTrack[0] || compTrack[1]) iTrack++; // at least one compatible cluster in st. 0
259   if (compTrack[2] || compTrack[3]) iTrack++; // at least one compatible cluster in st. 1
260   if (compTrack[4] || compTrack[5]) iTrack++; // at least one compatible cluster in st. 2
261   for (Int_t ch = 6; ch < AliMUONConstants::NTrackingCh(); ch++) {
262     if (compTrack[ch]) nCompClustersInLastStations++; 
263   }
264   if (nCompClustersInLastStations > 2) iTrack++; // at least 3 compatible clusters in st. 3 & 4
265   
266   return iTrack;
267   
268 }
269