1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDgtuParam.cxx 28397 2008-09-02 09:33:00Z cblume $ */
18 ////////////////////////////////////////////////////////////////////////////
20 // Parameters for GTU simulation //
22 // Author: J. Klein (Jochen.Klein@cern.ch) //
24 ////////////////////////////////////////////////////////////////////////////
29 #include "TDecompLU.h"
30 #include "TGraphAsymmErrors.h"
34 #include "AliTRDgtuParam.h"
35 #include "AliTRDgeometry.h"
36 #include "AliTRDpadPlane.h"
38 ClassImp(AliTRDgtuParam)
40 AliTRDgtuParam *AliTRDgtuParam::fgInstance = 0;
42 // ----- Bin widths (granularity) -----
43 const Float_t AliTRDgtuParam::fgkBinWidthY = 160e-4;
44 const Float_t AliTRDgtuParam::fgkBinWidthdY = 140e-4;
46 // ----- Bit widths (used for internal representation) -----
47 const Int_t AliTRDgtuParam::fgkBitWidthY = 13;
48 const Int_t AliTRDgtuParam::fgkBitWidthdY = 7;
49 const Int_t AliTRDgtuParam::fgkBitWidthYProj = 10;
50 const Int_t AliTRDgtuParam::fgkBitExcessY = 4;
51 const Int_t AliTRDgtuParam::fgkBitExcessAlpha = 10;
52 const Int_t AliTRDgtuParam::fgkBitExcessYProj = 2;
54 // ----- Tracking parameters -----
56 const Int_t AliTRDgtuParam::fgkNZChannels = 3; // No. of z-channels
57 const Int_t AliTRDgtuParam::fgkNLinks = 12; // No. of links
58 const Int_t AliTRDgtuParam::fgkFixLayer = 2; // which layer is fixed for the generation of the z-channel map
59 const Int_t AliTRDgtuParam::fgkDeltaY = 39; // accepted deviation in y_proj, default: 9
60 const Int_t AliTRDgtuParam::fgkDeltaAlpha = 31; // accepted deviation in alpha, default: 11
61 const Int_t AliTRDgtuParam::fgkNRefLayers = 3; // no. of reference layers
64 AliTRDgtuParam::AliTRDgtuParam() :
72 fGeo = new AliTRDgeometry();
73 fRefLayers = new Int_t[fgkNRefLayers];
77 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
83 GenerateZChannelMap();
86 AliTRDgtuParam::~AliTRDgtuParam()
94 AliTRDgtuParam* AliTRDgtuParam::Instance()
96 // get (or create) the single instance
99 fgInstance = new AliTRDgtuParam();
104 void AliTRDgtuParam::Terminate()
106 // destruct the instance
108 if (fgInstance != 0) {
114 Bool_t AliTRDgtuParam::IsInZChannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const
116 return (fZSubChannel[stack][zchannel][layer][zpos] != 0);
119 Int_t AliTRDgtuParam::GetZSubchannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const
121 return fZSubChannel[stack][zchannel][layer][zpos];
124 Int_t AliTRDgtuParam::GetRefLayer(Int_t refLayerIdx) const
126 // returns the reference layer indexed by refLayerIdx
128 if (refLayerIdx >= 0 && refLayerIdx < fgkNRefLayers)
129 return fRefLayers[refLayerIdx];
134 Int_t AliTRDgtuParam::GenerateZChannelMap()
136 // generate the z-channel map
137 // assuming that the tracks come from the vertex
138 // +/- fVertexSize in z-direction
140 Int_t iSec = 0; // sector is irrelevant
141 Bool_t collision = kFALSE;
143 for (Int_t iStack = 0; iStack < fGeo->Nstack(); iStack++) {
145 Float_t x[6] = { 0 };
146 Float_t z[6][16] = {{ 0 }};
147 Float_t dZ[6][16] = {{ 0 }};
149 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
150 AliTRDpadPlane *pp = fGeo->GetPadPlane(iLayer, iStack);
151 x[iLayer] = fGeo->GetTime0(iLayer) - fGeo->CdrHght(); // ???
152 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, iSec); iRow++) {
153 z[iLayer][iRow] = pp->GetRowPos(iRow); // this is the right (pos. z-direction) border of the pad
154 dZ[iLayer][iRow] = pp->GetRowSize(iRow); // length of the pad in z-direction
155 for (Int_t i = 0; i < fgkNZChannels; i++)
156 fZSubChannel[iStack][i][iLayer][iRow] = 0;
160 for (Int_t fixRow = 0; fixRow < fGeo->GetRowMax(fgkFixLayer, iStack, iSec); fixRow++) {
162 Double_t fixZmin = z[fgkFixLayer][fixRow] - dZ[fgkFixLayer][fixRow];
163 Double_t fixZmax = z[fgkFixLayer][fixRow];
164 Double_t fixX = x[fgkFixLayer] + 1.5; // ??? 1.5 from where?
166 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
167 Double_t leftZ, rightZ;
169 if (iLayer <= fgkFixLayer) {
170 leftZ = (fixZmin + fVertexSize) * (x[iLayer] + 1.5) / fixX - fVertexSize;
171 rightZ = (fixZmax - fVertexSize) * (x[iLayer] + 1.5) / fixX + fVertexSize;
174 leftZ = (fixZmin - fVertexSize) * (x[iLayer] + 1.5) / fixX + fVertexSize;
175 rightZ = (fixZmax + fVertexSize) * (x[iLayer] + 1.5) / fixX - fVertexSize;
178 Double_t epsilon = 0.001;
179 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, iSec); iRow++) {
180 if ( (z[iLayer][iRow] ) > (leftZ + epsilon) &&
181 (z[iLayer][iRow] - dZ[iLayer][iRow] ) < (rightZ - epsilon) ) {
182 fZChannelMap[iStack][fixRow][iLayer][iRow] = 1;
183 if (fZSubChannel[iStack][fixRow % fgkNZChannels][iLayer][iRow] != 0) {
184 AliError("Collision in Z-Channel assignment occured! No reliable tracking!!!");
188 fZSubChannel[iStack][fixRow % fgkNZChannels][iLayer][iRow] = fixRow / fgkNZChannels + 1;
199 Bool_t AliTRDgtuParam::DisplayZChannelMap(Int_t zchannel, Int_t subchannel) const
201 // display the z-channel map
203 if (zchannel >= fgkNZChannels) {
204 AliError("Invalid Z channel!");
208 Int_t zchmin = zchannel >= 0 ? zchannel : 0;
209 Int_t zchmax = zchannel >= 0 ? zchannel + 1 : fgkNZChannels;
212 TCanvas *c = new TCanvas("zchmap", "Z-Chhannel Mapping");
214 TGraph **graphz = new TGraph*[fgkNZChannels];
215 for (Int_t zch = zchmin; zch < zchmax; zch++)
216 graphz[zch] = new TGraph;
217 TGraphAsymmErrors *graph = new TGraphAsymmErrors();
218 graph->SetTitle("Z-Channel Map");
219 graph->SetPoint(i, 0, 0); // vertex
220 graph->SetPointError(i++, 20, 20, 0, 0);
221 // graph->SetRange //????
222 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
223 for (Int_t iStack = 0; iStack < fGeo->Nstack(); iStack++) {
224 AliTRDpadPlane *pp = fGeo->GetPadPlane(iLayer, iStack);
225 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, 0); iRow++) {
226 graph->SetPoint(i, pp->GetRowPos(iRow), fGeo->GetTime0(iLayer) - fGeo->CdrHght());
227 graph->SetPointError(i++, pp->GetRowSize(iRow), 0, 0, 0);
228 for (Int_t zch = zchmin; zch < zchmax; zch++)
229 if (fZSubChannel[iStack][zch][iLayer][iRow] != 0)
230 if (subchannel == 0 || fZSubChannel[iStack][zch][iLayer][iRow] == subchannel)
231 graphz[zch]->SetPoint(j++, pp->GetRowPos(iRow) - pp->GetRowSize(iRow)/2, fGeo->GetTime0(iLayer) - fGeo->CdrHght());
235 graph->SetMarkerStyle(kDot);
238 for (Int_t zch = zchmin; zch < zchmax; zch++) {
239 graphz[zch]->SetMarkerStyle(kCircle);
240 graphz[zch]->SetMarkerColor(zch+2);
241 graphz[zch]->SetMarkerSize(0.3 + zch*0.2);
242 graphz[zch]->Draw("P");
243 gROOT->Add(graphz[zch]);
248 Int_t AliTRDgtuParam::GetCiAlpha(Int_t layer) const
250 // get the constant for the calculation of alpha
252 Int_t ci = (Int_t) (GetChamberThickness() / fGeo->GetTime0(layer) * GetBinWidthY() / GetBinWidthdY() * (1 << (GetBitExcessAlpha() + GetBitExcessY() + 1)) );
256 Int_t AliTRDgtuParam::GetCiYProj(Int_t layer) const
258 // get the constant for the calculation of y_proj
260 Float_t xmid = (fGeo->GetTime0(0) + fGeo->GetTime0(fGeo->Nlayer()-1)) / 2.;
261 Int_t ci = (Int_t) (- (fGeo->GetTime0(layer) - xmid) / GetChamberThickness() * GetBinWidthdY() / GetBinWidthY() * (1 << GetBitExcessYProj()) );
265 Int_t AliTRDgtuParam::GetYt(Int_t stack, Int_t layer, Int_t zrow) const
267 return (Int_t) (- ( (layer % 2 ? 1 : -1) *
268 (GetGeo()->GetPadPlane(layer, stack)->GetRowPos(zrow) - GetGeo()->GetPadPlane(layer, stack)->GetRowSize(zrow) / 2) *
269 TMath::Tan(- 2.0 / 180.0 * TMath::Pi()) ) / 0.016 );
272 Bool_t AliTRDgtuParam::GenerateRecoCoefficients(Int_t trackletMask)
274 // calculate the coefficients for the straight line fit
275 // depending on the mask of contributing tracklets
277 fCurrTrackletMask = trackletMask;
279 TMatrix a(GetNLayers(), 3);
280 TMatrix b(3, GetNLayers());
283 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
284 if ( (trackletMask & (1 << layer)) == 0) {
291 a(layer, 1) = fGeo->GetTime0(layer);
292 a(layer, 2) = (layer % 2 ? 1 : -1) * fGeo->GetTime0(layer);
301 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
302 fAki[layer] = b.GetMatrixArray()[layer];
303 fBki[layer] = b.GetMatrixArray()[GetNLayers() + layer];
304 fCki[layer] = b.GetMatrixArray()[2 * GetNLayers() + layer];
309 Float_t AliTRDgtuParam::GetAki(Int_t k, Int_t i)
311 // get A_ki for the calculation of the tracking parameters
312 if (fCurrTrackletMask != k)
313 GenerateRecoCoefficients(k);
318 Float_t AliTRDgtuParam::GetBki(Int_t k, Int_t i)
320 // get B_ki for the calculation of the tracking parameters
322 if (fCurrTrackletMask != k)
323 GenerateRecoCoefficients(k);
328 Float_t AliTRDgtuParam::GetCki(Int_t k, Int_t i)
330 // get B_ki for the calculation of the tracking parameters
332 if (fCurrTrackletMask != k)
333 GenerateRecoCoefficients(k);
339 Float_t AliTRDgtuParam::GetD(Int_t k) const
341 // get the determinant for the calculation of the tracking parameters
344 for (Int_t i = 0; i < GetNLayers(); i++) {
345 if ( !((k >> i) & 0x1) )
347 Float_t xi = fGeo->GetTime0(i);
350 t(2,0) += TMath::Power(-1, i) * xi;
352 t(1,1) += TMath::Power(xi, 2);
353 t(2,1) += TMath::Power(-1, i) * TMath::Power(xi, 2);
354 t(0,2) += TMath::Power(-1, i) * xi;
355 t(1,2) += TMath::Power(-1, i) * TMath::Power(xi, 2);
356 t(2,2) += TMath::Power(xi, 2);
358 return t.Determinant();
361 Bool_t AliTRDgtuParam::GetFitParams(TVectorD& rhs, Int_t k)
363 // calculate the fitting parameters
367 for (Int_t i = 0; i < GetNLayers(); i++) {
368 if ( !((k >> i) & 0x1) )
370 Float_t xi = fGeo->GetTime0(i);
373 t(2,0) += TMath::Power(-1, i) * xi;
375 t(1,1) += TMath::Power(xi, 2);
376 t(2,1) += TMath::Power(-1, i) * TMath::Power(xi, 2);
377 t(0,2) -= TMath::Power(-1, i) * xi;
378 t(1,2) -= TMath::Power(-1, i) * TMath::Power(xi, 2);
379 t(2,2) -= TMath::Power(xi, 2);
383 return lr.Decompose();
387 Bool_t AliTRDgtuParam::GetIntersectionPoints(Int_t k, Float_t &x1, Float_t &x2)
389 // get the x-coord. of the assumed circle/straight line intersection points
394 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
395 if ( (k >> layer) & 0x1 ) {
403 if ( (l1 >= 0) && (l2 >= 0) ) {
404 x1 = fGeo->GetTime0(l1) + 10./6 * (nHits -1);
405 x2 = fGeo->GetTime0(l2) - 10./6 * (nHits -1);
412 Float_t AliTRDgtuParam::GetPt(Int_t a, Float_t /* b */, Float_t x1, Float_t x2) const
414 // returns 0.3 * B * 1/a (1/128 GeV/c)
415 // a : offset, b : slope (not used)
417 Float_t c1 = x1 * x2 / 2. / 10000.; // conversion cm to m
420 r = (0.3 * fMagField / 2. / (fgkBinWidthY/100.)) * (((Int_t) c1) << 8) / (a >> 1); //??? why shift of a?