]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/MUONRecoCheck.C
removing TPC analysis classes from end 2006
[u/mrichter/AliRoot.git] / MUON / MUONRecoCheck.C
CommitLineData
b8dc484b 1/**************************************************************************
48217459 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**************************************************************************/
b8dc484b 15
e54bf126 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///
c687aa97 26/// \author Jean-Pierre Cussonneau, Philippe Pillot, Subatech
e54bf126 27
b8dc484b 28// ROOT includes
e1fe81b8 29#include <Riostream.h>
f202486b 30#include "TMath.h"
b8dc484b 31#include "TClonesArray.h"
32#include "TH1.h"
f202486b 33#include "TH2.h"
c687aa97 34#include "TH3.h"
f202486b 35#include "TGraphErrors.h"
e1fe81b8 36#include "TGraphAsymmErrors.h"
f202486b 37#include "TF1.h"
b8dc484b 38#include "TFile.h"
c687aa97 39#include "TCanvas.h"
40#include "TLegend.h"
e1fe81b8 41#include "TGeoManager.h"
b8dc484b 42
43// STEER includes
a99c3449 44#include "AliCDBManager.h"
e1fe81b8 45#include "AliGeomManager.h"
c687aa97 46#include "AliLog.h"
b8dc484b 47
48// MUON includes
a99c3449 49#include "AliMUONCDB.h"
b8dc484b 50#include "AliMUONConstants.h"
51#include "AliMUONTrack.h"
52#include "AliMUONRecoCheck.h"
53#include "AliMUONTrackParam.h"
a99c3449 54#include "AliMUONRecoParam.h"
48217459 55#include "AliMUONVTrackStore.h"
f202486b 56#include "AliMUONVCluster.h"
c687aa97 57#include "AliMUONTrackExtrap.h"
e1fe81b8 58#include "AliMUONESDInterface.h"
9f164762 59#include "AliMUONVTriggerTrackStore.h"
60#include "AliMUONTriggerTrack.h"
b8dc484b 61
c687aa97 62Double_t langaufun(Double_t *x, Double_t *par);
e1fe81b8 63void FitGausResVsMom(TH2* h, Int_t nBins, const Double_t mean0, const Double_t sigma0, const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma);
64void FitPDCAVsMom(TH2* h, Int_t nBins, const char* fitting, TGraphAsymmErrors* gSigma);
65TCanvas* DrawVsAng(const char* name, const char* title, TH1* h1, TH2* h2);
66TCanvas* DrawVsPos(const char* name, const char* title, TH2* h1, TH2* h2, TH2* h3);
67TCanvas* DrawResMomVsMom(const char* name, const char* title, TH2* h, Int_t nBins, TF1* f2 = 0x0, const char* fitting = "");
c687aa97 68
69//------------------------------------------------------------------------------------
70void MUONRecoCheck (Int_t nEvent = -1, const char* pathSim="./generated/", const char* esdFileName="AliESDs.root",
e1fe81b8 71 const char* ocdbPath = "local://$ALICE_ROOT/OCDB", Int_t absorberRegion = -1)
e54bf126 72{
e1fe81b8 73 /// Associate the reconstructed tracks with the simulated ones and check the quality of the reconstruction
74 /// (tracking/trigger efficiency; momentum, slope,... resolutions at first cluster and at vertex; cluster resolution).
75 /// You can limit the calculation of track resolution at vertex to the tracks crossing the absorber in a given region
76 /// with the flag "absorberRegion": -1=all, 1=[2,3]deg, 2=[3,10]deg.
77
78 Double_t aAbsLimits[2];
79 if (absorberRegion > -1) {
80 if (absorberRegion == 1) {
81 aAbsLimits[0] = 2.;
82 aAbsLimits[1] = 3.;
83 } else if (absorberRegion == 2) {
84 aAbsLimits[0] = 3.;
85 aAbsLimits[1] = 10.;
86 } else {
87 cout<<"Unknown absorber region. Valid choices are: -1=all, 1=[2,3]deg, 2=[3,10]deg"<<endl;
88 return;
89 }
90 } else {
91 aAbsLimits[0] = 0.;
92 aAbsLimits[1] = 90.;
93 }
b8dc484b 94
c687aa97 95 AliLog::SetClassDebugLevel("AliMCEvent",-1);
96
97 // ###################################### define histograms ###################################### //
b8dc484b 98 // File for histograms and histogram booking
99 TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
48217459 100
9f164762 101 TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks / evt",15,-0.5,14.5);
b8dc484b 102 TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
96ebe67e 103 TH1F *hNClusterComp = new TH1F("hNClusterComp"," Nb of compatible clusters / track ",15,-0.5,14.5);
b8dc484b 104 TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
9f164762 105 TH1F *hTriggerable = new TH1F("hTriggerable"," Nb of triggerable tracks / evt",15,-0.5,14.5);
106 TH1F *hTriggered = new TH1F("hTriggered"," Nb of triggered tracks / evt",15,-0.5,14.5);
b8dc484b 107
c687aa97 108 // momentum resolution at vertex
109 histoFile->mkdir("momentumAtVertex","momentumAtVertex");
110 histoFile->cd("momentumAtVertex");
111
112 const Int_t pNBins = 30;
113 const Double_t pEdges[2] = {0., 300.};
114 const Int_t deltaPAtVtxNBins = 250;
115 const Double_t deltaPAtVtxEdges[2] = {-35., 15.};
116
117 TH1F *hResMomVertex = new TH1F("hResMomVertex"," delta P at vertex;#Delta_{p} (GeV/c)",deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
118
119 TH2D *hResMomVertexVsMom = new TH2D("hResMomVertexVsMom","#Delta_{p} at vertex versus p;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
120 TH2D *hResMomVertexVsMom_2_3_Deg = new TH2D("hResMomVertexVsMom_2_3_Deg","#Delta_{p} at vertex versus p for tracks between 2 and 3 degrees at absorber end;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
121 TH2D *hResMomVertexVsMom_3_10_Deg = new TH2D("hResMomVertexVsMom_3_10_Deg","#Delta_{p} at vertex versus p for tracks between 3 and 10 degrees at absorber end;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
122 TH2D *hResMomVertexVsMom_0_2_DegMC = new TH2D("hResMomVertexVsMom_0_2_DegMC","#Delta_{p} at vertex versus p for tracks with MC angle below 2 degrees;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins/10,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
123
e1fe81b8 124 TH2D *hResMomVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_0_2_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
125 TH2D *hResMomVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_2_3_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
126 TH2D *hResMomVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_3_10_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
c687aa97 127
128 TH2D *hResMomVertexVsAngle = new TH2D("hResMomVertexVsAngle","#Delta_{p} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{p} (GeV/c)",10,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
129 TH2D *hResMomVertexVsMCAngle = new TH2D("hResMomVertexVsMCAngle","#Delta_{p} at vertex versus MC angle;MC angle (Deg);#Delta_{p} (GeV/c)",10,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
130 TH3D *hResMomVertexVsAngleVsMom = new TH3D("hResMomVertexVsAngleVsMom","#Delta_{p} at vertex versus track position at absorber end converted to degrees versus momentum;p (GeV/c);angle (Deg);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],100,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
131
e1fe81b8 132 TGraphAsymmErrors* gMeanResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
c687aa97 133 gMeanResMomVertexVsMom->SetName("gMeanResMomVertexVsMom");
134 gMeanResMomVertexVsMom->SetTitle("<#Delta_{p}> at vertex versus p;p (GeV/c);<#Delta_{p}> (GeV/c)");
e1fe81b8 135 TGraphAsymmErrors* gMostProbResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
c687aa97 136 gMostProbResMomVertexVsMom->SetName("gMostProbResMomVertexVsMom");
137 gMostProbResMomVertexVsMom->SetTitle("Most probable #Delta_{p} at vertex versus p;p (GeV/c);Most prob. #Delta_{p} (GeV/c)");
e1fe81b8 138 TGraphAsymmErrors* gSigmaResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
c687aa97 139 gSigmaResMomVertexVsMom->SetName("gSigmaResMomVertexVsMom");
e1fe81b8 140 gSigmaResMomVertexVsMom->SetTitle("#sigma_{p}/p at vertex versus p;p (GeV/c);#sigma_{p}/p (%)");
c687aa97 141
142 // momentum resolution at first cluster
143 histoFile->mkdir("momentumAtFirstCluster","momentumAtFirstCluster");
144 histoFile->cd("momentumAtFirstCluster");
145
146 const Int_t deltaPAtFirstClNBins = 500;
147 const Double_t deltaPAtFirstClEdges[2] = {-25., 25.};
148
149 TH1F *hResMomFirstCluster = new TH1F("hResMomFirstCluster"," delta P at first cluster;#Delta_{p} (GeV/c)",deltaPAtFirstClNBins,deltaPAtFirstClEdges[0],deltaPAtFirstClEdges[1]);
150 TH2D *hResMomFirstClusterVsMom = new TH2D("hResMomFirstClusterVsMom","#Delta_{p} at first cluster versus p;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtFirstClNBins,deltaPAtFirstClEdges[0],deltaPAtFirstClEdges[1]);
151
e1fe81b8 152 TGraphAsymmErrors* gMeanResMomFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
c687aa97 153 gMeanResMomFirstClusterVsMom->SetName("gMeanResMomFirstClusterVsMom");
154 gMeanResMomFirstClusterVsMom->SetTitle("<#Delta_{p}> at first cluster versus p;p (GeV/c);<#Delta_{p}> (GeV/c)");
e1fe81b8 155 TGraphAsymmErrors* gSigmaResMomFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
c687aa97 156 gSigmaResMomFirstClusterVsMom->SetName("gSigmaResMomFirstClusterVsMom");
e1fe81b8 157 gSigmaResMomFirstClusterVsMom->SetTitle("#sigma_{p}/p at first cluster versus p;p (GeV/c);#sigma_{p}/p (%)");
158
159 // angular resolution at vertex
160 histoFile->mkdir("slopesAtVertex","slopesAtVertex");
161 histoFile->cd("slopesAtVertex");
162
163 const Int_t deltaSlopeAtVtxNBins = 500;
164 const Double_t deltaSlopeAtVtxEdges[2] = {-0.05, 0.05};
165
166 TH1F *hResSlopeXVertex = new TH1F("hResSlopeXVertex","#Delta_{slope_{X}} at vertex;#Delta_{slope_{X}}", deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
167 TH1F *hResSlopeYVertex = new TH1F("hResSlopeYVertex","#Delta_{slope_{Y}} at vertex;#Delta_{slope_{Y}}", deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
168 TH2D *hResSlopeXVertexVsMom = new TH2D("hResSlopeXVertexVsMom","#Delta_{slope_{X}} at vertex versus p;p (GeV/c);#Delta_{slope_{X}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
169 TH2D *hResSlopeYVertexVsMom = new TH2D("hResSlopeYVertexVsMom","#Delta_{slope_{Y}} at vertex versus p;p (GeV/c);#Delta_{slope_{Y}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
170
171 TH2D *hResSlopeXVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_0_2_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
172 TH2D *hResSlopeYVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_0_2_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
173 TH2D *hResSlopeXVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_2_3_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
174 TH2D *hResSlopeYVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_2_3_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
175 TH2D *hResSlopeXVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_3_10_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
176 TH2D *hResSlopeYVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_3_10_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
177
178 TH2D *hResSlopeXVertexVsAngle = new TH2D("hResSlopeXVertexVsAngle","#Delta_{slope_{X}} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{slope_{X}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
179 TH2D *hResSlopeYVertexVsAngle = new TH2D("hResSlopeYVertexVsAngle","#Delta_{slope_{Y}} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{slope_{Y}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
180 TH2D *hResSlopeXVertexVsMCAngle = new TH2D("hResSlopeXVertexVsMCAngle","#Delta_{slope_{X}} at vertex versus MC angle;MC angle (Deg);#Delta_{slope_{X}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
181 TH2D *hResSlopeYVertexVsMCAngle = new TH2D("hResSlopeYVertexVsMCAngle","#Delta_{slope_{Y}} at vertex versus MC angle;MC angle (Deg);#Delta_{slope_{Y}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
182
183 TGraphAsymmErrors* gMeanResSlopeXVertexVsMom = new TGraphAsymmErrors(pNBins);
184 gMeanResSlopeXVertexVsMom->SetName("gMeanResSlopeXVertexVsMom");
185 gMeanResSlopeXVertexVsMom->SetTitle("<#Delta_{slope_{X}}> at vertex versus p;p (GeV/c);<#Delta_{slope_{X}}>");
186 TGraphAsymmErrors* gSigmaResSlopeXVertexVsMom = new TGraphAsymmErrors(pNBins);
187 gSigmaResSlopeXVertexVsMom->SetName("gSigmaResSlopeXVertexVsMom");
188 gSigmaResSlopeXVertexVsMom->SetTitle("#sigma_{slope_{X}} at vertex versus p;p (GeV/c);#sigma_{slope_{X}}");
189 TGraphAsymmErrors* gMeanResSlopeYVertexVsMom = new TGraphAsymmErrors(pNBins);
190 gMeanResSlopeYVertexVsMom->SetName("gMeanResSlopeYVertexVsMom");
191 gMeanResSlopeYVertexVsMom->SetTitle("<#Delta_{slope_{Y}}> at vertex versus p;p (GeV/c);<#Delta_{slope_{Y}}>");
192 TGraphAsymmErrors* gSigmaResSlopeYVertexVsMom = new TGraphAsymmErrors(pNBins);
193 gSigmaResSlopeYVertexVsMom->SetName("gSigmaResSlopeYVertexVsMom");
194 gSigmaResSlopeYVertexVsMom->SetTitle("#sigma_{slope_{Y}} at vertex versus p;p (GeV/c);#sigma_{slope_{Y}}");
195
196 // angular resolution at first cluster
197 histoFile->mkdir("slopesAtFirstCluster","slopesAtFirstCluster");
198 histoFile->cd("slopesAtFirstCluster");
199
200 const Int_t deltaSlopeAtFirstClNBins = 500;
201 const Double_t deltaSlopeAtFirstClEdges[2] = {-0.01, 0.01};
202
203 TH1F *hResSlopeXFirstCluster = new TH1F("hResSlopeXFirstCluster","#Delta_{slope_{X}} at first cluster;#Delta_{slope_{X}}", deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
204 TH2D *hResSlopeXFirstClusterVsMom = new TH2D("hResSlopeXFirstClusterVsMom","#Delta_{slope_{X}} at first cluster versus p;p (GeV/c);#Delta_{slope_{X}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
205 TH1F *hResSlopeYFirstCluster = new TH1F("hResSlopeYFirstCluster","#Delta_{slope_{Y}} at first cluster;#Delta_{slope_{Y}}", deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
206 TH2D *hResSlopeYFirstClusterVsMom = new TH2D("hResSlopeYFirstClusterVsMom","#Delta_{slope_{Y}} at first cluster versus p;p (GeV/c);#Delta_{slope_{Y}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
207
208 TGraphAsymmErrors* gMeanResSlopeXFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
209 gMeanResSlopeXFirstClusterVsMom->SetName("gMeanResSlopeXFirstClusterVsMom");
210 gMeanResSlopeXFirstClusterVsMom->SetTitle("<#Delta_{slope_{X}}> at first cluster versus p;p (GeV/c);<#Delta_{slope_{X}}>");
211 TGraphAsymmErrors* gSigmaResSlopeXFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
212 gSigmaResSlopeXFirstClusterVsMom->SetName("gSigmaResSlopeXFirstClusterVsMom");
213 gSigmaResSlopeXFirstClusterVsMom->SetTitle("#sigma_{slope_{X}} at first cluster versus p;p (GeV/c);#sigma_{slope_{X}}");
214 TGraphAsymmErrors* gMeanResSlopeYFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
215 gMeanResSlopeYFirstClusterVsMom->SetName("gMeanResSlopeYFirstClusterVsMom");
216 gMeanResSlopeYFirstClusterVsMom->SetTitle("<#Delta_{slope_{Y}}> at first cluster versus p;p (GeV/c);<#Delta_{slope_{Y}}>");
217 TGraphAsymmErrors* gSigmaResSlopeYFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
218 gSigmaResSlopeYFirstClusterVsMom->SetName("gSigmaResSlopeYFirstClusterVsMom");
219 gSigmaResSlopeYFirstClusterVsMom->SetTitle("#sigma_{slope_{Y}} at first cluster versus p;p (GeV/c);#sigma_{slope_{Y}}");
220
221 // DCA resolution and MCS angular dispersion
222 histoFile->mkdir("DCA","DCA");
223 histoFile->cd("DCA");
224
225 const Int_t deltaPDCANBins = 500;
226 const Double_t deltaPDCAEdges[2] = {0., 1000.};
227 const Double_t deltaPMCSAngEdges[2] = {-0.5, 0.5};
228
229 TH1F *hPDCA = new TH1F("hPDCA","p #times DCA at vertex;p #times DCA (GeV #times cm)", deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
230 TH2D *hPDCAVsMom_2_3_Deg = new TH2D("hPDCAVsMom_2_3_Deg","p #times DCA versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times DCA (GeV #times cm)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
231 TH2D *hPDCAVsMom_3_10_Deg = new TH2D("hPDCAVsMom_3_10_Deg","p #times DCA versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);p #times DCA (GeV #times cm)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
232 TH2D *hPMCSAngVsMom_2_3_Deg = new TH2D("hPMCSAngVsMom_2_3_Deg","p #times #Delta#theta_{MCS} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times #Delta#theta_{MCS} (GeV)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPMCSAngEdges[0], deltaPMCSAngEdges[1]);
233 TH2D *hPMCSAngVsMom_3_10_Deg = new TH2D("hPMCSAngVsMom_3_10_Deg","p #times #Delta#theta_{MCS} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times #Delta#theta_{MCS} (GeV)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPMCSAngEdges[0], deltaPMCSAngEdges[1]);
234
235 TH2D *hPDCAVsPosAbsEnd_0_2_DegMC = new TH2D("hPDCAVsPosAbsEnd_0_2_DegMC","p #times DCA versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
236 TH2D *hPDCAVsPosAbsEnd_2_3_DegMC = new TH2D("hPDCAVsPosAbsEnd_2_3_DegMC","p #times DCA}versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
237 TH2D *hPDCAVsPosAbsEnd_3_10_DegMC = new TH2D("hPDCAVsPosAbsEnd_3_10_DegMC","p #times DCA versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
238
239 TH2D *hPDCAVsAngle = new TH2D("hPDCAVsAngle","p #times DCA versus track position at absorber end converted to degrees;angle (Deg);p #times DCA (GeV #times cm)",10,0.,10.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
240 TH2D *hPDCAVsMCAngle = new TH2D("hPDCAVsMCAngle","p #times DCA versus MC angle;MC angle (Deg);p #times DCA (GeV #times cm)",10,0.,10.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
241
242 TGraphAsymmErrors* gSigmaPDCAVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
243 gSigmaPDCAVsMom_2_3_Deg->SetName("gSigmaPDCAVsMom_2_3_Deg");
244 gSigmaPDCAVsMom_2_3_Deg->SetTitle("#sigma_{p #times DCA} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);#sigma_{p #times DCA} (GeV #times cm)");
245 TGraphAsymmErrors* gSigmaPDCAVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
246 gSigmaPDCAVsMom_3_10_Deg->SetName("gSigmaPDCAVsMom_3_10_Deg");
247 gSigmaPDCAVsMom_3_10_Deg->SetTitle("#sigma_{p #times DCA} versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);#sigma_{p #times DCA} (GeV #times cm)");
248 TGraphAsymmErrors* gMeanPMCSAngVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
249 gMeanPMCSAngVsMom_2_3_Deg->SetName("gMeanPMCSAngVsMom_2_3_Deg");
250 gMeanPMCSAngVsMom_2_3_Deg->SetTitle("<p #times #Delta#theta_{MCS}> versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);<p #times #Delta#theta_{MCS}> (GeV)");
251 TGraphAsymmErrors* gSigmaPMCSAngVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
252 gSigmaPMCSAngVsMom_2_3_Deg->SetName("gSigmaPMCSAngVsMom_2_3_Deg");
253 gSigmaPMCSAngVsMom_2_3_Deg->SetTitle("#sigma_{p #times #Delta#theta_{MCS}} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);#sigma_{p #times #Delta#theta_{MCS}} (GeV)");
254 TGraphAsymmErrors* gMeanPMCSAngVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
255 gMeanPMCSAngVsMom_3_10_Deg->SetName("gMeanPMCSAngVsMom_3_10_Deg");
256 gMeanPMCSAngVsMom_3_10_Deg->SetTitle("<p #times #Delta#theta_{MCS}> versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);<p #times #Delta#theta_{MCS}> (GeV)");
257 TGraphAsymmErrors* gSigmaPMCSAngVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
258 gSigmaPMCSAngVsMom_3_10_Deg->SetName("gSigmaPMCSAngVsMom_3_10_Deg");
259 gSigmaPMCSAngVsMom_3_10_Deg->SetTitle("#sigma_{p #times #Delta#theta_{MCS}} versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);#sigma_{p #times #Delta#theta_{MCS}} (GeV)");
260
261 // eta resolution at vertex
262 histoFile->mkdir("etaAtVertex","etaAtVertex");
263 histoFile->cd("etaAtVertex");
264
265 const Int_t deltaEtaAtVtxNBins = 500;
266 const Double_t deltaEtaAtVtxEdges[2] = {-0.5, 0.5};
267
268 TH1F *hResEtaVertex = new TH1F("hResEtaVertex","#Delta_{eta} at vertex;#Delta_{eta}", deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
269 TH2D *hResEtaVertexVsMom = new TH2D("hResEtaVertexVsMom","#Delta_{eta} at vertex versus p;p (GeV/c);#Delta_{eta}",2*pNBins,pEdges[0],pEdges[1], deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
270
271 TH2D *hResEtaVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_0_2_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
272 TH2D *hResEtaVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_2_3_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
273 TH2D *hResEtaVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_3_10_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
274
275 TH2D *hResEtaVertexVsAngle = new TH2D("hResEtaVertexVsAngle","#Delta_{eta} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{eta}",10,0.,10.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
276 TH2D *hResEtaVertexVsMCAngle = new TH2D("hResEtaVertexVsMCAngle","#Delta_{eta} at vertex versus MC angle;MC angle (Deg);#Delta_{eta}",10,0.,10.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
277
278 TGraphAsymmErrors* gMeanResEtaVertexVsMom = new TGraphAsymmErrors(pNBins);
279 gMeanResEtaVertexVsMom->SetName("gMeanResEtaVertexVsMom");
280 gMeanResEtaVertexVsMom->SetTitle("<#Delta_{eta}> at vertex versus p;p (GeV/c);<#Delta_{eta}>");
281 TGraphAsymmErrors* gSigmaResEtaVertexVsMom = new TGraphAsymmErrors(pNBins);
282 gSigmaResEtaVertexVsMom->SetName("gSigmaResEtaVertexVsMom");
283 gSigmaResEtaVertexVsMom->SetTitle("#sigma_{eta} at vertex versus p;p (GeV/c);#sigma_{eta}");
284
285 // phi resolution at vertex
286 histoFile->mkdir("phiAtVertex","phiAtVertex");
287 histoFile->cd("phiAtVertex");
288
289 const Int_t deltaPhiAtVtxNBins = 500;
290 const Double_t deltaPhiAtVtxEdges[2] = {-0.5, 0.5};
291
292 TH1F *hResPhiVertex = new TH1F("hResPhiVertex","#Delta_{phi} at vertex;#Delta_{phi}", deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
293 TH2D *hResPhiVertexVsMom = new TH2D("hResPhiVertexVsMom","#Delta_{phi} at vertex versus p;p (GeV/c);#Delta_{phi}",2*pNBins,pEdges[0],pEdges[1], deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
294
295 TH2D *hResPhiVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_0_2_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
296 TH2D *hResPhiVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_2_3_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
297 TH2D *hResPhiVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_3_10_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
298
299 TH2D *hResPhiVertexVsAngle = new TH2D("hResPhiVertexVsAngle","#Delta_{phi} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{phi}",10,0.,10.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
300 TH2D *hResPhiVertexVsMCAngle = new TH2D("hResPhiVertexVsMCAngle","#Delta_{phi} at vertex versus MC angle;MC angle (Deg);#Delta_{phi}",10,0.,10.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
301
302 TGraphAsymmErrors* gMeanResPhiVertexVsMom = new TGraphAsymmErrors(pNBins);
303 gMeanResPhiVertexVsMom->SetName("gMeanResPhiVertexVsMom");
304 gMeanResPhiVertexVsMom->SetTitle("<#Delta_{phi}> at vertex versus p;p (GeV/c);<#Delta_{phi}>");
305 TGraphAsymmErrors* gSigmaResPhiVertexVsMom = new TGraphAsymmErrors(pNBins);
306 gSigmaResPhiVertexVsMom->SetName("gSigmaResPhiVertexVsMom");
307 gSigmaResPhiVertexVsMom->SetTitle("#sigma_{phi} at vertex versus p;p (GeV/c);#sigma_{phi}");
c687aa97 308
309 // cluster resolution
310 histoFile->mkdir("clusters","clusters");
311 histoFile->cd("clusters");
312
f202486b 313 TH1F* hResidualXInCh[AliMUONConstants::NTrackingCh()];
314 TH1F* hResidualYInCh[AliMUONConstants::NTrackingCh()];
315 for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) {
316 hResidualXInCh[i] = new TH1F(Form("hResidualXInCh%d",i+1), Form("cluster-track residual-X distribution in chamber %d;#Delta_{X} (cm)",i+1), 1000, -1., 1.);
317 hResidualYInCh[i] = new TH1F(Form("hResidualYInCh%d",i+1), Form("cluster-track residual-Y distribution in chamber %d;#Delta_{Y} (cm)",i+1), 1000, -0.5, 0.5);
318 }
c687aa97 319
f202486b 320 TGraphErrors* gResidualXPerChMean = new TGraphErrors(AliMUONConstants::NTrackingCh());
321 gResidualXPerChMean->SetName("gResidualXPerChMean");
e1fe81b8 322 gResidualXPerChMean->SetTitle("cluster-trackRef residual-X per Ch: mean;chamber ID;<#Delta_{X}> (cm)");
f202486b 323 gResidualXPerChMean->SetMarkerStyle(kFullDotLarge);
324 TGraphErrors* gResidualYPerChMean = new TGraphErrors(AliMUONConstants::NTrackingCh());
325 gResidualYPerChMean->SetName("gResidualYPerChMean");
e1fe81b8 326 gResidualYPerChMean->SetTitle("cluster-trackRef residual-Y per Ch: mean;chamber ID;<#Delta_{Y}> (cm)");
f202486b 327 gResidualYPerChMean->SetMarkerStyle(kFullDotLarge);
328 TGraphErrors* gResidualXPerChSigma = new TGraphErrors(AliMUONConstants::NTrackingCh());
329 gResidualXPerChSigma->SetName("gResidualXPerChSigma");
e1fe81b8 330 gResidualXPerChSigma->SetTitle("cluster-trackRef residual-X per Ch: sigma;chamber ID;#sigma_{X} (cm)");
f202486b 331 gResidualXPerChSigma->SetMarkerStyle(kFullDotLarge);
332 TGraphErrors* gResidualYPerChSigma = new TGraphErrors(AliMUONConstants::NTrackingCh());
333 gResidualYPerChSigma->SetName("gResidualYPerChSigma");
e1fe81b8 334 gResidualYPerChSigma->SetTitle("cluster-trackRef residual-Y per Ch: sigma;chamber ID;#sigma_{Y} (cm)");
f202486b 335 gResidualYPerChSigma->SetMarkerStyle(kFullDotLarge);
9f164762 336
337 histoFile->mkdir("trigger");
338 histoFile->cd("trigger");
339 TH1F* hResidualTrigX11 = new TH1F("hResiudalTrigX11", "Residual X11", 100, -10., 10.);
340 TH1F* hResidualTrigY11 = new TH1F("hResiudalTrigY11", "Residual Y11", 100, -10., 10.);
341 TH1F* hResidualTrigSlopeY = new TH1F("hResiudalTrigSlopeY", "Residual Y slope", 100, -0.1, 0.1);
342 TH1F* hTriggerableMatchFailed = new TH1F("hTriggerableMatchFailed", "Triggerable multiplicity for events with no match", 15, -0.5, 14.5);
343
48217459 344
c687aa97 345 // ###################################### initialize ###################################### //
a99c3449 346 AliMUONRecoCheck rc(esdFileName, pathSim);
347
348 // load necessary data from OCDB
349 AliCDBManager::Instance()->SetDefaultStorage(ocdbPath);
350 AliCDBManager::Instance()->SetRun(rc.GetRunNumber());
c687aa97 351 if (!AliMUONCDB::LoadField()) return;
352 AliMUONTrackExtrap::SetField();
e1fe81b8 353 AliGeomManager::LoadGeometry();
354 if (!AliGeomManager::GetGeometry()) return;
a99c3449 355 AliMUONRecoParam* recoParam = AliMUONCDB::LoadRecoParam();
356 if (!recoParam) return;
e1fe81b8 357 AliMUONESDInterface::ResetTracker(recoParam);
a0bfad0e 358
f202486b 359 // get sigma cut from recoParam to associate clusters with TrackRefs in case the label are not used
a99c3449 360 Double_t sigmaCut = (recoParam->ImproveTracks()) ? recoParam->GetSigmaCutForImprovement() : recoParam->GetSigmaCutForTracking();
f202486b 361 // compute the mask of requested stations from recoParam
362 UInt_t requestedStationMask = 0;
363 for (Int_t i = 0; i < 5; i++) if (recoParam->RequestStation(i)) requestedStationMask |= ( 1 << i );
364 // get from recoParam whether a track need 2 chambers hit in the same station (4 or 5) or not to be reconstructible
365 Bool_t request2ChInSameSt45 = !recoParam->MakeMoreTrackCandidates();
48217459 366
367 Int_t nevents = rc.NumberOfEvents();
b8dc484b 368
edcc5b56 369 if (nevents < nEvent || nEvent < 0) nEvent = nevents;
b8dc484b 370
371 Int_t ievent;
372 Int_t nReconstructibleTracks = 0;
373 Int_t nReconstructedTracks = 0;
374 Int_t nReconstructibleTracksCheck = 0;
f202486b 375 AliMUONTrackParam *trackParam;
e1fe81b8 376 Double_t x1,y1,z1,slopex1,slopey1,pX1,pY1,pZ1,p1,pT1,eta1,phi1;
377 Double_t x2,y2,z2,slopex2,slopey2,pX2,pY2,pZ2,p2,pT2,eta2,phi2;
378 Double_t xAbs,yAbs,dAbs,aAbs,aMCS,aMC;
379 Double_t xDCA,yDCA,dca,pU;
380 Double_t aMCSMoy = 0., aMCS2Moy = 0., dMCSMoy = 0., dMCS2Moy = 0., adMCSMoy = 0.;
381 Int_t nMCS = 0;
48217459 382
c687aa97 383 // ###################################### fill histograms ###################################### //
48217459 384 for (ievent=0; ievent<nEvent; ievent++)
385 {
c687aa97 386 if ((ievent+1)%100 == 0) cout<<"\rEvent processing... "<<ievent+1<<flush;
b8dc484b 387
a99c3449 388 AliMUONVTrackStore* trackStore = rc.ReconstructedTracks(ievent, kFALSE);
f202486b 389 AliMUONVTrackStore* trackRefStore = rc.ReconstructibleTracks(ievent, requestedStationMask, request2ChInSameSt45);
b8dc484b 390
48217459 391 hReconstructible->Fill(trackRefStore->GetSize());
392 hReco->Fill(trackStore->GetSize());
393
394 nReconstructibleTracks += trackRefStore->GetSize();
395 nReconstructedTracks += trackStore->GetSize();
9f164762 396
397 AliMUONVTriggerTrackStore* triggerTrackRefStore = rc.TriggerableTracks(ievent);
398 AliMUONVTriggerTrackStore* triggerTrackStore = rc.TriggeredTracks(ievent);
399
400 hTriggerable->Fill(triggerTrackRefStore->GetSize());
401 hTriggered->Fill(triggerTrackStore->GetSize());
402
403 // loop over trigger trackRef
404 TIter nextTrig(triggerTrackRefStore->CreateIterator());
405 AliMUONTriggerTrack* triggerTrackRef;
406 Int_t nTriggerMatches = 0;
407 while ( ( triggerTrackRef = static_cast<AliMUONTriggerTrack*>(nextTrig()) ) )
408 {
409
410 AliMUONTriggerTrack* triggerTrackMatched = 0x0;
411
412 // loop over trackReco and look for compatible track
413 TIter nextTrig2(triggerTrackStore->CreateIterator());
414 AliMUONTriggerTrack* triggerTrackReco;
415 while ( ( triggerTrackReco = static_cast<AliMUONTriggerTrack*>(nextTrig2()) ) )
416 {
417
418 // check if trackReco is compatible with trackRef
419 if (triggerTrackReco->Match(*triggerTrackRef, sigmaCut)) {
420 triggerTrackMatched = triggerTrackReco;
421 nTriggerMatches++;
422 break;
423 }
424 }
425
426 if (triggerTrackMatched) { // tracking requirements verified, track is found
427 hResidualTrigX11->Fill( triggerTrackMatched->GetX11() - triggerTrackRef->GetX11() );
428 hResidualTrigY11->Fill( triggerTrackMatched->GetY11() - triggerTrackRef->GetY11() );
429 hResidualTrigSlopeY->Fill( triggerTrackMatched->GetSlopeY() - triggerTrackRef->GetSlopeY() );
430 }
431 } // loop on trigger track ref
432
433 if ( nTriggerMatches != triggerTrackStore->GetSize() )
434 hTriggerableMatchFailed->Fill(triggerTrackRefStore->GetSize());
48217459 435
f202486b 436 // loop over trackRef
48217459 437 TIter next(trackRefStore->CreateIterator());
438 AliMUONTrack* trackRef;
48217459 439 while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) )
440 {
48217459 441
9f164762 442 hTrackRefID->Fill(trackRef->GetUniqueID());
f202486b 443
444 AliMUONTrack* trackMatched = 0x0;
445 Int_t nMatchClusters = 0;
446
447 // loop over trackReco and look for compatible track
48217459 448 TIter next2(trackStore->CreateIterator());
449 AliMUONTrack* trackReco;
48217459 450 while ( ( trackReco = static_cast<AliMUONTrack*>(next2()) ) )
451 {
f202486b 452
453 // check if trackReco is compatible with trackRef
f202486b 454 if (trackReco->Match(*trackRef, sigmaCut, nMatchClusters)) {
455 trackMatched = trackReco;
f202486b 456 break;
61fed964 457 }
f202486b 458
48217459 459 }
460
f202486b 461 if (trackMatched) { // tracking requirements verified, track is found
48217459 462 nReconstructibleTracksCheck++;
f202486b 463 hNClusterComp->Fill(nMatchClusters);
c687aa97 464
465 // compute track position at the end of the absorber
466 AliMUONTrackParam trackParamAtAbsEnd(*((AliMUONTrackParam*)trackMatched->GetTrackParamAtCluster()->First()));
467 AliMUONTrackExtrap::ExtrapToZ(&trackParamAtAbsEnd, AliMUONConstants::AbsZEnd());
468 xAbs = trackParamAtAbsEnd.GetNonBendingCoor();
469 yAbs = trackParamAtAbsEnd.GetBendingCoor();
470 dAbs = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
471 aAbs = TMath::ATan(-dAbs/AliMUONConstants::AbsZEnd()) * TMath::RadToDeg();
e1fe81b8 472 pX2 = trackParamAtAbsEnd.Px();
473 pY2 = trackParamAtAbsEnd.Py();
474 pZ2 = trackParamAtAbsEnd.Pz();
475 pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
476 aMCS = TMath::ATan(-pT2/pZ2) * TMath::RadToDeg();
c687aa97 477
48217459 478 trackParam = trackRef->GetTrackParamAtVertex();
479 x1 = trackParam->GetNonBendingCoor();
480 y1 = trackParam->GetBendingCoor();
481 z1 = trackParam->GetZ();
e1fe81b8 482 slopex1 = trackParam->GetNonBendingSlope();
483 slopey1 = trackParam->GetBendingSlope();
48217459 484 pX1 = trackParam->Px();
485 pY1 = trackParam->Py();
486 pZ1 = trackParam->Pz();
487 p1 = trackParam->P();
f202486b 488 pT1 = TMath::Sqrt(pX1*pX1 + pY1*pY1);
c687aa97 489 aMC = TMath::ATan(-pT1/pZ1) * TMath::RadToDeg();
e1fe81b8 490 eta1 = TMath::Log(TMath::Tan(0.5*TMath::ATan(-pT1/pZ1)));
491 phi1 = TMath::Pi()+TMath::ATan2(-pY1, -pX1);
c687aa97 492
f202486b 493 trackParam = trackMatched->GetTrackParamAtVertex();
48217459 494 x2 = trackParam->GetNonBendingCoor();
495 y2 = trackParam->GetBendingCoor();
496 z2 = trackParam->GetZ();
e1fe81b8 497 slopex2 = trackParam->GetNonBendingSlope();
498 slopey2 = trackParam->GetBendingSlope();
48217459 499 pX2 = trackParam->Px();
500 pY2 = trackParam->Py();
501 pZ2 = trackParam->Pz();
502 p2 = trackParam->P();
f202486b 503 pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
e1fe81b8 504 eta2 = TMath::Log(TMath::Tan(0.5*TMath::ATan(-pT2/pZ2)));
505 phi2 = TMath::Pi()+TMath::ATan2(-pY2, -pX2);
48217459 506
e1fe81b8 507 AliMUONTrackParam trackParamAtDCA(*((AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First()));
508 pU = trackParamAtDCA.P();
509 AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(&trackParamAtDCA, z2);
510 xDCA = trackParamAtDCA.GetNonBendingCoor();
511 yDCA = trackParamAtDCA.GetBendingCoor();
512 dca = TMath::Sqrt(xDCA*xDCA + yDCA*yDCA);
513
48217459 514 hResMomVertex->Fill(p2-p1);
e1fe81b8 515 hResSlopeXVertex->Fill(slopex2-slopex1);
516 hResSlopeYVertex->Fill(slopey2-slopey1);
517 hPDCA->Fill(0.5*(p2+pU)*dca);
518 hResEtaVertex->Fill(eta2-eta1);
519 hResPhiVertex->Fill(phi2-phi1);
520 if (aMC >= aAbsLimits[0] && aMC <= aAbsLimits[1]) {
521 hResMomVertexVsMom->Fill(p1,p2-p1);
522 hResSlopeXVertexVsMom->Fill(p1,slopex2-slopex1);
523 hResSlopeYVertexVsMom->Fill(p1,slopey2-slopey1);
524 hResEtaVertexVsMom->Fill(p1,eta2-eta1);
525 hResPhiVertexVsMom->Fill(p1,phi2-phi1);
526 }
c687aa97 527 hResMomVertexVsAngleVsMom->Fill(p1,aAbs,p2-p1);
e1fe81b8 528 if (aAbs > 2. && aAbs < 3.) {
529 hResMomVertexVsMom_2_3_Deg->Fill(p1,p2-p1);
530 hPDCAVsMom_2_3_Deg->Fill(p1,0.5*(p2+pU)*dca);
531 hPMCSAngVsMom_2_3_Deg->Fill(p1,0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
532 }
533 else if (aAbs >= 3. && aAbs < 10.) {
534 hResMomVertexVsMom_3_10_Deg->Fill(p1,p2-p1);
535 hPDCAVsMom_3_10_Deg->Fill(p1,0.5*(p2+pU)*dca);
536 hPMCSAngVsMom_3_10_Deg->Fill(p1,0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
537 aMCSMoy += 0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad();
538 aMCS2Moy += (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad()) * (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
539 dMCSMoy += 0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd());
540 dMCS2Moy += (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd())) * (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd()));
541 adMCSMoy += (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad()) * (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd()));
542 nMCS++;
543 }
c687aa97 544 if (aMC < 2.) {
545 hResMomVertexVsMom_0_2_DegMC->Fill(p1,p2-p1);
546 hResMomVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,p2-p1);
e1fe81b8 547 hResSlopeXVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,slopex2-slopex1);
548 hResSlopeYVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,slopey2-slopey1);
549 hPDCAVsPosAbsEnd_0_2_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
550 hResEtaVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,eta2-eta1);
551 hResPhiVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,phi2-phi1);
552 }
553 else if (aMC >= 2. && aMC < 3) {
554 hResMomVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,p2-p1);
555 hResSlopeXVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,slopex2-slopex1);
556 hResSlopeYVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,slopey2-slopey1);
557 hPDCAVsPosAbsEnd_2_3_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
558 hResEtaVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,eta2-eta1);
559 hResPhiVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,phi2-phi1);
560 }
561 else if (aMC >= 3. && aMC < 10.) {
562 hResMomVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,p2-p1);
563 hResSlopeXVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,slopex2-slopex1);
564 hResSlopeYVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,slopey2-slopey1);
565 hPDCAVsPosAbsEnd_3_10_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
566 hResEtaVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,eta2-eta1);
567 hResPhiVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,phi2-phi1);
c687aa97 568 }
c687aa97 569 hResMomVertexVsAngle->Fill(aAbs,p2-p1);
e1fe81b8 570 hResSlopeXVertexVsAngle->Fill(aAbs,slopex2-slopex1);
571 hResSlopeYVertexVsAngle->Fill(aAbs,slopey2-slopey1);
572 hPDCAVsAngle->Fill(aAbs,0.5*(p2+pU)*dca);
573 hResEtaVertexVsAngle->Fill(aAbs,eta2-eta1);
574 hResPhiVertexVsAngle->Fill(aAbs,phi2-phi1);
c687aa97 575 hResMomVertexVsMCAngle->Fill(aMC,p2-p1);
e1fe81b8 576 hResSlopeXVertexVsMCAngle->Fill(aMC,slopex2-slopex1);
577 hResSlopeYVertexVsMCAngle->Fill(aMC,slopey2-slopey1);
578 hPDCAVsMCAngle->Fill(aMC,0.5*(p2+pU)*dca);
579 hResEtaVertexVsMCAngle->Fill(aMC,eta2-eta1);
580 hResPhiVertexVsMCAngle->Fill(aMC,phi2-phi1);
f202486b 581
96ebe67e 582 trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
48217459 583 x1 = trackParam->GetNonBendingCoor();
584 y1 = trackParam->GetBendingCoor();
585 z1 = trackParam->GetZ();
e1fe81b8 586 slopex1 = trackParam->GetNonBendingSlope();
587 slopey1 = trackParam->GetBendingSlope();
48217459 588 pX1 = trackParam->Px();
589 pY1 = trackParam->Py();
590 pZ1 = trackParam->Pz();
591 p1 = trackParam->P();
f202486b 592 pT1 = TMath::Sqrt(pX1*pX1 + pY1*pY1);
593
f202486b 594 trackParam = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First();
48217459 595 x2 = trackParam->GetNonBendingCoor();
596 y2 = trackParam->GetBendingCoor();
597 z2 = trackParam->GetZ();
e1fe81b8 598 slopex2 = trackParam->GetNonBendingSlope();
599 slopey2 = trackParam->GetBendingSlope();
48217459 600 pX2 = trackParam->Px();
601 pY2 = trackParam->Py();
602 pZ2 = trackParam->Pz();
603 p2 = trackParam->P();
f202486b 604 pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
48217459 605
96ebe67e 606 hResMomFirstCluster->Fill(p2-p1);
f202486b 607 hResMomFirstClusterVsMom->Fill(p1,p2-p1);
608
e1fe81b8 609 hResSlopeXFirstCluster->Fill(slopex2-slopex1);
610 hResSlopeYFirstCluster->Fill(slopey2-slopey1);
611 hResSlopeXFirstClusterVsMom->Fill(p1,slopex2-slopex1);
612 hResSlopeYFirstClusterVsMom->Fill(p1,slopey2-slopey1);
613
f202486b 614 // Fill residuals
615 // Loop over clusters of first track
616 AliMUONTrackParam* trackParamAtCluster1 = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First();
617 while (trackParamAtCluster1) {
618 AliMUONVCluster* cluster1 = trackParamAtCluster1->GetClusterPtr();
619 AliMUONTrackParam* trackParamAtCluster2 = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
620 while (trackParamAtCluster2) {
621 AliMUONVCluster* cluster2 = trackParamAtCluster2->GetClusterPtr();
622 if (cluster1->GetDetElemId() == cluster2->GetDetElemId()) {
623 hResidualXInCh[cluster1->GetChamberId()]->Fill(cluster1->GetX() - cluster2->GetX());
624 hResidualYInCh[cluster1->GetChamberId()]->Fill(cluster1->GetY() - cluster2->GetY());
625 break;
626 }
627 trackParamAtCluster2 = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->After(trackParamAtCluster2);
628 }
629 trackParamAtCluster1 = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->After(trackParamAtCluster1);
630 }
631
b8dc484b 632 }
f202486b 633
b8dc484b 634 } // end loop track ref.
b8dc484b 635
48217459 636 } // end loop on event
c687aa97 637 cout<<"\rEvent processing... "<<nevents<<" done"<<endl;
48217459 638
c687aa97 639 // ###################################### compute stuff ###################################### //
e1fe81b8 640 cout<<"\nWhen not specified, resolution at vertex is computed for ";
641 if (absorberRegion == 1) cout<<"tracks in the absorber region [2,3] deg."<<endl;
642 else if (absorberRegion == 2) cout<<"tracks in the absorber region [3,10] deg."<<endl;
643 else cout<<"all tracks"<<endl;
644
c687aa97 645 // compute momentum resolution at vertex versus p
e1fe81b8 646 TF1 *f2 = new TF1("f2",langaufun,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1],4);
c687aa97 647 Int_t rebinFactorX = TMath::Max(hResMomVertexVsMom->GetNbinsX()/pNBins, 1);
648 for (Int_t i = rebinFactorX; i <= hResMomVertexVsMom->GetNbinsX(); i+=rebinFactorX) {
649 cout<<"\rFitting momentum residuals at vertex... "<<i/rebinFactorX<<"/"<<pNBins<<flush;
650 TH1D *tmp = hResMomVertexVsMom->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
c687aa97 651 f2->SetParameters(0.2,0.,(Double_t)tmp->GetEntries(),1.);
652 tmp->Fit("f2","WWNQ");
653 Double_t fwhm = f2->GetParameter(0);
654 Double_t sigma = f2->GetParameter(3);
655 Double_t sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
656 Int_t rebin = TMath::Max(Int_t(0.5*sigmaP/tmp->GetBinWidth(1)),1);
657 while (deltaPAtVtxNBins%rebin!=0) rebin--;
658 tmp->Rebin(rebin);
659 tmp->Fit("f2","NQ");
660 fwhm = f2->GetParameter(0);
661 sigma = f2->GetParameter(3);
662 sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
663 Double_t fwhmErr = f2->GetParError(0);
664 Double_t sigmaErr = f2->GetParError(3);
665 Double_t sigmaPErr = TMath::Sqrt(sigma*sigma*sigmaErr*sigmaErr + fwhm*fwhm*fwhmErr*fwhmErr/(64.*log(2.)*log(2.))) / sigmaP;
e1fe81b8 666 hResMomVertexVsMom->GetXaxis()->SetRange(i-rebinFactorX+1,i);
667 Double_t p = hResMomVertexVsMom->GetMean();
668 hResMomVertexVsMom->GetXaxis()->SetRange();
669 Double_t pErr[2] = {p-hResMomVertexVsMom->GetBinLowEdge(i-rebinFactorX+1), hResMomVertexVsMom->GetBinLowEdge(i+1)-p};
670 gMeanResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, tmp->GetMean());
671 gMeanResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], tmp->GetMeanError(), tmp->GetMeanError());
672 gMostProbResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, -f2->GetParameter(1));
673 gMostProbResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], f2->GetParError(1), f2->GetParError(1));
674 gSigmaResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, 100.*sigmaP/p);
675 gSigmaResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], 100.*sigmaPErr/p, 100.*sigmaPErr/p);
f202486b 676 delete tmp;
677 }
c687aa97 678 cout<<"\rFitting momentum residuals at vertex... "<<pNBins<<"/"<<pNBins<<endl;
679
e1fe81b8 680 // compute momentum relative resolution at first cluster versus p
681 FitGausResVsMom(hResMomFirstClusterVsMom, pNBins, 0., 1., "momentum residuals at first cluster", gMeanResMomFirstClusterVsMom, gSigmaResMomFirstClusterVsMom);
c687aa97 682 rebinFactorX = TMath::Max(hResMomFirstClusterVsMom->GetNbinsX()/pNBins, 1);
683 for (Int_t i = rebinFactorX; i <= hResMomFirstClusterVsMom->GetNbinsX(); i+=rebinFactorX) {
e1fe81b8 684 Double_t x,y;
685 gSigmaResMomFirstClusterVsMom->GetPoint(i/rebinFactorX-1, x, y);
686 gSigmaResMomFirstClusterVsMom->SetPoint(i/rebinFactorX-1, x, 100.*y/x);
687 gSigmaResMomFirstClusterVsMom->SetPointEYlow(i/rebinFactorX-1, 100.*gSigmaResMomFirstClusterVsMom->GetErrorYlow(i/rebinFactorX-1)/x);
688 gSigmaResMomFirstClusterVsMom->SetPointEYhigh(i/rebinFactorX-1, 100.*gSigmaResMomFirstClusterVsMom->GetErrorYhigh(i/rebinFactorX-1)/x);
f202486b 689 }
690
e1fe81b8 691 // compute slopeX resolution at vertex versus p
692 FitGausResVsMom(hResSlopeXVertexVsMom, pNBins, 0., 2.e-3, "slopeX residuals at vertex", gMeanResSlopeXVertexVsMom, gSigmaResSlopeXVertexVsMom);
693
694 // compute slopeY resolution at vertex versus p
695 FitGausResVsMom(hResSlopeYVertexVsMom, pNBins, 0., 2.e-3, "slopeY residuals at vertex", gMeanResSlopeYVertexVsMom, gSigmaResSlopeYVertexVsMom);
696
697 // compute slopeX resolution at first cluster versus p
698 FitGausResVsMom(hResSlopeXFirstClusterVsMom, pNBins, 0., 3.e-4, "slopeX residuals at first cluster", gMeanResSlopeXFirstClusterVsMom, gSigmaResSlopeXFirstClusterVsMom);
699
700 // compute slopeY resolution at first cluster versus p
701 FitGausResVsMom(hResSlopeYFirstClusterVsMom, pNBins, 0., 2.e-4, "slopeY residuals at first cluster", gMeanResSlopeYFirstClusterVsMom, gSigmaResSlopeYFirstClusterVsMom);
702
703 // compute p*DCA resolution in the region [2,3] deg at absorber end
704 FitPDCAVsMom(hPDCAVsMom_2_3_Deg, pNBins, "p*DCA (tracks in [2,3] deg.)", gSigmaPDCAVsMom_2_3_Deg);
705
706 // compute p*DCA resolution in the region [3,10] deg at absorber end
707 FitPDCAVsMom(hPDCAVsMom_3_10_Deg, pNBins, "p*DCA (tracks in [3,10] deg.)", gSigmaPDCAVsMom_3_10_Deg);
708
709 // compute MCS angular dispersion in the region [2,3] deg at absorber end
710 FitGausResVsMom(hPMCSAngVsMom_2_3_Deg, pNBins, 0., 2.e-3, "p*MCSAngle (tracks in [2,3] deg.)", gMeanPMCSAngVsMom_2_3_Deg, gSigmaPMCSAngVsMom_2_3_Deg);
711
712 // compute MCS angular dispersion in the region [3,10] deg at absorber end
713 FitGausResVsMom(hPMCSAngVsMom_3_10_Deg, pNBins, 0., 2.e-3, "p*MCSAngle (tracks in [3,10] deg.)", gMeanPMCSAngVsMom_3_10_Deg, gSigmaPMCSAngVsMom_3_10_Deg);
714
715 // compute eta resolution at vertex versus p
716 FitGausResVsMom(hResEtaVertexVsMom, pNBins, 0., 0.1, "eta residuals at vertex", gMeanResEtaVertexVsMom, gSigmaResEtaVertexVsMom);
717
718 // compute phi resolution at vertex versus p
719 FitGausResVsMom(hResPhiVertexVsMom, pNBins, 0., 0.01, "phi residuals at vertex", gMeanResPhiVertexVsMom, gSigmaResPhiVertexVsMom);
720
721 // compute cluster-track residual mean and dispersion
f202486b 722 for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) {
723 hResidualXInCh[i]->GetXaxis()->SetRangeUser(-3.*hResidualXInCh[i]->GetRMS(), 3.*hResidualXInCh[i]->GetRMS());
724 gResidualXPerChMean->SetPoint(i, i+1, hResidualXInCh[i]->GetMean());
725 gResidualXPerChMean->SetPointError(i, 0., hResidualXInCh[i]->GetMeanError());
726 gResidualXPerChSigma->SetPoint(i, i+1, hResidualXInCh[i]->GetRMS());
727 gResidualXPerChSigma->SetPointError(i, 0., hResidualXInCh[i]->GetRMSError());
728 hResidualXInCh[i]->GetXaxis()->SetRange(0,0);
729 hResidualYInCh[i]->GetXaxis()->SetRangeUser(-3.*hResidualYInCh[i]->GetRMS(), 3.*hResidualYInCh[i]->GetRMS());
730 gResidualYPerChMean->SetPoint(i, i+1, hResidualYInCh[i]->GetMean());
731 gResidualYPerChMean->SetPointError(i, 0., hResidualYInCh[i]->GetMeanError());
732 gResidualYPerChSigma->SetPoint(i, i+1, hResidualYInCh[i]->GetRMS());
733 gResidualYPerChSigma->SetPointError(i, 0., hResidualYInCh[i]->GetRMSError());
734 hResidualYInCh[i]->GetXaxis()->SetRange(0,0);
735 }
736
c687aa97 737 // ###################################### display histograms ###################################### //
e1fe81b8 738 // diplay momentum residuals
739 TCanvas* cResMom = DrawVsAng("cResMom", "momentum residual at vertex in 3 angular regions", hResMomVertex, hResMomVertexVsAngle);
740 TCanvas* cResMomMC = DrawVsAng("cResMomMC", "momentum residual at vertex in 3 MC angular regions", hResMomVertex, hResMomVertexVsMCAngle);
741 TCanvas* cResMomVsPos = DrawVsPos("cResMomVsPos", "momentum residual at vertex versus position at absorber end in 3 MC angular regions",
742 hResMomVertexVsPosAbsEnd_0_2_DegMC, hResMomVertexVsPosAbsEnd_2_3_DegMC, hResMomVertexVsPosAbsEnd_3_10_DegMC);
743 TCanvas* cResMom_2_3_Deg = DrawResMomVsMom("cResMom_2_3_Deg", "momentum residual for tracks between 2 and 3 degrees",
744 hResMomVertexVsMom_2_3_Deg, 10, f2, "momentum residuals at vertex (tracks in [2,3] deg.)");
745 TCanvas* cResMom_3_10_Deg = DrawResMomVsMom("cResMom_3_10_Deg", "momentum residual for tracks between 3 and 10 degrees",
746 hResMomVertexVsMom_3_10_Deg, 10, f2, "momentum residuals at vertex (tracks in [3,10] deg.)");
747 TCanvas* cResMom_0_2_DegMC = DrawResMomVsMom("cResMom_0_2_DegMC", "momentum residuals for tracks with MC angle < 2 degrees", hResMomVertexVsMom_0_2_DegMC, 5);
748
749 // diplay slopeX residuals
750 TCanvas* cResSlopeX = DrawVsAng("cResSlopeX", "slope_{X} residual at vertex in 3 angular regions", hResSlopeXVertex, hResSlopeXVertexVsAngle);
751 TCanvas* cResSlopeXMC = DrawVsAng("cResSlopeXMC", "slope_{X} residual at vertex in 3 MC angular regions", hResSlopeXVertex, hResSlopeXVertexVsMCAngle);
752 TCanvas* cResSlopeXVsPos = DrawVsPos("cResSlopeXVsPos", "slope_{X} residual at vertex versus position at absorber end in 3 MC angular regions",
753 hResSlopeXVertexVsPosAbsEnd_0_2_DegMC, hResSlopeXVertexVsPosAbsEnd_2_3_DegMC, hResSlopeXVertexVsPosAbsEnd_3_10_DegMC);
754
755 // diplay slopeY residuals
756 TCanvas* cResSlopeY = DrawVsAng("cResSlopeY", "slope_{Y} residual at vertex in 3 angular regions", hResSlopeYVertex, hResSlopeYVertexVsAngle);
757 TCanvas* cResSlopeYMC = DrawVsAng("cResSlopeYMC", "slope_{Y} residual at vertex in 3 MC angular regions", hResSlopeYVertex, hResSlopeYVertexVsMCAngle);
758 TCanvas* cResSlopeYVsPos = DrawVsPos("cResSlopeYVsPos", "slope_{Y} residual at vertex versus position at absorber end in 3 MC angular regions",
759 hResSlopeYVertexVsPosAbsEnd_0_2_DegMC, hResSlopeYVertexVsPosAbsEnd_2_3_DegMC, hResSlopeYVertexVsPosAbsEnd_3_10_DegMC);
760
761 // diplay P*DCA
762 TCanvas* cPDCA = DrawVsAng("cPDCA", "p #times DCA in 3 angular regions", hPDCA, hPDCAVsAngle);
763 TCanvas* cPDCAMC = DrawVsAng("cPDCAMC", "p #times DCA in 3 MC angular regions", hPDCA, hPDCAVsMCAngle);
764 TCanvas* cPDCAVsPos = DrawVsPos("cPDCAVsPos", "p #times DCA versus position at absorber end in 3 MC angular regions",
765 hPDCAVsPosAbsEnd_0_2_DegMC, hPDCAVsPosAbsEnd_2_3_DegMC, hPDCAVsPosAbsEnd_3_10_DegMC);
766
767 // diplay eta residuals
768 TCanvas* cResEta = DrawVsAng("cResEta", "eta residual at vertex in 3 angular regions", hResEtaVertex, hResEtaVertexVsAngle);
769 TCanvas* cResEtaMC = DrawVsAng("cResEtaMC", "eta residual at vertex in 3 MC angular regions", hResEtaVertex, hResEtaVertexVsMCAngle);
770 TCanvas* cResEtaVsPos = DrawVsPos("cResEtaVsPos", "eta residual at vertex versus position at absorber end in 3 MC angular regions",
771 hResEtaVertexVsPosAbsEnd_0_2_DegMC, hResEtaVertexVsPosAbsEnd_2_3_DegMC, hResEtaVertexVsPosAbsEnd_3_10_DegMC);
772
773 // diplay phi residuals
774 TCanvas* cResPhi = DrawVsAng("cResPhi", "phi residual at vertex in 3 angular regions", hResPhiVertex, hResPhiVertexVsAngle);
775 TCanvas* cResPhiMC = DrawVsAng("cResPhiMC", "phi residual at vertex in 3 MC angular regions", hResPhiVertex, hResPhiVertexVsMCAngle);
776 TCanvas* cResPhiVsPos = DrawVsPos("cResPhiVsPos", "phi residual at vertex versus position at absorber end in 3 MC angular regions",
777 hResPhiVertexVsPosAbsEnd_0_2_DegMC, hResPhiVertexVsPosAbsEnd_2_3_DegMC, hResPhiVertexVsPosAbsEnd_3_10_DegMC);
c687aa97 778
779 // ###################################### save histogram ###################################### //
b8dc484b 780 histoFile->Write();
c687aa97 781
782 histoFile->cd("momentumAtVertex");
783 gMeanResMomVertexVsMom->Write();
784 gMostProbResMomVertexVsMom->Write();
785 gSigmaResMomVertexVsMom->Write();
e1fe81b8 786 cResMom->Write();
787 cResMomMC->Write();
788 cResMomVsPos->Write();
789 cResMom_2_3_Deg->Write();
790 cResMom_3_10_Deg->Write();
791 cResMom_0_2_DegMC->Write();
792
793 histoFile->cd("slopesAtVertex");
794 gMeanResSlopeXVertexVsMom->Write();
795 gMeanResSlopeYVertexVsMom->Write();
796 gSigmaResSlopeXVertexVsMom->Write();
797 gSigmaResSlopeYVertexVsMom->Write();
798 cResSlopeX->Write();
799 cResSlopeY->Write();
800 cResSlopeXMC->Write();
801 cResSlopeYMC->Write();
802 cResSlopeXVsPos->Write();
803 cResSlopeYVsPos->Write();
804
805 histoFile->cd("DCA");
806 gSigmaPDCAVsMom_2_3_Deg->Write();
807 gSigmaPDCAVsMom_3_10_Deg->Write();
808 gMeanPMCSAngVsMom_2_3_Deg->Write();
809 gSigmaPMCSAngVsMom_2_3_Deg->Write();
810 gMeanPMCSAngVsMom_3_10_Deg->Write();
811 gSigmaPMCSAngVsMom_3_10_Deg->Write();
812 cPDCA->Write();
813 cPDCAMC->Write();
814 cPDCAVsPos->Write();
815
816 histoFile->cd("etaAtVertex");
817 gMeanResEtaVertexVsMom->Write();
818 gSigmaResEtaVertexVsMom->Write();
819 cResEta->Write();
820 cResEtaMC->Write();
821 cResEtaVsPos->Write();
822
823 histoFile->cd("phiAtVertex");
824 gMeanResPhiVertexVsMom->Write();
825 gSigmaResPhiVertexVsMom->Write();
826 cResPhi->Write();
827 cResPhiMC->Write();
828 cResPhiVsPos->Write();
c687aa97 829
830 histoFile->cd("momentumAtFirstCluster");
831 gMeanResMomFirstClusterVsMom->Write();
832 gSigmaResMomFirstClusterVsMom->Write();
833
e1fe81b8 834 histoFile->cd("slopesAtFirstCluster");
835 gMeanResSlopeXFirstClusterVsMom->Write();
836 gMeanResSlopeYFirstClusterVsMom->Write();
837 gSigmaResSlopeXFirstClusterVsMom->Write();
838 gSigmaResSlopeYFirstClusterVsMom->Write();
839
c687aa97 840 histoFile->cd("clusters");
f202486b 841 gResidualXPerChMean->Write();
842 gResidualXPerChSigma->Write();
843 gResidualYPerChMean->Write();
844 gResidualYPerChSigma->Write();
c687aa97 845
b8dc484b 846 histoFile->Close();
c687aa97 847
e1fe81b8 848 // ###################################### clean memory ###################################### //
849 delete cResMom;
850 delete cResMomMC;
851 delete cResMomVsPos;
852 delete cResMom_2_3_Deg;
853 delete cResMom_3_10_Deg;
854 delete cResMom_0_2_DegMC;
855 delete cResSlopeX;
856 delete cResSlopeY;
857 delete cResSlopeXMC;
858 delete cResSlopeYMC;
859 delete cResSlopeXVsPos;
860 delete cResSlopeYVsPos;
861 delete cPDCA;
862 delete cPDCAMC;
863 delete cPDCAVsPos;
864 delete cResEta;
865 delete cResEtaMC;
866 delete cResEtaVsPos;
867 delete cResPhi;
868 delete cResPhiMC;
869 delete cResPhiVsPos;
c687aa97 870
e1fe81b8 871 // ###################################### print statistics ###################################### //
872 printf("\n");
873 printf("nb of reconstructible tracks: %d \n", nReconstructibleTracks);
874 printf("nb of reconstructed tracks: %d \n", nReconstructedTracks);
875 printf("nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
876
877 aMCSMoy /= (Double_t) nMCS;
878 aMCS2Moy /= (Double_t) nMCS;
879 dMCSMoy /= (Double_t) nMCS;
880 dMCS2Moy /= (Double_t) nMCS;
881 adMCSMoy /= (Double_t) nMCS;
882 Double_t sigma2_ThetaMCS = aMCS2Moy - aMCSMoy*aMCSMoy;
883 Double_t sigma2_PosMCS = dMCS2Moy - dMCSMoy*dMCSMoy;
884 Double_t cov_ThetaPosMCS = - (adMCSMoy - aMCSMoy*dMCSMoy);
885 printf("\nmultiple scattering of tracks between 3 and 10 deg. at absorber end:\n");
886 printf(" sigma_ThetaMCS = %f\n", TMath::Sqrt(sigma2_ThetaMCS));
887 printf(" sigma_PosMCS = %f\n", TMath::Sqrt(sigma2_PosMCS));
888 printf(" cov_ThetaPosMCS = %f\n", cov_ThetaPosMCS);
889 printf(" --> sigma_DCA = %f\n", TMath::Sqrt(AliMUONConstants::AbsZEnd()*AliMUONConstants::AbsZEnd()*sigma2_ThetaMCS
890 - 2.*AliMUONConstants::AbsZEnd()*cov_ThetaPosMCS + sigma2_PosMCS));
891 printf("\n");
c687aa97 892}
893
894//------------------------------------------------------------------------------------
895Double_t langaufun(Double_t *x, Double_t *par) {
896
897 //Fit parameters:
898 //par[0]=Width (scale) parameter of Landau density
899 //par[1]=Most Probable (MP, location) parameter of Landau density
900 //par[2]=Total area (integral -inf to inf, normalization constant)
901 //par[3]=Width (sigma) of convoluted Gaussian function
902 //
903 //In the Landau distribution (represented by the CERNLIB approximation),
904 //the maximum is located at x=-0.22278298 with the location parameter=0.
905 //This shift is corrected within this function, so that the actual
906 //maximum is identical to the MP parameter.
907
908 // Numeric constants
909 Double_t invsq2pi = 0.3989422804014; // (2 pi)^(-1/2)
910 Double_t mpshift = -0.22278298; // Landau maximum location
911
912 // Control constants
913 Double_t np = 100.0; // number of convolution steps
914 Double_t sc = 5.0; // convolution extends to +-sc Gaussian sigmas
915
916 // Variables
917 Double_t xx;
918 Double_t mpc;
919 Double_t fland;
920 Double_t sum = 0.0;
921 Double_t xlow,xupp;
922 Double_t step;
923 Double_t i;
924
925
926 // MP shift correction
927 mpc = par[1] - mpshift * par[0];
928
929 // Range of convolution integral
930 xlow = x[0] - sc * par[3];
931 xupp = x[0] + sc * par[3];
932
933 step = (xupp-xlow) / np;
934
935 // Convolution integral of Landau and Gaussian by sum
936 for(i=1.0; i<=np/2; i++) {
937 xx = xlow + (i-.5) * step;
938 //change x -> -x because the tail of the Landau is at the left here...
939 fland = TMath::Landau(-xx,mpc,par[0]) / par[0];
940 sum += fland * TMath::Gaus(x[0],xx,par[3]);
941
942 //change x -> -x because the tail of the Landau is at the left here...
943 xx = xupp - (i-.5) * step;
944 fland = TMath::Landau(-xx,mpc,par[0]) / par[0];
945 sum += fland * TMath::Gaus(x[0],xx,par[3]);
946 }
947
948 return (par[2] * step * sum * invsq2pi / par[3]);
b8dc484b 949}
950
e1fe81b8 951//------------------------------------------------------------------------------------
952void FitGausResVsMom(TH2* h, Int_t nBins, const Double_t mean0, const Double_t sigma0,
953 const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma)
954{
955 /// generic function to fit residuals versus momentum with a gaussian
956 static TF1* fGaus = 0x0;
957 if (!fGaus) fGaus = new TF1("fGaus","gaus");
958
959 Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
960 for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
961 cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
962 TH1D *tmp = h->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
963 fGaus->SetParameters(tmp->GetEntries(), mean0, sigma0);
964 tmp->Fit("fGaus","WWNQ");
965 Int_t rebin = TMath::Max(Int_t(0.5*fGaus->GetParameter(2)/tmp->GetBinWidth(1)),1);
966 while (tmp->GetNbinsX()%rebin!=0) rebin--;
967 tmp->Rebin(rebin);
968 tmp->Fit("fGaus","NQ");
969 h->GetXaxis()->SetRange(i-rebinFactorX+1,i);
970 Double_t p = h->GetMean();
971 h->GetXaxis()->SetRange();
972 Double_t pErr[2] = {p-h->GetBinLowEdge(i-rebinFactorX+1), h->GetBinLowEdge(i+1)-p};
973 gMean->SetPoint(i/rebinFactorX-1, p, fGaus->GetParameter(1));
974 gMean->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fGaus->GetParError(1), fGaus->GetParError(1));
975 gSigma->SetPoint(i/rebinFactorX-1, p, fGaus->GetParameter(2));
976 gSigma->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fGaus->GetParError(2), fGaus->GetParError(2));
977 delete tmp;
978 }
979 cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
980}
981
982//------------------------------------------------------------------------------------
983void FitPDCAVsMom(TH2* h, Int_t nBins, const char* fitting, TGraphAsymmErrors* gSigma)
984{
985 /// generic function to fit p*DCA distributions
986 static TF1* fPGaus = 0x0;
987 if (!fPGaus) fPGaus = new TF1("fPGaus","x*gaus");
988
989 Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
990 for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
991 cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
992 TH1D *tmp = h->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
993 fPGaus->SetParameters(1.,-100.,100.);
994 Int_t rebin = 50.*(tmp->GetNbinsX()/(tmp->GetBinLowEdge(tmp->GetNbinsX()+1)-tmp->GetBinLowEdge(1)));
995 while (tmp->GetNbinsX()%rebin!=0) rebin--;
996 tmp->Rebin(rebin);
997 tmp->Fit("fPGaus","NQ");
998 h->GetXaxis()->SetRange(i-rebinFactorX+1,i);
999 Double_t p = h->GetMean();
1000 h->GetXaxis()->SetRange();
1001 Double_t pErr[2] = {p-h->GetBinLowEdge(i-rebinFactorX+1), h->GetBinLowEdge(i+1)-p};
1002 gSigma->SetPoint(i/rebinFactorX-1, p, fPGaus->GetParameter(2));
1003 gSigma->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fPGaus->GetParError(2), fPGaus->GetParError(2));
1004 delete tmp;
1005 }
1006 cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
1007}
1008
1009//------------------------------------------------------------------------------------
1010TCanvas* DrawVsAng(const char* name, const char* title, TH1* h1, TH2* h2)
1011{
1012 /// generic function to draw histograms versus absorber angular region
1013 TCanvas* c = new TCanvas(name, title);
1014 c->cd();
1015 h1->Draw();
1016 TH1D *proj1 = h2->ProjectionY(Form("%s_proj_0_2",h2->GetName()),1,2);
1017 proj1->Draw("sames");
1018 proj1->SetLineColor(2);
1019 TH1D *proj2 = h2->ProjectionY(Form("%s_proj_2_3",h2->GetName()),3,3);
1020 proj2->Draw("sames");
1021 proj2->SetLineColor(4);
1022 TH1D *proj3 = h2->ProjectionY(Form("%s__proj_3_10",h2->GetName()),4,10);
1023 proj3->Draw("sames");
1024 proj3->SetLineColor(3);
1025 return c;
1026}
1027
1028//------------------------------------------------------------------------------------
1029TCanvas* DrawVsPos(const char* name, const char* title, TH2* h1, TH2* h2, TH2* h3)
1030{
1031 /// generic function to draw histograms versus position at absorber end
1032 TCanvas* c = new TCanvas(name, title);
1033 c->cd();
1034 h1->Draw();
1035 h1->SetMarkerColor(2);
1036 h2->Draw("sames");
1037 h2->SetMarkerColor(4);
1038 h3->Draw("sames");
1039 h3->SetMarkerColor(3);
1040 return c;
1041}
1042
1043//------------------------------------------------------------------------------------
1044TCanvas* DrawResMomVsMom(const char* name, const char* title, TH2* h, Int_t nBins, TF1* f2, const char* fitting)
1045{
1046 /// generic function to draw and eventually fit momentum residuals versus momentum
1047 TLegend* l = new TLegend(0.15,0.25,0.3,0.85);
1048 TCanvas* c = new TCanvas(name, title);
1049 c->cd();
1050 TH1D* proj = 0x0;
1051 h->Sumw2();
1052 Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
1053 for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
1054 if (f2) cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
1055 proj = h->ProjectionY(Form("%s_%d",h->GetName(),i/rebinFactorX),i-rebinFactorX+1,i);
1056 if (proj->GetEntries() > 0) proj->Scale(1./proj->GetEntries());
1057 proj->Draw((i==rebinFactorX)?"hist":"histsames");
1058 proj->SetLineColor(i/rebinFactorX);
1059 if (f2) {
1060 f2->SetParameters(0.2,0.,1.,1.);
1061 f2->SetLineColor(i/rebinFactorX);
1062 proj->Fit("f2","WWNQ","sames");
1063 Double_t fwhm = f2->GetParameter(0);
1064 Double_t sigma = f2->GetParameter(3);
1065 Double_t sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
1066 Int_t rebin = TMath::Max(Int_t(0.5*sigmaP/proj->GetBinWidth(1)),1);
1067 while (proj->GetNbinsX()%rebin!=0) rebin--;
1068 proj->Rebin(rebin);
1069 proj->Scale(1./rebin);
1070 proj->Fit("f2","Q","sames");
1071 } else proj->SetLineWidth(2);
1072 Double_t p = 0.5 * (h->GetBinLowEdge(i-rebinFactorX+1) + h->GetBinLowEdge(i+1));
1073 l->AddEntry(proj,Form("%5.1f GeV",p));
1074 }
1075 if (f2) cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
1076 l->Draw("same");
1077 return c;
1078}
1079