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 ////////////////////////////////////////////////////////////////////////////
28 #include "TDecompLU.h"
29 #include "TGraphAsymmErrors.h"
33 #include "AliTRDgtuParam.h"
34 #include "AliTRDgeometry.h"
35 #include "AliTRDpadPlane.h"
37 ClassImp(AliTRDgtuParam)
39 AliTRDgtuParam *AliTRDgtuParam::fgInstance = 0;
41 // ----- Bin widths (granularity) -----
42 const Float_t AliTRDgtuParam::fgkBinWidthY = 160e-4;
43 const Float_t AliTRDgtuParam::fgkBinWidthdY = 140e-4;
45 // ----- Bit widths (used for internal representation) -----
46 const Int_t AliTRDgtuParam::fgkBitWidthY = 13;
47 const Int_t AliTRDgtuParam::fgkBitWidthdY = 7;
48 const Int_t AliTRDgtuParam::fgkBitWidthYProj = 10;
49 const Int_t AliTRDgtuParam::fgkBitExcessY = 4;
50 const Int_t AliTRDgtuParam::fgkBitExcessAlpha = 10;
51 const Int_t AliTRDgtuParam::fgkBitExcessYProj = 2;
53 // ----- Tracking parameters -----
55 const Int_t AliTRDgtuParam::fgkNZChannels = 3; // No. of z-channels
56 const Int_t AliTRDgtuParam::fgkNLinks = 12; // No. of links
57 const Int_t AliTRDgtuParam::fgkFixLayer = 2; // which layer is fixed for the generation of the z-channel map
58 const Int_t AliTRDgtuParam::fgkDeltaY = 39; // accepted deviation in y_proj, default: 9
59 const Int_t AliTRDgtuParam::fgkDeltaAlpha = 31; // accepted deviation in alpha, default: 11
60 const Int_t AliTRDgtuParam::fgkNRefLayers = 3; // no. of reference layers
63 AliTRDgtuParam::AliTRDgtuParam() :
71 fGeo = new AliTRDgeometry();
72 fRefLayers = new Int_t[fgkNRefLayers];
76 GenerateZChannelMap();
79 AliTRDgtuParam::~AliTRDgtuParam()
87 AliTRDgtuParam* AliTRDgtuParam::Instance()
89 // get (or create) the single instance
92 fgInstance = new AliTRDgtuParam();
97 void AliTRDgtuParam::Terminate()
99 // destruct the instance
101 if (fgInstance != 0) {
107 Bool_t AliTRDgtuParam::IsInZChannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const
109 return (fZSubChannel[stack][zchannel][layer][zpos] != 0);
112 Int_t AliTRDgtuParam::GetZSubchannel(Int_t stack, Int_t layer, Int_t zchannel, Int_t zpos) const
114 return fZSubChannel[stack][zchannel][layer][zpos];
117 Int_t AliTRDgtuParam::GetRefLayer(Int_t refLayerIdx) const
119 // returns the reference layer indexed by refLayerIdx
121 if (refLayerIdx >= 0 && refLayerIdx < fgkNRefLayers)
122 return fRefLayers[refLayerIdx];
127 Int_t AliTRDgtuParam::GenerateZChannelMap()
129 // generate the z-channel map
130 // assuming that the tracks come from the vertex
131 // +/- fVertexSize in z-direction
133 Int_t iSec = 0; // sector is irrelevant
134 Bool_t collision = kFALSE;
136 for (Int_t iStack = 0; iStack < fGeo->Nstack(); iStack++) {
138 Float_t x[6] = { 0 };
139 Float_t z[6][16] = {{ 0 }};
140 Float_t dZ[6][16] = {{ 0 }};
142 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
143 AliTRDpadPlane *pp = fGeo->GetPadPlane(iLayer, iStack);
144 x[iLayer] = fGeo->GetTime0(iLayer) - fGeo->CdrHght(); // ???
145 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, iSec); iRow++) {
146 z[iLayer][iRow] = pp->GetRowPos(iRow); // this is the right (pos. z-direction) border of the pad
147 dZ[iLayer][iRow] = pp->GetRowSize(iRow); // length of the pad in z-direction
148 for (Int_t i = 0; i < fgkNZChannels; i++)
149 fZSubChannel[iStack][i][iLayer][iRow] = 0;
153 for (Int_t fixRow = 0; fixRow < fGeo->GetRowMax(fgkFixLayer, iStack, iSec); fixRow++) {
155 Double_t fixZmin = z[fgkFixLayer][fixRow] - dZ[fgkFixLayer][fixRow];
156 Double_t fixZmax = z[fgkFixLayer][fixRow];
157 Double_t fixX = x[fgkFixLayer] + 1.5; // ??? 1.5 from where?
159 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
160 Double_t leftZ, rightZ;
162 if (iLayer <= fgkFixLayer) {
163 leftZ = (fixZmin + fVertexSize) * (x[iLayer] + 1.5) / fixX - fVertexSize;
164 rightZ = (fixZmax - fVertexSize) * (x[iLayer] + 1.5) / fixX + fVertexSize;
167 leftZ = (fixZmin - fVertexSize) * (x[iLayer] + 1.5) / fixX + fVertexSize;
168 rightZ = (fixZmax + fVertexSize) * (x[iLayer] + 1.5) / fixX - fVertexSize;
171 Double_t epsilon = 0.001;
172 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, iSec); iRow++) {
173 if ( (z[iLayer][iRow] ) > (leftZ + epsilon) &&
174 (z[iLayer][iRow] - dZ[iLayer][iRow] ) < (rightZ - epsilon) ) {
175 fZChannelMap[iStack][fixRow][iLayer][iRow] = 1;
176 if (fZSubChannel[iStack][fixRow % fgkNZChannels][iLayer][iRow] != 0) {
177 AliError("Collision in Z-Channel assignment occured! No reliable tracking!!!");
181 fZSubChannel[iStack][fixRow % fgkNZChannels][iLayer][iRow] = fixRow / fgkNZChannels + 1;
192 Bool_t AliTRDgtuParam::DisplayZChannelMap(Int_t zchannel, Int_t subchannel) const
194 // display the z-channel map
196 if (zchannel > fgkNZChannels) {
197 AliError("Invalid Z channel!");
201 Int_t zchmin = zchannel >= 0 ? zchannel : 0;
202 Int_t zchmax = zchannel >= 0 ? zchannel + 1 : fgkNZChannels;
205 TCanvas *c = new TCanvas("zchmap", "Z-Chhannel Mapping");
207 TGraph **graphz = new TGraph*[fgkNZChannels];
208 for (Int_t zch = zchmin; zch < zchmax; zch++)
209 graphz[zch] = new TGraph;
210 TGraphAsymmErrors *graph = new TGraphAsymmErrors();
211 graph->SetTitle("Z-Channel Map");
212 graph->SetPoint(i, 0, 0); // vertex
213 graph->SetPointError(i++, 20, 20, 0, 0);
214 // graph->SetRange //????
215 for (Int_t iLayer = 0; iLayer < fGeo->Nlayer(); iLayer++) {
216 for (Int_t iStack = 0; iStack < fGeo->Nstack(); iStack++) {
217 AliTRDpadPlane *pp = fGeo->GetPadPlane(iLayer, iStack);
218 for (Int_t iRow = 0; iRow < fGeo->GetRowMax(iLayer, iStack, 0); iRow++) {
219 graph->SetPoint(i, pp->GetRowPos(iRow), fGeo->GetTime0(iLayer) - fGeo->CdrHght());
220 graph->SetPointError(i++, pp->GetRowSize(iRow), 0, 0, 0);
221 for (Int_t zch = zchmin; zch < zchmax; zch++)
222 if (fZSubChannel[iStack][zch][iLayer][iRow] != 0)
223 if (subchannel == 0 || fZSubChannel[iStack][zch][iLayer][iRow] == subchannel)
224 graphz[zch]->SetPoint(j++, pp->GetRowPos(iRow) - pp->GetRowSize(iRow)/2, fGeo->GetTime0(iLayer) - fGeo->CdrHght());
228 graph->SetMarkerStyle(kDot);
230 for (Int_t zch = zchmin; zch < zchmax; zch++) {
231 graphz[zch]->SetMarkerStyle(kCircle);
232 graphz[zch]->SetMarkerColor(zch+2);
233 graphz[zch]->SetMarkerSize(0.3 + zch*0.2);
234 graphz[zch]->Draw("P");
239 Int_t AliTRDgtuParam::GetCiAlpha(Int_t layer) const
241 // get the constant for the calculation of alpha
243 Int_t ci = (Int_t) (GetChamberThickness() / fGeo->GetTime0(layer) * GetBinWidthY() / GetBinWidthdY() * (1 << (GetBitExcessAlpha() + GetBitExcessY() + 1)) );
247 Int_t AliTRDgtuParam::GetCiYProj(Int_t layer) const
249 // get the constant for the calculation of y_proj
251 Float_t xmid = (fGeo->GetTime0(0) + fGeo->GetTime0(fGeo->Nlayer()-1)) / 2.;
252 Int_t ci = (Int_t) (- (fGeo->GetTime0(layer) - xmid) / GetChamberThickness() * GetBinWidthdY() / GetBinWidthY() * (1 << GetBitExcessYProj()) );
256 Int_t AliTRDgtuParam::GetYt(Int_t stack, Int_t layer, Int_t zrow) const
258 return (Int_t) (- ( (layer % 2 ? 1 : -1) *
259 (GetGeo()->GetPadPlane(layer, stack)->GetRowPos(zrow) - GetGeo()->GetPadPlane(layer, stack)->GetRowSize(zrow) / 2) *
260 TMath::Tan(- 2.0 / 180.0 * TMath::Pi()) ) / 0.016 );
263 Bool_t AliTRDgtuParam::GenerateRecoCoefficients(Int_t trackletMask)
265 // calculate the coefficients for the straight line fit
266 // depending on the mask of contributing tracklets
268 fCurrTrackletMask = trackletMask;
270 TMatrix a(GetNLayers(), 3);
271 TMatrix b(3, GetNLayers());
274 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
275 if ( (trackletMask & (1 << layer)) == 0) {
282 a(layer, 1) = fGeo->GetTime0(layer);
283 a(layer, 2) = (layer % 2 ? 1 : -1) * fGeo->GetTime0(layer);
292 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
293 fAki[layer] = b.GetMatrixArray()[layer];
294 fBki[layer] = b.GetMatrixArray()[GetNLayers() + layer];
295 fCki[layer] = b.GetMatrixArray()[2 * GetNLayers() + layer];
300 Float_t AliTRDgtuParam::GetAki(Int_t k, Int_t i)
302 // get A_ki for the calculation of the tracking parameters
303 if (fCurrTrackletMask != k)
304 GenerateRecoCoefficients(k);
309 Float_t AliTRDgtuParam::GetBki(Int_t k, Int_t i)
311 // get B_ki for the calculation of the tracking parameters
313 if (fCurrTrackletMask != k)
314 GenerateRecoCoefficients(k);
319 Float_t AliTRDgtuParam::GetCki(Int_t k, Int_t i)
321 // get B_ki for the calculation of the tracking parameters
323 if (fCurrTrackletMask != k)
324 GenerateRecoCoefficients(k);
330 Float_t AliTRDgtuParam::GetD(Int_t k) const
332 // get the determinant for the calculation of the tracking parameters
335 for (Int_t i = 0; i < GetNLayers(); i++) {
336 if ( !((k >> i) & 0x1) )
338 Float_t xi = fGeo->GetTime0(i);
341 t(2,0) += TMath::Power(-1, i) * xi;
343 t(1,1) += TMath::Power(xi, 2);
344 t(2,1) += TMath::Power(-1, i) * TMath::Power(xi, 2);
345 t(0,2) += TMath::Power(-1, i) * xi;
346 t(1,2) += TMath::Power(-1, i) * TMath::Power(xi, 2);
347 t(2,2) += TMath::Power(xi, 2);
349 return t.Determinant();
352 Bool_t AliTRDgtuParam::GetFitParams(TVectorD& rhs, Int_t k)
354 // calculate the fitting parameters
358 for (Int_t i = 0; i < GetNLayers(); i++) {
359 if ( !((k >> i) & 0x1) )
361 Float_t xi = fGeo->GetTime0(i);
364 t(2,0) += TMath::Power(-1, i) * xi;
366 t(1,1) += TMath::Power(xi, 2);
367 t(2,1) += TMath::Power(-1, i) * TMath::Power(xi, 2);
368 t(0,2) -= TMath::Power(-1, i) * xi;
369 t(1,2) -= TMath::Power(-1, i) * TMath::Power(xi, 2);
370 t(2,2) -= TMath::Power(xi, 2);
374 return lr.Decompose();
378 Bool_t AliTRDgtuParam::GetIntersectionPoints(Int_t k, Float_t &x1, Float_t &x2)
380 // get the x-coord. of the assumed circle/straight line intersection points
385 for (Int_t layer = 0; layer < GetNLayers(); layer++) {
386 if ( (k >> layer) & 0x1 ) {
394 x1 = fGeo->GetTime0(l1) + 10./6 * (nHits -1);
395 x2 = fGeo->GetTime0(l2) - 10./6 * (nHits -1);
397 return ( (l1 >= 0) && (l2 >= 0) );
400 Float_t AliTRDgtuParam::GetPt(Int_t a, Float_t /* b */, Float_t x1, Float_t x2) const
402 // returns 0.3 * B * 1/a (1/128 GeV/c)
403 // a : offset, b : slope (not used)
405 Float_t c1 = x1 * x2 / 2. / 10000.; // conversion cm to m
408 r = (0.3 * fMagField / 2. / (fgkBinWidthY/100.)) * (((Int_t) c1) << 8) / (a >> 1); //??? why shift of a?