AliRoot Core  3dc7879 (3dc7879)
AliTPCParamSR.cxx
Go to the documentation of this file.
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 
27 
28 //#include <Riostream.h>
29 #include <TMath.h>
30 #include "TBuffer.h"
31 
32 #include "AliTPCPRF2D.h"
33 #include "AliTPCParamSR.h"
34 #include "AliTPCRF1D.h"
35 #include "TH1.h"
36 #include "AliTPCROC.h"
37 #include "TGeoManager.h"
38 
40 ClassImp(AliTPCParamSR)
42 static const Int_t kMaxRows=600;
43 static const Float_t kEdgeSectorSpace = 2.5;
44 static const Float_t kFacSigmaPadRow=3.;
45 static const Float_t kFacSigmaPad=3.;
46 static const Float_t kFacSigmaTime=3.;
47 
48 
50  :AliTPCParam(),
51  fInnerPRF(0),
52  fOuter1PRF(0),
53  fOuter2PRF(0),
54  fTimeRF(0),
55  fFacSigmaPadRow(0),
56  fFacSigmaPad(0),
57  fFacSigmaTime(0)
58 {
60 
61  fFacSigmaPadRow = Float_t(kFacSigmaPadRow);
62  fFacSigmaPad = Float_t(kFacSigmaPad);
63  fFacSigmaTime = Float_t(kFacSigmaTime);
64  SetDefault();
65  Update();
66 }
67 
69 {
71 
72  if (fInnerPRF != 0) delete fInnerPRF;
73  if (fOuter1PRF != 0) delete fOuter1PRF;
74  if (fOuter2PRF != 0) delete fOuter2PRF;
75  if (fTimeRF != 0) delete fTimeRF;
76 }
77 
79 {
81 
82  fbStatus = kFALSE;
84 }
85 
86 Int_t AliTPCParamSR::CalcResponse(Float_t* xyz, Int_t * index, Int_t row)
87 {
95 
96  if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
97  Error("AliTPCParamSR", "response function was not adjusted");
98  return -1;
99  }
100 
101  Float_t sfpadrow; // sigma of response function
102  Float_t sfpad; // sigma of
103  Float_t sftime= fFacSigmaTime*fTimeRF->GetSigma()/fZWidth; //3 sigma of time response
104  if (index[1]<fNInnerSector){
107  }
108  else{
109  if(row<fNRowUp1){
112  else{
115  }
116  }
117 
118  Int_t fpadrow = TMath::Max(TMath::Nint(index[2]+xyz[0]-sfpadrow),0); //"first" padrow
119  Int_t fpad = TMath::Nint(xyz[1]-sfpad); //first pad
120  Int_t ftime = TMath::Max(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()-sftime),0); // first time
121  Int_t lpadrow = TMath::Min(TMath::Nint(index[2]+xyz[0]+sfpadrow),fpadrow+19); //"last" padrow
122  lpadrow = TMath::Min(GetNRow(index[1])-1,lpadrow);
123  Int_t lpad = TMath::Min(TMath::Nint(xyz[1]+sfpad),fpad+19); //last pad
124  Int_t ltime = TMath::Min(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()+sftime),ftime+19); // last time
125  ltime = TMath::Min(ltime,GetMaxTBin()-1);
126  //
127  Int_t npads = GetNPads(index[1],row);
128  if (fpad<-npads/2)
129  fpad = -npads/2;
130  if (lpad>npads/2)
131  lpad= npads/2;
132  if (ftime<0) ftime=0;
133  //
134  if (row>=0) { //if we are interesting about given pad row
135  if (fpadrow<=row) fpadrow =row;
136  else
137  return 0;
138  if (lpadrow>=row) lpadrow = row;
139  else
140  return 0;
141  }
142 
143 
144  Float_t padres[20][20]; //I don't expect bigger number of bins
145  Float_t timeres[20];
146  Int_t cindex3=0;
147  Int_t cindex=0;
148  Float_t cweight = 0;
149  if (fpadrow>=0) {
150  //calculate padresponse function
151  Int_t padrow, pad;
152  for (padrow = fpadrow;padrow<=lpadrow;padrow++)
153  for (pad = fpad;pad<=lpad;pad++){
154  Float_t dy = (xyz[0]+Float_t(index[2]-padrow));
155  Float_t dx = (xyz[1]+Float_t(pad));
156  if (index[1]<fNInnerSector)
157  padres[padrow-fpadrow][pad-fpad]=fInnerPRF->GetPRF(dx*fInnerPadPitchWidth,dy*fInnerPadPitchLength);
158  else{
159  if(row<fNRowUp1){
160  padres[padrow-fpadrow][pad-fpad]=fOuter1PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter1PadPitchLength);}
161  else{
162  padres[padrow-fpadrow][pad-fpad]=fOuter2PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter2PadPitchLength);}}}
163  //calculate time response function
164  Int_t time;
165  for (time = ftime;time<=ltime;time++)
166  timeres[time-ftime]= fTimeRF->GetRF((-xyz[2]-xyz[3]+Float_t(time))*fZWidth);
167  //write over threshold values to stack
168  for (padrow = fpadrow;padrow<=lpadrow;padrow++)
169  for (pad = fpad;pad<=lpad;pad++)
170  for (time = ftime;time<=ltime;time++){
171  cweight = timeres[time-ftime]*padres[padrow-fpadrow][pad-fpad];
172  if (cweight>fResponseThreshold) {
173  fResponseBin[cindex3]=padrow;
174  fResponseBin[cindex3+1]=pad;
175  fResponseBin[cindex3+2]=time;
176  cindex3+=3;
177  fResponseWeight[cindex]=cweight;
178  cindex++;
179  }
180  }
181  }
182  fCurrentMax=cindex;
183  return fCurrentMax;
184 }
185 
186 void AliTPCParamSR::TransformTo8(Float_t *xyz, Int_t *index) const
187 {
189 
190  if (index[0]==0) Transform0to1(xyz,index);
191  if (index[0]==1) Transform1to2(xyz,index);
192  if (index[0]==2) Transform2to3(xyz,index);
193  if (index[0]==3) Transform3to4(xyz,index);
194  if (index[0]==4) Transform4to8(xyz,index);
195 }
196 
197 void AliTPCParamSR::TransformTo2(Float_t *xyz, Int_t *index) const
198 {
202 
203  if (index[0]==0) Transform0to1(xyz,index);
204  if (index[0]==1) Transform1to2(xyz,index);
205  if (index[0]==4) Transform4to3(xyz,index);
206  if (index[0]==8) { //if we are in digit coordinate system transform to global
207  Transform8to4(xyz,index);
208  Transform4to3(xyz,index);
209  }
210 }
211 
212 void AliTPCParamSR::CRXYZtoXYZ(Float_t *xyz,
213  const Int_t &sector, const Int_t & padrow, Int_t option) const
214 {
216 
217  Bool_t rel = ( (option&2)!=0);
218  Int_t index[3]={sector,padrow,0};
219  if (rel==kTRUE) Transform4to3(xyz,index);//if the position is relative to pad row
220  Transform2to1(xyz,index);
221 }
222 
223 void AliTPCParamSR::XYZtoCRXYZ(Float_t *xyz,
224  Int_t &sector, Int_t & padrow, Int_t option) const
225 {
232 
233  Int_t index[3];
234  Bool_t rel = ( (option&2)!=0);
235 
236  //option 0 and 2 means that we don't have information about sector
237  if ((option&1)==0) Transform0to1(xyz,index); //we calculate sector number
238  else
239  index[0]=sector;
240  Transform1to2(xyz,index);
241  Transform2to3(xyz,index);
242  //if we store relative position calculate position relative to pad row
243  if (rel==kTRUE) Transform3to4(xyz,index);
244  sector = index[0];
245  padrow = index[1];
246 }
247 
248 Float_t AliTPCParamSR::GetPrimaryLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
249 {
251 
252  Float_t padlength=GetPadPitchLength(index[1]);
253  Float_t a1=TMath::Sin(angle[0]);
254  a1*=a1;
255  Float_t a2=TMath::Sin(angle[1]);
256  a2*=a2;
257  Float_t length =padlength*TMath::Sqrt(1+a1+a2);
258  return length*fNPrimLoss;
259 }
260 
261 Float_t AliTPCParamSR::GetTotalLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
262 {
264 
265  Float_t padlength=GetPadPitchLength(index[1]);
266  Float_t a1=TMath::Sin(angle[0]);
267  a1*=a1;
268  Float_t a2=TMath::Sin(angle[1]);
269  a2*=a2;
270  Float_t length =padlength*TMath::Sqrt(1+a1+a2);
271  return length*fNTotalLoss;
272 
273 }
274 
275 
276 void AliTPCParamSR::GetClusterSize(Float_t *x, Int_t *index, Float_t */*angle*/, Int_t /*mode*/, Float_t *sigma)
277 {
282 
283  Float_t xx;
284  Float_t lx[3] = {x[0],x[1],x[2]};
285  Int_t li[3] = {index[0],index[1],index[2]};
286  TransformTo2(lx,li);
287  // Float_t sigmadiff;
288  sigma[0]=0;
289  sigma[1]=0;
290 
291  xx = lx[2]; //calculate drift length in cm
292  if (xx>0) {
293  sigma[0]+= xx*GetDiffL()*GetDiffL();
294  sigma[1]+= xx*GetDiffT()*GetDiffT();
295  }
296 
297 
298  //sigma[0]=sigma[1]=0;
299  if (GetTimeRF()!=0) sigma[0]+=GetTimeRF()->GetSigma()*GetTimeRF()->GetSigma();
300  if ( (index[1]<fNInnerSector) &&(GetInnerPRF()!=0))
301  sigma[1]+=GetInnerPRF()->GetSigmaX()*GetInnerPRF()->GetSigmaX();
302  if ( (index[1]>=fNInnerSector) &&(index[2]<fNRowUp1) && (GetOuter1PRF()!=0))
303  sigma[1]+=GetOuter1PRF()->GetSigmaX()*GetOuter1PRF()->GetSigmaX();
304  if( (index[1]>=fNInnerSector) &&(index[2]>=fNRowUp1) && (GetOuter2PRF()!=0))
305  sigma[1]+=GetOuter2PRF()->GetSigmaX()*GetOuter2PRF()->GetSigmaX();
306 
307 
308  sigma[0]/= GetZWidth()*GetZWidth();
309  sigma[1]/=GetPadPitchWidth(index[0])*GetPadPitchWidth(index[0]);
310 }
311 
312 
313 
314 
315 void AliTPCParamSR::GetSpaceResolution(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/,
316  Float_t /*amplitude*/, Int_t /*mode*/, Float_t */*sigma*/)
317 {
319 
320 }
321 Float_t AliTPCParamSR::GetAmp(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/)
322 {
324 
325  return 0;
326 }
327 
328 Float_t * AliTPCParamSR::GetAnglesAccMomentum(Float_t *x, Int_t * index, Float_t* momentum, Float_t *angle)
329 {
332 
333  TransformTo2(x,index);
334  AliDetectorParam::GetAnglesAccMomentum(x,index,momentum,angle);
335  Float_t addangle = TMath::ASin(x[1]/GetPadRowRadii(index[1],index[2]));
336  angle[1] +=addangle;
337  return angle;
338 }
339 
340 
342 {
343  Int_t i;
344  if (AliTPCParam::Update()==kFALSE) return kFALSE;
345  fbStatus = kFALSE;
346 
347  Float_t firstrow = fInnerRadiusLow + 1.575;
348  for( i= 0;i<fNRowLow;i++)
349  {
350  Float_t x = firstrow + fInnerPadPitchLength*(Float_t)i;
351  fPadRowLow[i]=x;
352  // number of pads per row
353  // Float_t y = (x-0.5*fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount-
354  // fInnerPadPitchWidth/2.;
355  // 0 and fNRowLow+1 reserved for cross talk rows
356  fYInner[i+1] = x*tan(fInnerAngle/2.)-fInnerWireMount;
357  //fNPadsLow[i] = 1+2*(Int_t)(y/fInnerPadPitchWidth) ;
358  fNPadsLow[i] = AliTPCROC::Instance()->GetNPads(0,i) ; // ROC implement
359  }
360  // cross talk rows
362  fYInner[fNRowLow+1]=(fPadRowLow[fNRowLow-1]+fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
363  firstrow = fOuterRadiusLow + 1.6;
364  for(i=0;i<fNRowUp;i++)
365  {
366  if(i<fNRowUp1){
367  Float_t x = firstrow + fOuter1PadPitchLength*(Float_t)i;
368  fPadRowUp[i]=x;
369 // Float_t y =(x-0.5*fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
370 // fOuterPadPitchWidth/2.;
371  fYOuter[i+1]= x*tan(fOuterAngle/2.)-fOuterWireMount;
372  //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
373  fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
374  if(i==fNRowUp1-1) {
375  fLastWireUp1=fPadRowUp[i] +0.625;
376  firstrow = fPadRowUp[i] + 0.5*(fOuter1PadPitchLength+fOuter2PadPitchLength);
377  }
378  }
379  else
380  {
381  Float_t x = firstrow + fOuter2PadPitchLength*(Float_t)(i-64);
382  fPadRowUp[i]=x;
383  //Float_t y =(x-0.5*fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
384  // fOuterPadPitchWidth/2.;
385  //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
386  fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
387  }
388  fYOuter[i+1] = fPadRowUp[i]*tan(fOuterAngle/2.)-fOuterWireMount;
389  }
390  // cross talk rows
392  fYOuter[fNRowUp+1]=(fPadRowUp[fNRowUp-1]+fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
394  fbStatus = kTRUE;
395  return kTRUE;
396 }
397 Float_t AliTPCParamSR::GetYInner(Int_t irow) const
398 {
399  return fYInner[irow];
400 }
401 Float_t AliTPCParamSR::GetYOuter(Int_t irow) const
402 {
403  return fYOuter[irow];
404 }
405 
406 void AliTPCParamSR::Streamer(TBuffer &R__b)
407 {
409 
410  if (R__b.IsReading()) {
411  Version_t R__v = R__b.ReadVersion(); if (R__v) { }
412  // TObject::Streamer(R__b);
413  AliTPCParam::Streamer(R__b);
414  // if (R__v < 2) return;
415  Update();
417  } else {
418  R__b.WriteVersion(AliTPCParamSR::IsA());
419  //TObject::Streamer(R__b);
420  AliTPCParam::Streamer(R__b);
421  }
422 }
423 Int_t AliTPCParamSR::CalcResponseFast(Float_t* xyz, Int_t * index, Int_t row, Float_t phase)
424 {
433 
434  if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
435  Error("AliTPCParamSR", "response function was not adjusted");
436  return -1;
437  }
438 
439  const Int_t kpadn = 500;
440  const Float_t kfpadn = 500.;
441  const Int_t ktimen = 500;
442  const Float_t kftimen = 500.;
443  const Int_t kpadrn = 500;
444  const Float_t kfpadrn = 500.;
445 
446 
447 
448  static Float_t prfinner[2*kpadrn][5*kpadn]; //pad divided by 50
449  static Float_t prfouter1[2*kpadrn][5*kpadn]; //prfouter division
450  static Float_t prfouter2[2*kpadrn][5*kpadn];
451  static Float_t kTanMax =0;
452 
453  static Float_t rftime[5*ktimen]; //time division
454  static Int_t blabla=0;
455  static Float_t zoffset=0;
456  static Float_t zwidth=0;
457  static Float_t zoffset2=0;
458  static TH1F * hdiff=0;
459  static TH1F * hdiff1=0;
460  static TH1F * hdiff2=0;
461 
462  if (blabla==0) { //calculate Response function - only at the begginning
463  kTanMax = TMath::ATan(10.*TMath::DegToRad());
464  hdiff =new TH1F("prf_diff","prf_diff",10000,-1,1);
465  hdiff1 =new TH1F("no_repsonse1","no_response1",10000,-1,1);
466  hdiff2 =new TH1F("no_response2","no_response2",10000,-1,1);
467 
468  blabla=1;
469  zoffset = GetZOffset();
470  zwidth = fZWidth;
471  zoffset2 = zoffset/zwidth;
472  for (Int_t i=0;i<5*ktimen;i++){
473  rftime[i] = fTimeRF->GetRF(((i-2.5*kftimen)/kftimen)*zwidth+zoffset);
474  }
475  for (Int_t i=0;i<5*kpadn;i++){
476  for (Int_t j=0;j<2*kpadrn;j++){
477  prfinner[j][i] =
478  fInnerPRF->GetPRF((i-2.5*kfpadn)/kfpadn
479  *fInnerPadPitchWidth,(j-kfpadrn)/kfpadrn*fInnerPadPitchLength);
480  prfouter1[j][i] =
481  fOuter1PRF->GetPRF((i-2.5*kfpadn)/kfpadn
482  *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter1PadPitchLength);
483 
484  //
485  prfouter2[j][i] =
486  fOuter2PRF->GetPRF((i-2.5*kfpadn)/kfpadn
487  *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter2PadPitchLength);
488  }
489  }
490  } // the above is calculated only once
491 
492  // calculate central padrow, pad, time
493  Int_t npads = GetNPads(index[1],index[3]-1);
494  Int_t cpadrow = index[2]; // electrons are here
495  Int_t cpad = TMath::Nint(xyz[1]);
496  Int_t ctime = TMath::Nint(xyz[2]+zoffset2+xyz[3]);
497  //calulate deviation
498  Float_t dpadrow = xyz[0];
499  Float_t dpad = xyz[1]-cpad;
500  Float_t dtime = xyz[2]+zoffset2+xyz[3]-ctime+phase*0.25;
501  Int_t cindex =0;
502  Int_t cindex3 =0;
503  Int_t maxt =GetMaxTBin();
504 
505  Int_t fpadrow;
506  Int_t lpadrow;
507 
508  if (row>=0) { //if we are interesting about given pad row
509  fpadrow = row-cpadrow;
510  lpadrow = row-cpadrow;
511  }else{
512  fpadrow = (index[2]>1) ? -1 :0;
513  lpadrow = (index[2]<GetNRow(index[1])-1) ? 1:0;
514  }
515 
516  Int_t fpad = (cpad > -npads/2+1) ? -2: -npads/2-cpad;
517  Int_t lpad = (cpad < npads/2-2) ? 2: npads/2-1-cpad;
518  Int_t ftime = (ctime>1) ? -2: -ctime;
519  Int_t ltime = (ctime<maxt-2) ? 2: maxt-ctime-1;
520 
521  // cross talk from long pad to short one
522  if(row==fNRowUp1 && fpadrow==-1) {
523  dpadrow *= fOuter2PadPitchLength;
524  dpadrow += fOuterWWPitch;
525  dpadrow /= fOuter1PadPitchLength;
526  }
527  // cross talk from short pad to long one
528  if(row==fNRowUp1+1 && fpadrow==1){
529  dpadrow *= fOuter1PadPitchLength;
530  if(dpadrow < 0.) dpadrow = -1.; //protection against 3rd wire
531  dpadrow += fOuterWWPitch;
532  dpadrow /= fOuter2PadPitchLength;
533 
534  }
535 
536  // "normal"
537  Int_t apadrow = TMath::Nint((dpadrow-fpadrow)*kfpadrn+kfpadrn);
538  for (Int_t ipadrow = fpadrow; ipadrow<=lpadrow;ipadrow++){
539  if ( (apadrow<0) || (apadrow>=2*kpadrn))
540  continue;
541  // pad angular correction
542  Float_t angle = 0.;
543  if (npads != 0)
544  angle = kTanMax*2.*(cpad+0.5)/Float_t(npads);
545  Float_t dpadangle =0;
546  if (index[1]<fNInnerSector){
547  dpadangle = angle*dpadrow*fInnerPadPitchLength/fInnerPadPitchWidth;
548  }
549  else{
550  if(row < fNRowUp1+1){
551  dpadangle = angle*dpadrow*fOuter1PadPitchLength/fOuterPadPitchWidth;
552  }
553  else {
554  dpadangle = angle*dpadrow*fOuter2PadPitchLength/fOuterPadPitchWidth;
555  }
556  }
557  if (ipadrow==0) dpadangle *=-1;
558  //
559  // Int_t apad= TMath::Nint((dpad-fpad)*kfpadn+2.5*kfpadn);
560  Int_t apad= TMath::Nint((dpad+dpadangle-fpad)*kfpadn+2.5*kfpadn);
561  for (Int_t ipad = fpad; ipad<=lpad;ipad++){
562  Float_t cweight;
563  if (index[1]<fNInnerSector){
564  cweight=prfinner[apadrow][apad];
565  }
566  else{
567  if(row < fNRowUp1+1){
568  cweight=prfouter1[apadrow][apad];
569  }
570  else {
571  cweight=prfouter2[apadrow][apad];
572  }
573  }
574  // if (cweight<fResponseThreshold) continue;
575  Int_t atime = TMath::Nint((dtime-ftime)*kftimen+2.5*kftimen);
576  for (Int_t itime = ftime;itime<=ltime;itime++){
577  Float_t cweight2 = cweight*rftime[atime];
578  if (cweight2>fResponseThreshold) {
579  fResponseBin[cindex3++]=cpadrow+ipadrow;
580  fResponseBin[cindex3++]=cpad+ipad;
581  fResponseBin[cindex3++]=ctime+itime;
582  fResponseWeight[cindex++]=cweight2;
583  }
584  atime-=ktimen;
585  }
586  apad-= kpadn;
587  }
588  apadrow-=kpadrn;
589  }
590  fCurrentMax=cindex;
591  return fCurrentMax;
592 
593 }
594 
595 
596 
597 
598 
599 
600 
Float_t fInnerWireMount
space for wire mount, inner sector
Definition: AliTPCParam.h:417
Int_t * fResponseBin
! array with bins -calulated
Definition: AliTPCParam.h:542
Float_t GetPadPitchLength(Int_t isector=0, Int_t padrow=0) const
Definition: AliTPCParam.h:302
virtual Float_t GetAmp(Float_t *x, Int_t *index, Float_t *angle)
Float_t fLastWireUp1
position of the last wire in outer1 sector
Definition: AliTPCParam.h:449
Int_t Transform0to1(Float_t *xyz, Int_t *index) const
Float_t fInnerPadPitchWidth
Inner pad pitch width.
Definition: AliTPCParam.h:461
UInt_t GetNPads(UInt_t sector, UInt_t row) const
Definition: AliTPCROC.h:30
Float_t GetDiffL() const
Definition: AliTPCParam.h:338
Float_t fOuterWWPitch
pitch between wires in outer sector -calculated
Definition: AliTPCParam.h:452
void CRXYZtoXYZ(Float_t *xyz, const Int_t &sector, const Int_t &padrow, Int_t option=3) const
Float_t fPadRowLow[600]
Lower sector, pad row radii -calculated.
Definition: AliTPCParam.h:478
Manager and of geomety classes for set: TPC.
Definition: AliTPCParamSR.h:15
Bool_t fbStatus
indicates consistency of the data
Definition: AliTPCParam.h:403
AliTPCPRF2D * GetOuter2PRF() const
Definition: AliTPCParamSR.h:48
AliTPCPRF2D * fOuter1PRF
pad response function for outer sector
Definition: AliTPCParamSR.h:69
Manager and of geomety classes for set: TPC.
Definition: AliTPCParam.h:18
virtual Float_t * GetAnglesAccMomentum(Float_t *x, Int_t *index, Float_t *momentum, Float_t *angle)
AliTPCPRF2D * GetInnerPRF() const
Definition: AliTPCParamSR.h:46
void Transform4to3(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:714
virtual Bool_t ReadGeoMatrices()
Float_t fNPrimLoss
number of produced primary electrons per cm
AliTPCPRF2D * fInnerPRF
pad response function for inner sector
Definition: AliTPCParamSR.h:68
AliTPCPRF2D * GetOuter1PRF() const
Definition: AliTPCParamSR.h:47
Float_t GetPadPitchWidth(Int_t isector=0) const
Definition: AliTPCParam.h:300
virtual Float_t GetSigmaX() const
Definition: AliTPCPRF2D.h:67
AliTPCRF1D * fTimeRF
time response function object
Definition: AliTPCParamSR.h:71
Float_t fOuterRadiusLow
lower radius of outer sector-IP
Definition: AliTPCParam.h:410
Float_t fInnerPadPitchLength
Inner pad pitch length.
Definition: AliTPCParam.h:460
void TransformTo8(Float_t *xyz, Int_t *index) const
void Transform4to8(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:759
Int_t CalcResponse(Float_t *x, Int_t *index, Int_t row)
Int_t Transform2to3(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:682
Int_t fNInnerSector
number of inner sectors -calculated
Definition: AliTPCParam.h:419
virtual void GetClusterSize(Float_t *x, Int_t *index, Float_t *angle, Int_t mode, Float_t *sigma)
void XYZtoCRXYZ(Float_t *xyz, Int_t &sector, Int_t &padrow, Int_t option=3) const
static const Float_t kFacSigmaPad
Float_t fPadRowUp[600]
Upper sector, pad row radii -calculated.
Definition: AliTPCParam.h:479
virtual Float_t GetTotalLoss(Float_t *x, Int_t *index, Float_t *angle)
Float_t fOuterWireMount
space for wire mount, outer sector
Definition: AliTPCParam.h:418
Float_t fNTotalLoss
total number of produced electrons per cm
Float_t fYInner[600]
Inner sector, wire-length.
Definition: AliTPCParam.h:482
Int_t fNRowUp
number of pad rows per sector up -calculated
Definition: AliTPCParam.h:476
virtual Float_t GetZOffset() const
Definition: AliTPCParam.h:370
Float_t fOuterAngle
opening angle of outer sector
Definition: AliTPCParam.h:413
Float_t fFacSigmaPad
factor-how many sigma of response I accept
Definition: AliTPCParamSR.h:73
Float_t GetZWidth() const
Definition: AliTPCParam.h:367
Float_t fResponseThreshold
threshold for accepted response
Definition: AliTPCParam.h:540
Float_t GetYOuter(Int_t irow) const
Float_t * fResponseWeight
! array with response -calulated
Definition: AliTPCParam.h:543
Int_t fNPadsLow[600]
Lower sector, number of pads per row -calculated.
Definition: AliTPCParam.h:480
Int_t fNPadsUp[600]
Upper sector, number of pads per row -calculated.
Definition: AliTPCParam.h:481
TGeoManager * gGeoManager
virtual void SetDefault()
Float_t fInnerAngle
opening angle of Inner sector
Definition: AliTPCParam.h:411
Int_t CalcResponseFast(Float_t *x, Int_t *index, Int_t row, Float_t phase)
Float_t GetDiffT() const
Definition: AliTPCParam.h:337
Int_t fCurrentMax
! current maximal dimension -calulated
Definition: AliTPCParam.h:541
virtual Bool_t Update()
virtual Float_t GetPRF(Float_t xin, Float_t yin)
virtual void GetSpaceResolution(Float_t *x, Int_t *index, Float_t *angle, Float_t amplitude, Int_t mode, Float_t *sigma)
static const Int_t kMaxRows
virtual Float_t * GetAnglesAccMomentum(Float_t *x, Int_t *index, Float_t *momentum, Float_t *angle)
Int_t GetNRow(Int_t isec) const
Definition: AliTPCParam.h:309
void Transform8to4(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:784
Float_t fYOuter[600]
Outer sector, wire-length.
Definition: AliTPCParam.h:483
Float_t fZWidth
derived value calculated using TSample and driftw -computed
Definition: AliTPCParam.h:517
Float_t fOuterPadPitchWidth
Outer pad pitch width.
Definition: AliTPCParam.h:466
Int_t GetNPads(Int_t isector, Int_t irow) const
Definition: AliTPCParam.h:318
Float_t fOuter1PadPitchLength
Outer pad pitch length.
Definition: AliTPCParam.h:464
virtual ~AliTPCParamSR()
virtual Float_t GetPrimaryLoss(Float_t *x, Int_t *index, Float_t *angle)
Float_t GetYInner(Int_t irow) const
Int_t fNOuterSector
number of outer sectors -calculated
Definition: AliTPCParam.h:420
static const Float_t kFacSigmaTime
Float_t GetRF(Float_t xin)
Definition: AliTPCRF1D.cxx:174
void Transform2to1(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:636
virtual Float_t GetSigmaY() const
Definition: AliTPCPRF2D.h:68
Int_t fNRowLow
number of pad rows per low sector -set
Definition: AliTPCParam.h:473
static AliTPCROC * Instance()
Definition: AliTPCROC.cxx:34
AliTPCPRF2D * fOuter2PRF
Definition: AliTPCParamSR.h:70
Float_t fInnerRadiusLow
lower radius of inner sector-IP
Definition: AliTPCParam.h:407
Int_t fNRowUp1
number of short pad rows per sector up -set
Definition: AliTPCParam.h:474
void Transform3to4(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:699
Float_t GetPadRowRadii(Int_t isec, Int_t irow) const
Definition: AliTPCParam.h:313
Float_t fOuter2PadPitchLength
Outer pad pitch length.
Definition: AliTPCParam.h:465
Int_t GetMaxTBin() const
Definition: AliTPCParam.h:371
AliTPCRF1D * GetTimeRF() const
Definition: AliTPCParamSR.h:49
void Transform1to2(Float_t *xyz, Int_t *index) const
Definition: AliTPCParam.h:617
void TransformTo2(Float_t *xyz, Int_t *index) const
static const Float_t kEdgeSectorSpace
Float_t fFacSigmaTime
factor-how many sigma of response I accept
Definition: AliTPCParamSR.h:74
Float_t fFacSigmaPadRow
factor-how many sigma of response I accept
Definition: AliTPCParamSR.h:72
Float_t GetSigma()
Definition: AliTPCRF1D.h:42
Int_t fNtRows
total number of rows in TPC -calculated
Definition: AliTPCParam.h:477
static const Float_t kFacSigmaPadRow