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()
59 //constructor set the default parameters
62 fFacSigmaPadRow = Float_t(kFacSigmaPadRow);
63 fFacSigmaPad = Float_t(kFacSigmaPad);
64 fFacSigmaTime = Float_t(kFacSigmaTime);
68 AliTPCParamSR::AliTPCParamSR(const AliTPCParamSR ¶m)
79 // copy constructor - dummy
81 fFacSigmaPadRow = param.fFacSigmaPadRow;
83 AliTPCParamSR & AliTPCParamSR::operator =(const AliTPCParamSR & param)
86 // assignment operator - dummy
88 fZLength=param.fZLength;
92 AliTPCParamSR::~AliTPCParamSR()
95 //destructor destroy some dynmicaly alocated variables
96 if (fInnerPRF != 0) delete fInnerPRF;
97 if (fOuter1PRF != 0) delete fOuter1PRF;
98 if (fOuter2PRF != 0) delete fOuter2PRF;
99 if (fTimeRF != 0) delete fTimeRF;
102 void AliTPCParamSR::SetDefault()
104 //set default TPC param
106 AliTPCParam::SetDefault();
109 Int_t AliTPCParamSR::CalcResponse(Float_t* xyz, Int_t * index, Int_t row)
112 //calculate bin response as function of the input position -x
113 //return number of valid response bin
115 //we suppose that coordinate is expressed in float digits
116 // it's mean coordinate system 8
117 //xyz[0] - float padrow xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
118 //xyz[3] - electron time in float time bin format
119 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
120 Error("AliTPCParamSR", "response function was not adjusted");
124 Float_t sfpadrow; // sigma of response function
125 Float_t sfpad; // sigma of
126 Float_t sftime= fFacSigmaTime*fTimeRF->GetSigma()/fZWidth; //3 sigma of time response
127 if (index[1]<fNInnerSector){
128 sfpadrow =fFacSigmaPadRow*fInnerPRF->GetSigmaY()/fInnerPadPitchLength;
129 sfpad =fFacSigmaPad*fInnerPRF->GetSigmaX()/fInnerPadPitchWidth;
133 sfpadrow =fFacSigmaPadRow*fOuter1PRF->GetSigmaY()/fOuter1PadPitchLength;
134 sfpad =fFacSigmaPad*fOuter1PRF->GetSigmaX()/fOuterPadPitchWidth;}
136 sfpadrow =fFacSigmaPadRow*fOuter2PRF->GetSigmaY()/fOuter2PadPitchLength;
137 sfpad =fFacSigmaPad*fOuter2PRF->GetSigmaX()/fOuterPadPitchWidth;
141 Int_t fpadrow = TMath::Max(TMath::Nint(index[2]+xyz[0]-sfpadrow),0); //"first" padrow
142 Int_t fpad = TMath::Nint(xyz[1]-sfpad); //first pad
143 Int_t ftime = TMath::Max(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()-sftime),0); // first time
144 Int_t lpadrow = TMath::Min(TMath::Nint(index[2]+xyz[0]+sfpadrow),fpadrow+19); //"last" padrow
145 lpadrow = TMath::Min(GetNRow(index[1])-1,lpadrow);
146 Int_t lpad = TMath::Min(TMath::Nint(xyz[1]+sfpad),fpad+19); //last pad
147 Int_t ltime = TMath::Min(TMath::Nint(xyz[2]+xyz[3]+GetZOffset()/GetZWidth()+sftime),ftime+19); // last time
148 ltime = TMath::Min(ltime,GetMaxTBin()-1);
150 Int_t npads = GetNPads(index[1],row);
155 if (ftime<0) ftime=0;
157 if (row>=0) { //if we are interesting about given pad row
158 if (fpadrow<=row) fpadrow =row;
161 if (lpadrow>=row) lpadrow = row;
167 Float_t padres[20][20]; //I don't expect bigger number of bins
173 //calculate padresponse function
175 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
176 for (pad = fpad;pad<=lpad;pad++){
177 Float_t dy = (xyz[0]+Float_t(index[2]-padrow));
178 Float_t dx = (xyz[1]+Float_t(pad));
179 if (index[1]<fNInnerSector)
180 padres[padrow-fpadrow][pad-fpad]=fInnerPRF->GetPRF(dx*fInnerPadPitchWidth,dy*fInnerPadPitchLength);
183 padres[padrow-fpadrow][pad-fpad]=fOuter1PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter1PadPitchLength);}
185 padres[padrow-fpadrow][pad-fpad]=fOuter2PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter2PadPitchLength);}}}
186 //calculate time response function
188 for (time = ftime;time<=ltime;time++)
189 timeres[time-ftime]= fTimeRF->GetRF((-xyz[2]-xyz[3]+Float_t(time))*fZWidth);
190 //write over threshold values to stack
191 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
192 for (pad = fpad;pad<=lpad;pad++)
193 for (time = ftime;time<=ltime;time++){
194 cweight = timeres[time-ftime]*padres[padrow-fpadrow][pad-fpad];
195 if (cweight>fResponseThreshold) {
196 fResponseBin[cindex3]=padrow;
197 fResponseBin[cindex3+1]=pad;
198 fResponseBin[cindex3+2]=time;
200 fResponseWeight[cindex]=cweight;
209 void AliTPCParamSR::TransformTo8(Float_t *xyz, Int_t *index) const
212 // transformate point to digit coordinate
214 if (index[0]==0) Transform0to1(xyz,index);
215 if (index[0]==1) Transform1to2(xyz,index);
216 if (index[0]==2) Transform2to3(xyz,index);
217 if (index[0]==3) Transform3to4(xyz,index);
218 if (index[0]==4) Transform4to8(xyz,index);
221 void AliTPCParamSR::TransformTo2(Float_t *xyz, Int_t *index) const
224 //transformate point to rotated coordinate
227 if (index[0]==0) Transform0to1(xyz,index);
228 if (index[0]==1) Transform1to2(xyz,index);
229 if (index[0]==4) Transform4to3(xyz,index);
230 if (index[0]==8) { //if we are in digit coordinate system transform to global
231 Transform8to4(xyz,index);
232 Transform4to3(xyz,index);
236 void AliTPCParamSR::CRXYZtoXYZ(Float_t *xyz,
237 const Int_t §or, const Int_t & padrow, Int_t option) const
239 //transform relative coordinates to absolute
240 Bool_t rel = ( (option&2)!=0);
241 Int_t index[2]={sector,padrow};
242 if (rel==kTRUE) Transform4to3(xyz,index);//if the position is relative to pad row
243 Transform2to1(xyz,index);
246 void AliTPCParamSR::XYZtoCRXYZ(Float_t *xyz,
247 Int_t §or, Int_t & padrow, Int_t option) const
249 //transform global position to the position relative to the sector padrow
250 //if option=0 X calculate absolute calculate sector
251 //if option=1 X absolute use input sector
252 //if option=2 X relative to pad row calculate sector
253 //if option=3 X relative use input sector
254 //!!!!!!!!! WE start to calculate rows from row = 0
256 Bool_t rel = ( (option&2)!=0);
258 //option 0 and 2 means that we don't have information about sector
259 if ((option&1)==0) Transform0to1(xyz,index); //we calculate sector number
262 Transform1to2(xyz,index);
263 Transform2to3(xyz,index);
264 //if we store relative position calculate position relative to pad row
265 if (rel==kTRUE) Transform3to4(xyz,index);
270 Float_t AliTPCParamSR::GetPrimaryLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
274 Float_t padlength=GetPadPitchLength(index[1]);
275 Float_t a1=TMath::Sin(angle[0]);
277 Float_t a2=TMath::Sin(angle[1]);
279 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
280 return length*fNPrimLoss;
283 Float_t AliTPCParamSR::GetTotalLoss(Float_t */*x*/, Int_t *index, Float_t *angle)
287 Float_t padlength=GetPadPitchLength(index[1]);
288 Float_t a1=TMath::Sin(angle[0]);
290 Float_t a2=TMath::Sin(angle[1]);
292 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
293 return length*fNTotalLoss;
298 void AliTPCParamSR::GetClusterSize(Float_t *x, Int_t *index, Float_t */*angle*/, Int_t /*mode*/, Float_t *sigma)
301 //return cluster sigma2 (x,y) for particle at position x
302 // in this case x coordinata is in drift direction
303 //and y in pad row direction
304 //we suppose that input coordinate system is digit system
307 Float_t lx[3] = {x[0],x[1],x[2]};
308 Int_t li[3] = {index[0],index[1],index[2]};
310 // Float_t sigmadiff;
314 xx = lx[2]; //calculate drift length in cm
316 sigma[0]+= xx*GetDiffL()*GetDiffL();
317 sigma[1]+= xx*GetDiffT()*GetDiffT();
321 //sigma[0]=sigma[1]=0;
322 if (GetTimeRF()!=0) sigma[0]+=GetTimeRF()->GetSigma()*GetTimeRF()->GetSigma();
323 if ( (index[1]<fNInnerSector) &&(GetInnerPRF()!=0))
324 sigma[1]+=GetInnerPRF()->GetSigmaX()*GetInnerPRF()->GetSigmaX();
325 if ( (index[1]>=fNInnerSector) &&(index[2]<fNRowUp1) && (GetOuter1PRF()!=0))
326 sigma[1]+=GetOuter1PRF()->GetSigmaX()*GetOuter1PRF()->GetSigmaX();
327 if( (index[1]>=fNInnerSector) &&(index[2]>=fNRowUp1) && (GetOuter2PRF()!=0))
328 sigma[1]+=GetOuter2PRF()->GetSigmaX()*GetOuter2PRF()->GetSigmaX();
331 sigma[0]/= GetZWidth()*GetZWidth();
332 sigma[1]/=GetPadPitchWidth(index[0])*GetPadPitchWidth(index[0]);
338 void AliTPCParamSR::GetSpaceResolution(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/,
339 Float_t /*amplitude*/, Int_t /*mode*/, Float_t */*sigma*/)
346 Float_t AliTPCParamSR::GetAmp(Float_t */*x*/, Int_t */*index*/, Float_t */*angle*/)
354 Float_t * AliTPCParamSR::GetAnglesAccMomentum(Float_t *x, Int_t * index, Float_t* momentum, Float_t *angle)
357 //calculate angle of track to padrow at given position
358 // for given magnetic field and momentum of the particle
361 TransformTo2(x,index);
362 AliDetectorParam::GetAnglesAccMomentum(x,index,momentum,angle);
363 Float_t addangle = TMath::ASin(x[1]/GetPadRowRadii(index[1],index[2]));
369 Bool_t AliTPCParamSR::Update()
372 if (AliTPCParam::Update()==kFALSE) return kFALSE;
375 Float_t firstrow = fInnerRadiusLow + 1.575;
376 for( i= 0;i<fNRowLow;i++)
378 Float_t x = firstrow + fInnerPadPitchLength*(Float_t)i;
380 // number of pads per row
381 // Float_t y = (x-0.5*fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount-
382 // fInnerPadPitchWidth/2.;
383 // 0 and fNRowLow+1 reserved for cross talk rows
384 fYInner[i+1] = x*tan(fInnerAngle/2.)-fInnerWireMount;
385 //fNPadsLow[i] = 1+2*(Int_t)(y/fInnerPadPitchWidth) ;
386 fNPadsLow[i] = AliTPCROC::Instance()->GetNPads(0,i) ; // ROC implement
389 fYInner[0]=(fPadRowLow[0]-fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
390 fYInner[fNRowLow+1]=(fPadRowLow[fNRowLow-1]+fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
391 firstrow = fOuterRadiusLow + 1.6;
392 for(i=0;i<fNRowUp;i++)
395 Float_t x = firstrow + fOuter1PadPitchLength*(Float_t)i;
397 // Float_t y =(x-0.5*fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
398 // fOuterPadPitchWidth/2.;
399 fYOuter[i+1]= x*tan(fOuterAngle/2.)-fOuterWireMount;
400 //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
401 fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
403 fLastWireUp1=fPadRowUp[i] +0.625;
404 firstrow = fPadRowUp[i] + 0.5*(fOuter1PadPitchLength+fOuter2PadPitchLength);
409 Float_t x = firstrow + fOuter2PadPitchLength*(Float_t)(i-64);
411 //Float_t y =(x-0.5*fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
412 // fOuterPadPitchWidth/2.;
413 //fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
414 fNPadsUp[i] = AliTPCROC::Instance()->GetNPads(36,i) ; // ROC implement
416 fYOuter[i+1] = fPadRowUp[i]*tan(fOuterAngle/2.)-fOuterWireMount;
419 fYOuter[0]=(fPadRowUp[0]-fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
420 fYOuter[fNRowUp+1]=(fPadRowUp[fNRowUp-1]+fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
421 fNtRows = fNInnerSector*fNRowLow+fNOuterSector*fNRowUp;
425 Float_t AliTPCParamSR::GetYInner(Int_t irow) const
427 return fYInner[irow];
429 Float_t AliTPCParamSR::GetYOuter(Int_t irow) const
431 return fYOuter[irow];
434 void AliTPCParamSR::Streamer(TBuffer &R__b)
436 // Stream an object of class AliTPC.
438 if (R__b.IsReading()) {
439 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
440 // TObject::Streamer(R__b);
441 AliTPCParam::Streamer(R__b);
442 // if (R__v < 2) return;
444 if (gGeoManager) ReadGeoMatrices();
446 R__b.WriteVersion(AliTPCParamSR::IsA());
447 //TObject::Streamer(R__b);
448 AliTPCParam::Streamer(R__b);
451 Int_t AliTPCParamSR::CalcResponseFast(Float_t* xyz, Int_t * index, Int_t row)
454 //calculate bin response as function of the input position -x
455 //return number of valid response bin
457 //we suppose that coordinate is expressed in float digits
458 // it's mean coordinate system 8
459 //xyz[0] - electron position w.r.t. pad center, normalized to pad length,
460 //xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
461 //xyz[3] - electron time in float time bin format
462 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
463 Error("AliTPCParamSR", "response function was not adjusted");
467 const Int_t kpadn = 500;
468 const Float_t kfpadn = 500.;
469 const Int_t ktimen = 500;
470 const Float_t kftimen = 500.;
471 const Int_t kpadrn = 500;
472 const Float_t kfpadrn = 500.;
476 static Float_t prfinner[2*kpadrn][5*kpadn]; //pad divided by 50
477 static Float_t prfouter1[2*kpadrn][5*kpadn]; //prfouter division
478 static Float_t prfouter2[2*kpadrn][5*kpadn];
479 static Float_t kTanMax =0;
481 static Float_t rftime[5*ktimen]; //time division
482 static Int_t blabla=0;
483 static Float_t zoffset=0;
484 static Float_t zwidth=0;
485 static Float_t zoffset2=0;
486 static TH1F * hdiff=0;
487 static TH1F * hdiff1=0;
488 static TH1F * hdiff2=0;
490 if (blabla==0) { //calculate Response function - only at the begginning
491 kTanMax = TMath::ATan(10.*TMath::DegToRad());
492 hdiff =new TH1F("prf_diff","prf_diff",10000,-1,1);
493 hdiff1 =new TH1F("no_repsonse1","no_response1",10000,-1,1);
494 hdiff2 =new TH1F("no_response2","no_response2",10000,-1,1);
497 zoffset = GetZOffset();
499 zoffset2 = zoffset/zwidth;
500 for (Int_t i=0;i<5*ktimen;i++){
501 rftime[i] = fTimeRF->GetRF(((i-2.5*kftimen)/kftimen)*zwidth+zoffset);
503 for (Int_t i=0;i<5*kpadn;i++){
504 for (Int_t j=0;j<2*kpadrn;j++){
506 fInnerPRF->GetPRF((i-2.5*kfpadn)/kfpadn
507 *fInnerPadPitchWidth,(j-kfpadrn)/kfpadrn*fInnerPadPitchLength);
509 fOuter1PRF->GetPRF((i-2.5*kfpadn)/kfpadn
510 *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter1PadPitchLength);
514 fOuter2PRF->GetPRF((i-2.5*kfpadn)/kfpadn
515 *fOuterPadPitchWidth,(j-kfpadrn)/kfpadrn*fOuter2PadPitchLength);
518 } // the above is calculated only once
520 // calculate central padrow, pad, time
521 Int_t npads = GetNPads(index[1],index[3]-1);
522 Int_t cpadrow = index[2]; // electrons are here
523 Int_t cpad = TMath::Nint(xyz[1]);
524 Int_t ctime = TMath::Nint(xyz[2]+zoffset2+xyz[3]);
526 Float_t dpadrow = xyz[0];
527 Float_t dpad = xyz[1]-cpad;
528 Float_t dtime = xyz[2]+zoffset2+xyz[3]-ctime;
531 Int_t maxt =GetMaxTBin();
536 if (row>=0) { //if we are interesting about given pad row
537 fpadrow = row-cpadrow;
538 lpadrow = row-cpadrow;
540 fpadrow = (index[2]>1) ? -1 :0;
541 lpadrow = (index[2]<GetNRow(index[1])-1) ? 1:0;
544 Int_t fpad = (cpad > -npads/2+1) ? -2: -npads/2-cpad;
545 Int_t lpad = (cpad < npads/2-2) ? 2: npads/2-1-cpad;
546 Int_t ftime = (ctime>1) ? -2: -ctime;
547 Int_t ltime = (ctime<maxt-2) ? 2: maxt-ctime-1;
549 // cross talk from long pad to short one
550 if(row==fNRowUp1 && fpadrow==-1) {
551 dpadrow *= fOuter2PadPitchLength;
552 dpadrow += fOuterWWPitch;
553 dpadrow /= fOuter1PadPitchLength;
555 // cross talk from short pad to long one
556 if(row==fNRowUp1+1 && fpadrow==1){
557 dpadrow *= fOuter1PadPitchLength;
558 if(dpadrow < 0.) dpadrow = -1.; //protection against 3rd wire
559 dpadrow += fOuterWWPitch;
560 dpadrow /= fOuter2PadPitchLength;
565 Int_t apadrow = TMath::Nint((dpadrow-fpadrow)*kfpadrn+kfpadrn);
566 for (Int_t ipadrow = fpadrow; ipadrow<=lpadrow;ipadrow++){
567 if ( (apadrow<0) || (apadrow>=2*kpadrn))
569 // pad angular correction
570 Float_t angle = kTanMax*2.*(cpad+0.5)/Float_t(npads);
571 Float_t dpadangle =0;
572 if (index[1]<fNInnerSector){
573 dpadangle = angle*dpadrow*fInnerPadPitchLength/fInnerPadPitchWidth;
576 if(row < fNRowUp1+1){
577 dpadangle = angle*dpadrow*fOuter1PadPitchLength/fOuterPadPitchWidth;
580 dpadangle = angle*dpadrow*fOuter2PadPitchLength/fOuterPadPitchWidth;
583 if (ipadrow==0) dpadangle *=-1;
585 // Int_t apad= TMath::Nint((dpad-fpad)*kfpadn+2.5*kfpadn);
586 Int_t apad= TMath::Nint((dpad+dpadangle-fpad)*kfpadn+2.5*kfpadn);
587 for (Int_t ipad = fpad; ipad<=lpad;ipad++){
589 if (index[1]<fNInnerSector){
590 cweight=prfinner[apadrow][apad];
593 if(row < fNRowUp1+1){
594 cweight=prfouter1[apadrow][apad];
597 cweight=prfouter2[apadrow][apad];
600 // if (cweight<fResponseThreshold) continue;
601 Int_t atime = TMath::Nint((dtime-ftime)*kftimen+2.5*kftimen);
602 for (Int_t itime = ftime;itime<=ltime;itime++){
603 Float_t cweight2 = cweight*rftime[atime];
604 if (cweight2>fResponseThreshold) {
605 fResponseBin[cindex3++]=cpadrow+ipadrow;
606 fResponseBin[cindex3++]=cpad+ipad;
607 fResponseBin[cindex3++]=ctime+itime;
608 fResponseWeight[cindex++]=cweight2;