]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCParamSR.cxx
Have some new set functions. Remove deconv=true setting from init.
[u/mrichter/AliRoot.git] / TPC / AliTPCParamSR.cxx
CommitLineData
cc80f89e 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/*
17$Log$
b584c7dd 18Revision 1.7 2002/03/18 17:59:13 kowal2
19Chnges in the pad geometry - 3 pad lengths introduced.
20
f03e3423 21Revision 1.6 2002/02/25 11:02:56 kowal2
22Changes towards speeding up the code. Thanks to Marian Ivanov.
23
de61d5d5 24Revision 1.5 2001/12/06 07:49:30 kowal2
25corrected number of pads calculation
26
73477b37 27Revision 1.4 2000/11/02 07:33:15 kowal2
28Improvements of the code.
29
c11cb93f 30Revision 1.3 2000/06/30 12:07:50 kowal2
31Updated from the TPC-PreRelease branch
32
73042f01 33Revision 1.2.4.2 2000/06/14 16:48:24 kowal2
34Parameter setting improved. Removed compiler warnings
35
36Revision 1.2.4.1 2000/06/09 07:55:39 kowal2
37
38Updated defaults
39
40Revision 1.2 2000/04/17 09:37:33 kowal2
41removed obsolete AliTPCDigitsDisplay.C
42
cc80f89e 43Revision 1.1.4.2 2000/04/10 11:36:13 kowal2
44
45New Detector parameters handling class
46
47*/
48
49///////////////////////////////////////////////////////////////////////
50// Manager and of geomety classes for set: TPC //
51// //
52// !sectors are numbered from 0 //
53// !pad rows are numbered from 0 //
54//
55// 27.7. - AliTPCPaaramSr object for TPC
56// TPC with straight pad rows
57// Origin: Marian Ivanov, Uni. of Bratislava, ivanov@fmph.uniba.sk //
58// //
59///////////////////////////////////////////////////////////////////////
60
61
62#include <iostream.h>
63#include <TMath.h>
64#include <TObject.h>
65#include <AliTPCParamSR.h>
73042f01 66#include "AliTPCPRF2D.h"
67#include "AliTPCRF1D.h"
de61d5d5 68#include "TH1.h"
cc80f89e 69
70
71ClassImp(AliTPCParamSR)
72const static Int_t kMaxRows=600;
73const static Float_t kEdgeSectorSpace = 2.5;
73042f01 74const static Float_t kFacSigmaPadRow=3.;
cc80f89e 75const static Float_t kFacSigmaPad=3.;
76const static Float_t kFacSigmaTime=3.;
77
78
79AliTPCParamSR::AliTPCParamSR()
80{
81 //
82 //constructor set the default parameters
83 fInnerPRF=0;
f03e3423 84 fOuter1PRF=0;
85 fOuter2PRF=0;
cc80f89e 86 fTimeRF = 0;
87 fFacSigmaPadRow = Float_t(kFacSigmaPadRow);
88 fFacSigmaPad = Float_t(kFacSigmaPad);
89 fFacSigmaTime = Float_t(kFacSigmaTime);
cc80f89e 90 SetDefault();
91 Update();
92}
93
94AliTPCParamSR::~AliTPCParamSR()
95{
96 //
97 //destructor destroy some dynmicaly alocated variables
98 if (fInnerPRF != 0) delete fInnerPRF;
f03e3423 99 if (fOuter1PRF != 0) delete fOuter1PRF;
100 if (fOuter2PRF != 0) delete fOuter2PRF;
cc80f89e 101 if (fTimeRF != 0) delete fTimeRF;
102}
103
104void AliTPCParamSR::SetDefault()
105{
106 //set default TPC param
107 fbStatus = kFALSE;
108 AliTPCParam::SetDefault();
109}
110
111Int_t AliTPCParamSR::CalcResponse(Float_t* xyz, Int_t * index, Int_t row)
112{
113 //
114 //calculate bin response as function of the input position -x
115 //return number of valid response bin
116 //
117 //we suppose that coordinate is expressed in float digits
118 // it's mean coordinate system 8
119 //xyz[0] - float padrow xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
f03e3423 120 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
cc80f89e 121 Error("AliTPCParamSR", "response function was not adjusted");
122 return -1;
123 }
124
125 Float_t sfpadrow; // sigma of response function
126 Float_t sfpad; // sigma of
127 Float_t sftime= fFacSigmaTime*fTimeRF->GetSigma()/fZWidth; //3 sigma of time response
128 if (index[1]<fNInnerSector){
129 sfpadrow =fFacSigmaPadRow*fInnerPRF->GetSigmaY()/fInnerPadPitchLength;
130 sfpad =fFacSigmaPad*fInnerPRF->GetSigmaX()/fInnerPadPitchWidth;
f03e3423 131 }
132 else{
133 if(row<fNRowUp1){
134 sfpadrow =fFacSigmaPadRow*fOuter1PRF->GetSigmaY()/fOuter1PadPitchLength;
135 sfpad =fFacSigmaPad*fOuter1PRF->GetSigmaX()/fOuterPadPitchWidth;}
136 else{
137 sfpadrow =fFacSigmaPadRow*fOuter2PRF->GetSigmaY()/fOuter2PadPitchLength;
138 sfpad =fFacSigmaPad*fOuter2PRF->GetSigmaX()/fOuterPadPitchWidth;
139 }
cc80f89e 140 }
141
142 Int_t fpadrow = TMath::Max(TMath::Nint(index[2]+xyz[0]-sfpadrow),0); //"first" padrow
143 Int_t fpad = TMath::Nint(xyz[1]-sfpad); //first pad
144 Int_t ftime = TMath::Max(TMath::Nint(xyz[2]+GetZOffset()/GetZWidth()-sftime),0); // first time
145 Int_t lpadrow = TMath::Min(TMath::Nint(index[2]+xyz[0]+sfpadrow),fpadrow+19); //"last" padrow
146 lpadrow = TMath::Min(GetNRow(index[1])-1,lpadrow);
147 Int_t lpad = TMath::Min(TMath::Nint(xyz[1]+sfpad),fpad+19); //last pad
148 Int_t ltime = TMath::Min(TMath::Nint(xyz[2]+GetZOffset()/GetZWidth()+sftime),ftime+19); // last time
149 ltime = TMath::Min(ltime,GetMaxTBin()-1);
de61d5d5 150 //
151 Int_t npads = GetNPads(index[1],row);
152 if (fpad<-npads/2)
153 fpad = -npads/2;
154 if (lpad>npads/2)
155 lpad= npads/2;
156 if (ftime<0) ftime=0;
157 //
cc80f89e 158 if (row>=0) { //if we are interesting about given pad row
159 if (fpadrow<=row) fpadrow =row;
160 else
161 return 0;
162 if (lpadrow>=row) lpadrow = row;
163 else
164 return 0;
165 }
166
167
168 Float_t padres[20][20]; //I don't expect bigger number of bins
169 Float_t timeres[20];
170 Int_t cindex3=0;
171 Int_t cindex=0;
172 Float_t cweight = 0;
173 if (fpadrow>=0) {
174 //calculate padresponse function
175 Int_t padrow, pad;
176 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
177 for (pad = fpad;pad<=lpad;pad++){
f03e3423 178 Float_t dy = (xyz[0]+Float_t(index[2]-padrow));
179 Float_t dx = (xyz[1]+Float_t(pad));
cc80f89e 180 if (index[1]<fNInnerSector)
181 padres[padrow-fpadrow][pad-fpad]=fInnerPRF->GetPRF(dx*fInnerPadPitchWidth,dy*fInnerPadPitchLength);
f03e3423 182 else{
183 if(row<fNRowUp1){
184 padres[padrow-fpadrow][pad-fpad]=fOuter1PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter1PadPitchLength);}
185 else{
186 padres[padrow-fpadrow][pad-fpad]=fOuter2PRF->GetPRF(dx*fOuterPadPitchWidth,dy*fOuter2PadPitchLength);}}}
cc80f89e 187 //calculate time response function
188 Int_t time;
189 for (time = ftime;time<=ltime;time++)
190 timeres[time-ftime]= fTimeRF->GetRF((-xyz[2]+Float_t(time))*fZWidth);
191 //write over threshold values to stack
192 for (padrow = fpadrow;padrow<=lpadrow;padrow++)
193 for (pad = fpad;pad<=lpad;pad++)
194 for (time = ftime;time<=ltime;time++){
195 cweight = timeres[time-ftime]*padres[padrow-fpadrow][pad-fpad];
196 if (cweight>fResponseThreshold) {
197 fResponseBin[cindex3]=padrow;
198 fResponseBin[cindex3+1]=pad;
199 fResponseBin[cindex3+2]=time;
200 cindex3+=3;
201 fResponseWeight[cindex]=cweight;
202 cindex++;
203 }
204 }
205 }
206 fCurrentMax=cindex;
207 return fCurrentMax;
208}
209
210void AliTPCParamSR::TransformTo8(Float_t *xyz, Int_t *index) const
211{
212 //
213 // transformate point to digit coordinate
214 //
215 if (index[0]==0) Transform0to1(xyz,index);
216 if (index[0]==1) Transform1to2(xyz,index);
217 if (index[0]==2) Transform2to3(xyz,index);
218 if (index[0]==3) Transform3to4(xyz,index);
219 if (index[0]==4) Transform4to8(xyz,index);
220}
221
222void AliTPCParamSR::TransformTo2(Float_t *xyz, Int_t *index) const
223{
224 //
225 //transformate point to rotated coordinate
226 //
227 //we suppose that
228 if (index[0]==0) Transform0to1(xyz,index);
229 if (index[0]==1) Transform1to2(xyz,index);
230 if (index[0]==4) Transform4to3(xyz,index);
231 if (index[0]==8) { //if we are in digit coordinate system transform to global
232 Transform8to4(xyz,index);
233 Transform4to3(xyz,index);
234 }
235}
236
237void AliTPCParamSR::CRXYZtoXYZ(Float_t *xyz,
238 const Int_t &sector, const Int_t & padrow, Int_t option) const
239{
240 //transform relative coordinates to absolute
241 Bool_t rel = ( (option&2)!=0);
242 Int_t index[2]={sector,padrow};
243 if (rel==kTRUE) Transform4to3(xyz,index);//if the position is relative to pad row
244 Transform2to1(xyz,index);
245}
246
247void AliTPCParamSR::XYZtoCRXYZ(Float_t *xyz,
248 Int_t &sector, Int_t & padrow, Int_t option) const
249{
250 //transform global position to the position relative to the sector padrow
251 //if option=0 X calculate absolute calculate sector
252 //if option=1 X absolute use input sector
253 //if option=2 X relative to pad row calculate sector
254 //if option=3 X relative use input sector
255 //!!!!!!!!! WE start to calculate rows from row = 0
256 Int_t index[2];
257 Bool_t rel = ( (option&2)!=0);
258
259 //option 0 and 2 means that we don't have information about sector
260 if ((option&1)==0) Transform0to1(xyz,index); //we calculate sector number
261 else
262 index[0]=sector;
263 Transform1to2(xyz,index);
264 Transform2to3(xyz,index);
265 //if we store relative position calculate position relative to pad row
266 if (rel==kTRUE) Transform3to4(xyz,index);
267 sector = index[0];
268 padrow = index[1];
269}
270
271Float_t AliTPCParamSR::GetPrimaryLoss(Float_t *x, Int_t *index, Float_t *angle)
272{
273 //
274 //
275 Float_t padlength=GetPadPitchLength(index[1]);
276 Float_t a1=TMath::Sin(angle[0]);
277 a1*=a1;
278 Float_t a2=TMath::Sin(angle[1]);
279 a2*=a2;
280 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
281 return length*fNPrimLoss;
282}
283
284Float_t AliTPCParamSR::GetTotalLoss(Float_t *x, Int_t *index, Float_t *angle)
285{
286 //
287 //
288 Float_t padlength=GetPadPitchLength(index[1]);
289 Float_t a1=TMath::Sin(angle[0]);
290 a1*=a1;
291 Float_t a2=TMath::Sin(angle[1]);
292 a2*=a2;
293 Float_t length =padlength*TMath::Sqrt(1+a1+a2);
294 return length*fNTotalLoss;
295
296}
297
298
299void AliTPCParamSR::GetClusterSize(Float_t *x, Int_t *index, Float_t *angle, Int_t mode, Float_t *sigma)
300{
301 //
302 //return cluster sigma2 (x,y) for particle at position x
303 // in this case x coordinata is in drift direction
304 //and y in pad row direction
305 //we suppose that input coordinate system is digit system
306
307 Float_t xx;
308 Float_t lx[3] = {x[0],x[1],x[2]};
309 Int_t li[3] = {index[0],index[1],index[2]};
310 TransformTo2(lx,li);
311 // Float_t sigmadiff;
312 sigma[0]=0;
313 sigma[1]=0;
314
315 xx = lx[2]; //calculate drift length in cm
316 if (xx>0) {
317 sigma[0]+= xx*GetDiffL()*GetDiffL();
318 sigma[1]+= xx*GetDiffT()*GetDiffT();
319 }
320
321
322 //sigma[0]=sigma[1]=0;
323 if (GetTimeRF()!=0) sigma[0]+=GetTimeRF()->GetSigma()*GetTimeRF()->GetSigma();
324 if ( (index[1]<fNInnerSector) &&(GetInnerPRF()!=0))
325 sigma[1]+=GetInnerPRF()->GetSigmaX()*GetInnerPRF()->GetSigmaX();
f03e3423 326 if ( (index[1]>=fNInnerSector) &&(index[2]<fNRowUp1) && (GetOuter1PRF()!=0))
327 sigma[1]+=GetOuter1PRF()->GetSigmaX()*GetOuter1PRF()->GetSigmaX();
328 if( (index[1]>=fNInnerSector) &&(index[2]>=fNRowUp1) && (GetOuter2PRF()!=0))
329 sigma[1]+=GetOuter2PRF()->GetSigmaX()*GetOuter2PRF()->GetSigmaX();
cc80f89e 330
331
332 sigma[0]/= GetZWidth()*GetZWidth();
333 sigma[1]/=GetPadPitchWidth(index[0])*GetPadPitchWidth(index[0]);
334}
335
336
337
338
339void AliTPCParamSR::GetSpaceResolution(Float_t *x, Int_t *index, Float_t *angle,
340 Float_t amplitude, Int_t mode, Float_t *sigma)
341{
342 //
343 //
344 //
345
346}
347Float_t AliTPCParamSR::GetAmp(Float_t *x, Int_t *index, Float_t *angle)
348{
349 //
350 //
351 //
352 return 0;
353}
354
355Float_t * AliTPCParamSR::GetAnglesAccMomentum(Float_t *x, Int_t * index, Float_t* momentum, Float_t *angle)
356{
357 //
358 //calculate angle of track to padrow at given position
359 // for given magnetic field and momentum of the particle
360 //
361
362 TransformTo2(x,index);
363 AliDetectorParam::GetAnglesAccMomentum(x,index,momentum,angle);
364 Float_t addangle = TMath::ASin(x[1]/GetPadRowRadii(index[1],index[2]));
365 angle[1] +=addangle;
366 return angle;
367}
368
369
370Bool_t AliTPCParamSR::Update()
371{
cc80f89e 372 Int_t i;
373 if (AliTPCParam::Update()==kFALSE) return kFALSE;
374 fbStatus = kFALSE;
375
f03e3423 376 Float_t firstrow = fInnerRadiusLow + 2.225 ;
377 for( i= 0;i<fNRowLow;i++)
378 {
379 Float_t x = firstrow + fInnerPadPitchLength*(Float_t)i;
380 fPadRowLow[i]=x;
381 // number of pads per row
382 Float_t y = (x-0.5*fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount-
383 fInnerPadPitchWidth/2.;
b584c7dd 384 // 0 and fNRowLow+1 reserved for cross talk rows
385 fYInner[i+1] = x*tan(fInnerAngle/2.)-fInnerWireMount;
f03e3423 386 fNPadsLow[i] = 1+2*(Int_t)(y/fInnerPadPitchWidth) ;
387 }
b584c7dd 388 // cross talk rows
389 fYInner[0]=(fPadRowLow[0]-fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
390 fYInner[fNRowLow+1]=(fPadRowLow[fNRowLow-1]+fInnerPadPitchLength)*tan(fInnerAngle/2.)-fInnerWireMount;
f03e3423 391 firstrow = fOuterRadiusLow + 1.6;
392 for(i=0;i<fNRowUp;i++)
393 {
394 if(i<fNRowUp1){
395 Float_t x = firstrow + fOuter1PadPitchLength*(Float_t)i;
396 fPadRowUp[i]=x;
397 Float_t y =(x-0.5*fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
398 fOuterPadPitchWidth/2.;
b584c7dd 399 fYOuter[i+1]= x*tan(fOuterAngle/2.)-fOuterWireMount;
f03e3423 400 fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
401 if(i==fNRowUp1-1) {
402 fLastWireUp1=fPadRowUp[i] +0.375;
403 firstrow = fPadRowUp[i] + 0.5*(fOuter1PadPitchLength+fOuter2PadPitchLength);
404 }
405 }
406 else
407 {
408 Float_t x = firstrow + fOuter2PadPitchLength*(Float_t)(i-64);
409 fPadRowUp[i]=x;
410Float_t y =(x-0.5*fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount-
411 fOuterPadPitchWidth/2.;
412 fNPadsUp[i] = 1+2*(Int_t)(y/fOuterPadPitchWidth) ;
cc80f89e 413 }
b584c7dd 414 fYOuter[i+1] = fPadRowUp[i]*tan(fOuterAngle/2.)-fOuterWireMount;
f03e3423 415 }
b584c7dd 416 // cross talk rows
417 fYOuter[0]=(fPadRowUp[0]-fOuter1PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
418 fYOuter[fNRowUp+1]=(fPadRowUp[fNRowUp-1]+fOuter2PadPitchLength)*tan(fOuterAngle/2.)-fOuterWireMount;
f03e3423 419 fNtRows = fNInnerSector*fNRowLow+fNOuterSector*fNRowUp;
420 fbStatus = kTRUE;
421 return kTRUE;
422}
423Float_t AliTPCParamSR::GetYInner(Int_t irow) const
424{
425 return fYInner[irow];
426}
427Float_t AliTPCParamSR::GetYOuter(Int_t irow) const
428{
429 return fYOuter[irow];
cc80f89e 430}
cc80f89e 431
cc80f89e 432void AliTPCParamSR::Streamer(TBuffer &R__b)
433{
434 // Stream an object of class AliTPC.
435
436 if (R__b.IsReading()) {
437 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
438 // TObject::Streamer(R__b);
439 AliTPCParam::Streamer(R__b);
440 // if (R__v < 2) return;
441 Update();
442 } else {
443 R__b.WriteVersion(AliTPCParamSR::IsA());
444 //TObject::Streamer(R__b);
445 AliTPCParam::Streamer(R__b);
446 }
447}
de61d5d5 448Int_t AliTPCParamSR::CalcResponseFast(Float_t* xyz, Int_t * index, Int_t row)
449{
450 //
451 //calculate bin response as function of the input position -x
452 //return number of valid response bin
453 //
454 //we suppose that coordinate is expressed in float digits
455 // it's mean coordinate system 8
b584c7dd 456 //xyz[0] - electron position w.r.t. pad center, normalized to pad length,
457 //xyz[1] is float pad (center pad is number 0) and xyz[2] is float time bin
f03e3423 458 if ( (fInnerPRF==0)||(fOuter1PRF==0)||(fOuter2PRF==0) ||(fTimeRF==0) ){
de61d5d5 459 Error("AliTPCParamSR", "response function was not adjusted");
460 return -1;
461 }
462
463 const Int_t padn = 500;
464 const Float_t fpadn = 500.;
465 const Int_t timen = 500;
466 const Float_t ftimen = 500.;
467 const Int_t padrn = 500;
468 const Float_t fpadrn = 500.;
469
470
cc80f89e 471
de61d5d5 472 static Float_t prfinner[2*padrn][5*padn]; //pad divided by 50
f03e3423 473 static Float_t prfouter1[2*padrn][5*padn]; //prfouter division
474 static Float_t prfouter2[2*padrn][5*padn];
475
de61d5d5 476 static Float_t rftime[5*timen]; //time division
477 static Int_t blabla=0;
478 static Float_t zoffset=0;
479 static Float_t zwidth=0;
480 static Float_t zoffset2=0;
481 static TH1F * hdiff=0;
482 static TH1F * hdiff1=0;
483 static TH1F * hdiff2=0;
484
485 if (blabla==0) { //calculate Response function - only at the begginning
486 hdiff =new TH1F("prf_diff","prf_diff",10000,-1,1);
487 hdiff1 =new TH1F("no_repsonse1","no_response1",10000,-1,1);
488 hdiff2 =new TH1F("no_response2","no_response2",10000,-1,1);
489
490 blabla=1;
491 zoffset = GetZOffset();
492 zwidth = fZWidth;
493 zoffset2 = zoffset/zwidth;
494 for (Int_t i=0;i<5*timen;i++){
495 rftime[i] = fTimeRF->GetRF(((i-2.5*ftimen)/ftimen)*zwidth+zoffset);
496 }
497 for (Int_t i=0;i<5*padn;i++){
498 for (Int_t j=0;j<2*padrn;j++){
499 prfinner[j][i] =
500 fInnerPRF->GetPRF((i-2.5*fpadn)/fpadn
501 *fInnerPadPitchWidth,(j-fpadrn)/fpadrn*fInnerPadPitchLength);
f03e3423 502 prfouter1[j][i] =
503 fOuter1PRF->GetPRF((i-2.5*fpadn)/fpadn
504 *fOuterPadPitchWidth,(j-fpadrn)/fpadrn*fOuter1PadPitchLength);
505
506 //
507 prfouter2[j][i] =
508 fOuter2PRF->GetPRF((i-2.5*fpadn)/fpadn
509 *fOuterPadPitchWidth,(j-fpadrn)/fpadrn*fOuter2PadPitchLength);
de61d5d5 510 }
511 }
f03e3423 512 } // the above is calculated only once
513
de61d5d5 514 // calculate central padrow, pad, time
b584c7dd 515 Int_t npads = GetNPads(index[1],index[3]-1);
f03e3423 516 Int_t cpadrow = index[2]; // electrons are here
de61d5d5 517 Int_t cpad = TMath::Nint(xyz[1]);
518 Int_t ctime = TMath::Nint(xyz[2]+zoffset2);
519 //calulate deviation
520 Float_t dpadrow = xyz[0];
521 Float_t dpad = xyz[1]-cpad;
522 Float_t dtime = xyz[2]+zoffset2-ctime;
523 Int_t cindex =0;
524 Int_t cindex3 =0;
525 Int_t maxt =GetMaxTBin();
526
527 Int_t fpadrow;
528 Int_t lpadrow;
529
530 if (row>=0) { //if we are interesting about given pad row
531 fpadrow = row-cpadrow;
532 lpadrow = row-cpadrow;
533 }else{
534 fpadrow = (index[2]>1) ? -1 :0;
535 lpadrow = (index[2]<GetNRow(index[1])-1) ? 1:0;
536 }
537 Int_t fpad = (cpad > -npads/2+1) ? -2: -npads/2-cpad;
538 Int_t lpad = (cpad < npads/2-1) ? 2: npads/2-cpad;
539 Int_t ftime = (ctime>1) ? -2: -ctime;
540 Int_t ltime = (ctime<maxt-2) ? 2: maxt-ctime-1;
541
f03e3423 542 // cross talk from long pad to short one
b584c7dd 543 if(row==fNRowUp1 && fpadrow==-1) {
f03e3423 544 dpadrow *= fOuter2PadPitchLength;
545 dpadrow += fOuterWWPitch;
546 dpadrow /= fOuter1PadPitchLength;
547 }
548 // cross talk from short pad to long one
b584c7dd 549 if(row==fNRowUp1+1 && fpadrow==1){
f03e3423 550 dpadrow *= fOuter1PadPitchLength;
b584c7dd 551 if(dpadrow < 0.) dpadrow = -1.; //protection against 3rd wire
f03e3423 552 dpadrow += fOuterWWPitch;
553 dpadrow /= fOuter2PadPitchLength;
554
555 }
556 // "normal"
557 Int_t apadrow = TMath::Nint((dpadrow-fpadrow)*fpadrn+fpadrn);
de61d5d5 558 for (Int_t ipadrow = fpadrow; ipadrow<=lpadrow;ipadrow++){
559 if ( (apadrow<0) || (apadrow>=2*padrn))
560 continue;
561 Int_t apad= TMath::Nint((dpad-fpad)*fpadn+2.5*fpadn);
562 for (Int_t ipad = fpad; ipad<=lpad;ipad++){
563 Float_t cweight;
564 if (index[1]<fNInnerSector)
565 cweight=prfinner[apadrow][apad];
f03e3423 566 else{
567 if(row < fNRowUp1)
568 cweight=prfouter1[apadrow][apad];
569 else cweight=prfouter2[apadrow][apad];
570 }
571
de61d5d5 572 // if (cweight<fResponseThreshold) continue;
573 Int_t atime = TMath::Nint((dtime-ftime)*ftimen+2.5*ftimen);
574 for (Int_t itime = ftime;itime<=ltime;itime++){
575 Float_t cweight2 = cweight*rftime[atime];
576 if (cweight2>fResponseThreshold) {
577 fResponseBin[cindex3++]=cpadrow+ipadrow;
578 fResponseBin[cindex3++]=cpad+ipad;
579 fResponseBin[cindex3++]=ctime+itime;
580 fResponseWeight[cindex++]=cweight2;
581
582 if (cweight2>100)
583 {
584 printf("Pici pici %d %f %d\n",ipad,dpad,apad);
585 }
586
587 }
588 atime-=timen;
589 }
590 apad-= padn;
591 }
592 apadrow-=padrn;
593 }
594 fCurrentMax=cindex;
595 return fCurrentMax;
596
597}
cc80f89e 598
599
600
c11cb93f 601
602
603
604