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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////
19 // Manager and of geomety classes for set: TPC //
21 // !sectors are numbered from 0 //
22 // !pad rows are numbered from 0 //
24 // 27.7. - AliTPCPaaramSr object for TPC
25 // TPC with straight pad rows
26 // Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
28 ///////////////////////////////////////////////////////////////////////
30 //#include <Riostream.h>
33 #include "AliTPCPRF2D.h"
34 #include "AliTPCParamSR.h"
35 #include "AliTPCRF1D.h"
37 #include "AliTPCROC.h"
38 #include "TGeoManager.h"
40 ClassImp(AliTPCParamSR)
41 static const Int_t kMaxRows=600;
42 static const Float_t kEdgeSectorSpace = 2.5;
43 static const Float_t kFacSigmaPadRow=3.;
44 static const Float_t kFacSigmaPad=3.;
45 static const Float_t kFacSigmaTime=3.;
48 AliTPCParamSR::AliTPCParamSR()
51 //constructor set the default parameters
56 fFacSigmaPadRow = Float_t(kFacSigmaPadRow);
57 fFacSigmaPad = Float_t(kFacSigmaPad);
58 fFacSigmaTime = Float_t(kFacSigmaTime);
63 AliTPCParamSR::~AliTPCParamSR()
66 //destructor destroy some dynmicaly alocated variables
67 if (fInnerPRF != 0) delete fInnerPRF;
68 if (fOuter1PRF != 0) delete fOuter1PRF;
69 if (fOuter2PRF != 0) delete fOuter2PRF;
70 if (fTimeRF != 0) delete fTimeRF;
73 void AliTPCParamSR::SetDefault()
75 //set default TPC param
77 AliTPCParam::SetDefault();
80 Int_t AliTPCParamSR::CalcResponse(Float_t* xyz, Int_t * index, Int_t row)
83 //calculate bin response as function of the input position -x
84 //return number of valid response bin
86 //we suppose that coordinate is expressed in float digits
87 // it's mean coordinate system 8
88 //xyz[0] - float padrow xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
89 //xyz[3] - electron time in float time bin format
90 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
91 Error("AliTPCParamSR", "response function was not adjusted");
95 Float_t sfpadrow; // sigma of response function
96 Float_t sfpad; // sigma of
97 Float_t sftime= fFacSigmaTime*fTimeRF->GetSigma()/fZWidth; //3 sigma of time response
98 if (index[1]<fNInnerSector){
99 sfpadrow =fFacSigmaPadRow*fInnerPRF->GetSigmaY()/fInnerPadPitchLength;
100 sfpad =fFacSigmaPad*fInnerPRF->GetSigmaX()/fInnerPadPitchWidth;
104 sfpadrow =fFacSigmaPadRow*fOuter1PRF->GetSigmaY()/fOuter1PadPitchLength;
105 sfpad =fFacSigmaPad*fOuter1PRF->GetSigmaX()/fOuterPadPitchWidth;}
107 sfpadrow =fFacSigmaPadRow*fOuter2PRF->GetSigmaY()/fOuter2PadPitchLength;
108 sfpad =fFacSigmaPad*fOuter2PRF->GetSigmaX()/fOuterPadPitchWidth;
112 Int_t fpadrow = TMath::Max(TMath::Nint(index[2]+xyz[0]-sfpadrow),0); //"first" padrow
113 Int_t fpad = TMath::Nint(xyz[1]-sfpad); //first pad
114 Int_t ftime = TMath::Max(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()-sftime-GetNTBinsL1()),0); // first time
115 Int_t lpadrow = TMath::Min(TMath::Nint(index[2]+xyz[0]+sfpadrow),fpadrow+19); //"last" padrow
116 lpadrow = TMath::Min(GetNRow(index[1])-1,lpadrow);
117 Int_t lpad = TMath::Min(TMath::Nint(xyz[1]+sfpad),fpad+19); //last pad
118 Int_t ltime = TMath::Min(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()+sftime-GetNTBinsL1()),ftime+19); // last time
119 ltime = TMath::Min(ltime,GetMaxTBin()-1);
121 Int_t npads = GetNPads(index[1],row);
126 if (ftime<0) ftime=0;
128 if (row>=0) { //if we are interesting about given pad row
129 if (fpadrow<=row) fpadrow =row;
132 if (lpadrow>=row) lpadrow = row;
138 Float_t padres[20][20]; //I don't expect bigger number of bins
144 //calculate padresponse function
146 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
147 for (pad = fpad;pad<=lpad;pad++){
148 Float_t dy = (xyz[0]+Float_t(index[2]-padrow));
149 Float_t dx = (xyz[1]+Float_t(pad));
150 if (index[1]<fNInnerSector)
151 padres[padrow-fpadrow][pad-fpad]=fInnerPRF->GetPRF(dx*fInnerPadPitchWidth,dy*fInnerPadPitchLength);
154 padres[padrow-fpadrow][pad-fpad]=fOuter1PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter1PadPitchLength);}
156 padres[padrow-fpadrow][pad-fpad]=fOuter2PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter2PadPitchLength);}}}
157 //calculate time response function
159 for (time = ftime;time<=ltime;time++)
160 timeres[time-ftime]= fTimeRF->GetRF((-xyz[2]-xyz[3]+Float_t(time))*fZWidth+GetNTBinsL1());
161 //write over threshold values to stack
162 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
163 for (pad = fpad;pad<=lpad;pad++)
164 for (time = ftime;time<=ltime;time++){
165 cweight = timeres[time-ftime]*padres[padrow-fpadrow][pad-fpad];
166 if (cweight>fResponseThreshold) {
167 fResponseBin[cindex3]=padrow;
168 fResponseBin[cindex3+1]=pad;
169 fResponseBin[cindex3+2]=time;
171 fResponseWeight[cindex]=cweight;
180 void AliTPCParamSR::TransformTo8(Float_t *xyz, Int_t *index) const
183 // transformate point to digit coordinate
185 if (index[0]==0) Transform0to1(xyz,index);
186 if (index[0]==1) Transform1to2(xyz,index);
187 if (index[0]==2) Transform2to3(xyz,index);
188 if (index[0]==3) Transform3to4(xyz,index);
189 if (index[0]==4) Transform4to8(xyz,index);
192 void AliTPCParamSR::TransformTo2(Float_t *xyz, Int_t *index) const
195 //transformate point to rotated coordinate
198 if (index[0]==0) Transform0to1(xyz,index);
199 if (index[0]==1) Transform1to2(xyz,index);
200 if (index[0]==4) Transform4to3(xyz,index);
201 if (index[0]==8) { //if we are in digit coordinate system transform to global
202 Transform8to4(xyz,index);
203 Transform4to3(xyz,index);
207 void AliTPCParamSR::CRXYZtoXYZ(Float_t *xyz,
208 const Int_t §or, const Int_t & padrow, Int_t option) const
210 //transform relative coordinates to absolute
211 Bool_t rel = ( (option&2)!=0);
212 Int_t index[2]={sector,padrow};
213 if (rel==kTRUE) Transform4to3(xyz,index);//if the position is relative to pad row
214 Transform2to1(xyz,index);
217 void AliTPCParamSR::XYZtoCRXYZ(Float_t *xyz,
218 Int_t §or, Int_t & padrow, Int_t option) const
220 //transform global position to the position relative to the sector padrow
221 //if option=0 X calculate absolute calculate sector
222 //if option=1 X absolute use input sector
223 //if option=2 X relative to pad row calculate sector
224 //if option=3 X relative use input sector
225 //!!!!!!!!! WE start to calculate rows from row = 0
227 Bool_t rel = ( (option&2)!=0);
229 //option 0 and 2 means that we don't have information about sector
230 if ((option&1)==0) Transform0to1(xyz,index); //we calculate sector number
233 Transform1to2(xyz,index);
234 Transform2to3(xyz,index);
235 //if we store relative position calculate position relative to pad row
236 if (rel==kTRUE) Transform3to4(xyz,index);
241 Float_t AliTPCParamSR::GetPrimaryLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
245 Float_t padlength=GetPadPitchLength(index[1]);
246 Float_t a1=TMath::Sin(angle[0]);
248 Float_t a2=TMath::Sin(angle[1]);
250 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
251 return length*fNPrimLoss;
254 Float_t AliTPCParamSR::GetTotalLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
258 Float_t padlength=GetPadPitchLength(index[1]);
259 Float_t a1=TMath::Sin(angle[0]);
261 Float_t a2=TMath::Sin(angle[1]);
263 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
264 return length*fNTotalLoss;
269 void AliTPCParamSR::GetClusterSize(Float_t *x, Int_t *index, Float_t */*angle*/, Int_t /*mode*/, Float_t *sigma)
272 //return cluster sigma2 (x,y) for particle at position x
273 // in this case x coordinata is in drift direction
274 //and y in pad row direction
275 //we suppose that input coordinate system is digit system
278 Float_t lx[3] = {x[0],x[1],x[2]};
279 Int_t li[3] = {index[0],index[1],index[2]};
281 // Float_t sigmadiff;
285 xx = lx[2]; //calculate drift length in cm
287 sigma[0]+= xx*GetDiffL()*GetDiffL();
288 sigma[1]+= xx*GetDiffT()*GetDiffT();
292 //sigma[0]=sigma[1]=0;
293 if (GetTimeRF()!=0) sigma[0]+=GetTimeRF()->GetSigma()*GetTimeRF()->GetSigma();
294 if ( (index[1]<fNInnerSector) &&(GetInnerPRF()!=0))
295 sigma[1]+=GetInnerPRF()->GetSigmaX()*GetInnerPRF()->GetSigmaX();
296 if ( (index[1]>=fNInnerSector) &&(index[2]<fNRowUp1) && (GetOuter1PRF()!=0))
297 sigma[1]+=GetOuter1PRF()->GetSigmaX()*GetOuter1PRF()->GetSigmaX();
298 if( (index[1]>=fNInnerSector) &&(index[2]>=fNRowUp1) && (GetOuter2PRF()!=0))
299 sigma[1]+=GetOuter2PRF()->GetSigmaX()*GetOuter2PRF()->GetSigmaX();
302 sigma[0]/= GetZWidth()*GetZWidth();
303 sigma[1]/=GetPadPitchWidth(index[0])*GetPadPitchWidth(index[0]);
309 void AliTPCParamSR::GetSpaceResolution(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/,
310 Float_t /*amplitude*/, Int_t /*mode*/, Float_t */*sigma*/)
317 Float_t AliTPCParamSR::GetAmp(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/)
325 Float_t * AliTPCParamSR::GetAnglesAccMomentum(Float_t *x, Int_t * index, Float_t* momentum, Float_t *angle)
328 //calculate angle of track to padrow at given position
329 // for given magnetic field and momentum of the particle
332 TransformTo2(x,index);
333 AliDetectorParam::GetAnglesAccMomentum(x,index,momentum,angle);
334 Float_t addangle = TMath::ASin(x[1]/GetPadRowRadii(index[1],index[2]));
340 Bool_t AliTPCParamSR::Update()
343 if (AliTPCParam::Update()==kFALSE) return kFALSE;
346 Float_t firstrow = fInnerRadiusLow + 1.575;
347 for( i= 0;i<fNRowLow;i++)
349 Float_t x = firstrow + fInnerPadPitchLength*(Float_t)i;
351 // number of pads per row
352 // Float_t y = (x-0.5*fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount-
353 // fInnerPadPitchWidth/2.;
354 // 0 and fNRowLow+1 reserved for cross talk rows
355 fYInner[i+1] = x*tan(fInnerAngle/2.)-fInnerWireMount;
356 //fNPadsLow[i] = 1+2*(Int_t)(y/fInnerPadPitchWidth) ;
357 fNPadsLow[i] = AliTPCROC::Instance()->GetNPads(0,i) ; // ROC implement
360 fYInner[0]=(fPadRowLow[0]-fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
361 fYInner[fNRowLow+1]=(fPadRowLow[fNRowLow-1]+fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
362 firstrow = fOuterRadiusLow + 1.6;
363 for(i=0;i<fNRowUp;i++)
366 Float_t x = firstrow + fOuter1PadPitchLength*(Float_t)i;
368 // Float_t y =(x-0.5*fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
369 // fOuterPadPitchWidth/2.;
370 fYOuter[i+1]= x*tan(fOuterAngle/2.)-fOuterWireMount;
371 //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
372 fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
374 fLastWireUp1=fPadRowUp[i] +0.625;
375 firstrow = fPadRowUp[i] + 0.5*(fOuter1PadPitchLength+fOuter2PadPitchLength);
380 Float_t x = firstrow + fOuter2PadPitchLength*(Float_t)(i-64);
382 //Float_t y =(x-0.5*fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
383 // fOuterPadPitchWidth/2.;
384 //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
385 fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
387 fYOuter[i+1] = fPadRowUp[i]*tan(fOuterAngle/2.)-fOuterWireMount;
390 fYOuter[0]=(fPadRowUp[0]-fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
391 fYOuter[fNRowUp+1]=(fPadRowUp[fNRowUp-1]+fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
392 fNtRows = fNInnerSector*fNRowLow+fNOuterSector*fNRowUp;
396 Float_t AliTPCParamSR::GetYInner(Int_t irow) const
398 return fYInner[irow];
400 Float_t AliTPCParamSR::GetYOuter(Int_t irow) const
402 return fYOuter[irow];
405 void AliTPCParamSR::Streamer(TBuffer &R__b)
407 // Stream an object of class AliTPC.
409 if (R__b.IsReading()) {
410 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
411 // TObject::Streamer(R__b);
412 AliTPCParam::Streamer(R__b);
413 // if (R__v < 2) return;
415 if (gGeoManager) ReadGeoMatrices();
417 R__b.WriteVersion(AliTPCParamSR::IsA());
418 //TObject::Streamer(R__b);
419 AliTPCParam::Streamer(R__b);
422 Int_t AliTPCParamSR::CalcResponseFast(Float_t* xyz, Int_t * index, Int_t row)
425 //calculate bin response as function of the input position -x
426 //return number of valid response bin
428 //we suppose that coordinate is expressed in float digits
429 // it's mean coordinate system 8
430 //xyz[0] - electron position w.r.t. pad center, normalized to pad length,
431 //xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
432 //xyz[3] - electron time in float time bin format
433 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
434 Error("AliTPCParamSR", "response function was not adjusted");
438 const Int_t kpadn = 500;
439 const Float_t kfpadn = 500.;
440 const Int_t ktimen = 500;
441 const Float_t kftimen = 500.;
442 const Int_t kpadrn = 500;
443 const Float_t kfpadrn = 500.;
447 static Float_t prfinner[2*kpadrn][5*kpadn]; //pad divided by 50
448 static Float_t prfouter1[2*kpadrn][5*kpadn]; //prfouter division
449 static Float_t prfouter2[2*kpadrn][5*kpadn];
450 static Float_t kTanMax =0;
452 static Float_t rftime[5*ktimen]; //time division
453 static Int_t blabla=0;
454 static Float_t zoffset=0;
455 static Float_t zwidth=0;
456 static Float_t zoffset2=0;
457 static TH1F * hdiff=0;
458 static TH1F * hdiff1=0;
459 static TH1F * hdiff2=0;
461 if (blabla==0) { //calculate Response function - only at the begginning
462 kTanMax = TMath::ATan(10.*TMath::DegToRad());
463 hdiff =new TH1F("prf_diff","prf_diff",10000,-1,1);
464 hdiff1 =new TH1F("no_repsonse1","no_response1",10000,-1,1);
465 hdiff2 =new TH1F("no_response2","no_response2",10000,-1,1);
468 zoffset = GetZOffset();
470 zoffset2 = zoffset/zwidth;
471 for (Int_t i=0;i<5*ktimen;i++){
472 rftime[i] = fTimeRF->GetRF(((i-2.5*kftimen)/kftimen)*zwidth+zoffset);
474 for (Int_t i=0;i<5*kpadn;i++){
475 for (Int_t j=0;j<2*kpadrn;j++){
477 fInnerPRF->GetPRF((i-2.5*kfpadn)/kfpadn
478 *fInnerPadPitchWidth,(j-kfpadrn)/kfpadrn*fInnerPadPitchLength);
480 fOuter1PRF->GetPRF((i-2.5*kfpadn)/kfpadn
481 *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter1PadPitchLength);
485 fOuter2PRF->GetPRF((i-2.5*kfpadn)/kfpadn
486 *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter2PadPitchLength);
489 } // the above is calculated only once
491 // calculate central padrow, pad, time
492 Int_t npads = GetNPads(index[1],index[3]-1);
493 Int_t cpadrow = index[2]; // electrons are here
494 Int_t cpad = TMath::Nint(xyz[1]);
495 Int_t ctime = TMath::Nint(xyz[2]+zoffset2+xyz[3]-GetNTBinsL1());
497 Float_t dpadrow = xyz[0];
498 Float_t dpad = xyz[1]-cpad;
499 Float_t dtime = xyz[2]+zoffset2+xyz[3]-ctime-GetNTBinsL1();
502 Int_t maxt =GetMaxTBin();
507 if (row>=0) { //if we are interesting about given pad row
508 fpadrow = row-cpadrow;
509 lpadrow = row-cpadrow;
511 fpadrow = (index[2]>1) ? -1 :0;
512 lpadrow = (index[2]<GetNRow(index[1])-1) ? 1:0;
515 Int_t fpad = (cpad > -npads/2+1) ? -2: -npads/2-cpad;
516 Int_t lpad = (cpad < npads/2-2) ? 2: npads/2-1-cpad;
517 Int_t ftime = (ctime>1) ? -2: -ctime;
518 Int_t ltime = (ctime<maxt-2) ? 2: maxt-ctime-1;
520 // cross talk from long pad to short one
521 if(row==fNRowUp1 && fpadrow==-1) {
522 dpadrow *= fOuter2PadPitchLength;
523 dpadrow += fOuterWWPitch;
524 dpadrow /= fOuter1PadPitchLength;
526 // cross talk from short pad to long one
527 if(row==fNRowUp1+1 && fpadrow==1){
528 dpadrow *= fOuter1PadPitchLength;
529 if(dpadrow < 0.) dpadrow = -1.; //protection against 3rd wire
530 dpadrow += fOuterWWPitch;
531 dpadrow /= fOuter2PadPitchLength;
536 Int_t apadrow = TMath::Nint((dpadrow-fpadrow)*kfpadrn+kfpadrn);
537 for (Int_t ipadrow = fpadrow; ipadrow<=lpadrow;ipadrow++){
538 if ( (apadrow<0) || (apadrow>=2*kpadrn))
540 // pad angular correction
541 Float_t angle = kTanMax*2.*(cpad+0.5)/Float_t(npads);
542 Float_t dpadangle =0;
543 if (index[1]<fNInnerSector){
544 dpadangle = angle*dpadrow*fInnerPadPitchLength/fInnerPadPitchWidth;
547 if(row < fNRowUp1+1){
548 dpadangle = angle*dpadrow*fOuter1PadPitchLength/fOuterPadPitchWidth;
551 dpadangle = angle*dpadrow*fOuter2PadPitchLength/fOuterPadPitchWidth;
554 if (ipadrow==0) dpadangle *=-1;
556 // Int_t apad= TMath::Nint((dpad-fpad)*kfpadn+2.5*kfpadn);
557 Int_t apad= TMath::Nint((dpad+dpadangle-fpad)*kfpadn+2.5*kfpadn);
558 for (Int_t ipad = fpad; ipad<=lpad;ipad++){
560 if (index[1]<fNInnerSector){
561 cweight=prfinner[apadrow][apad];
564 if(row < fNRowUp1+1){
565 cweight=prfouter1[apadrow][apad];
568 cweight=prfouter2[apadrow][apad];
571 // if (cweight<fResponseThreshold) continue;
572 Int_t atime = TMath::Nint((dtime-ftime)*kftimen+2.5*kftimen);
573 for (Int_t itime = ftime;itime<=ltime;itime++){
574 Float_t cweight2 = cweight*rftime[atime];
575 if (cweight2>fResponseThreshold) {
576 fResponseBin[cindex3++]=cpadrow+ipadrow;
577 fResponseBin[cindex3++]=cpad+ipad;
578 fResponseBin[cindex3++]=ctime+itime;
579 fResponseWeight[cindex++]=cweight2;