]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HLT/comp/AliL3Modeller.cxx
First version, stand alone detector
[u/mrichter/AliRoot.git] / HLT / comp / AliL3Modeller.cxx
CommitLineData
735e167e 1//$Id$
2
029912b7 3// Author: Anders Vestbo <mailto:vestbo@fi.uib.no>
735e167e 4//*-- Copyright &copy ASV
5
4994a25d 6#include "AliL3StandardIncludes.h"
735e167e 7
8#include "AliL3Modeller.h"
9#include "AliL3MemHandler.h"
8487f697 10#ifdef use_aliroot
11#include "AliL3FileHandler.h"
12#endif
735e167e 13#include "AliL3TrackArray.h"
14#include "AliL3ModelTrack.h"
15#include "AliL3DigitData.h"
16#include "AliL3Transform.h"
17
4994a25d 18#if GCCVERSION == 3
19using namespace std;
20#endif
21
4a838220 22//_____________________________________________________________
23// AliL3Modeller
24//
25// Class for modeling TPC data.
26//
27// This performs the cluster finding, based on track parameters.
28// Basically it propagates the tracks to all padrows, and looks
29// for a corresponding cluster. For the moment only cog is calculated,
30// and no deconvolution is done.
735e167e 31
32ClassImp(AliL3Modeller)
33
34AliL3Modeller::AliL3Modeller()
35{
36 fMemHandler=0;
37 fTracks=0;
4a838220 38 fTrackThreshold=0;
39 fPadOverlap=0;
2357bb38 40 SetOverlap();
41 SetTrackThreshold();
735e167e 42}
43
44
45AliL3Modeller::~AliL3Modeller()
46{
47 if(fMemHandler)
48 delete fMemHandler;
49 if(fTracks)
50 delete fTracks;
735e167e 51}
52
4994a25d 53void AliL3Modeller::Init(Int_t slice,Int_t patch,Char_t *trackdata,Char_t *path,Bool_t houghtracks,Bool_t binary)
735e167e 54{
55 fSlice = slice;
56 fPatch = patch;
2357bb38 57
be6ddb10 58 sprintf(fPath,"%s",path);
029912b7 59
735e167e 60 fTracks = new AliL3TrackArray("AliL3ModelTrack");
92a876e2 61
62 Char_t fname[100];
735e167e 63 AliL3MemHandler *file = new AliL3MemHandler();
8487f697 64 //sprintf(fname,"%s/tracks_tr_%d_0.raw",trackdata,fSlice); //output tracks from the tracker (no merging)
65 sprintf(fname,"%s/tracks_ho_%d.raw",trackdata,fSlice);
92a876e2 66 if(!file->SetBinaryInput(fname))
735e167e 67 {
92a876e2 68 cerr<<"AliL3Modeller::Init : Error opening trackfile: "<<fname<<endl;
735e167e 69 return;
70 }
71 file->Binary2TrackArray(fTracks);
72 file->CloseBinaryInput();
73 delete file;
74
4a838220 75 if(houghtracks)
76 cout<<"AliL3Modeller is assuming local hough tracksegments!"<<endl;
77 else
78 cout<<"AliL3Modeller is assuming global tracks!"<<endl;
79
80 if(!houghtracks)
81 fTracks->QSort();
82
735e167e 83 for(Int_t i=0; i<fTracks->GetNTracks(); i++)
84 {
85 AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i);
86 if(!track) continue;
87 track->Init(fSlice,fPatch);
4a838220 88
89 //Only if the tracks has been merged across sector boundaries:
90 //if(!houghtracks)
91 //track->Rotate(fSlice,kTRUE); //!!!!!!!!!!!!!!!!!!!
92
735e167e 93 track->CalculateHelix();
94 }
95
96 CalculateCrossingPoints();
4a838220 97
2357bb38 98 CheckForOverlaps();
99
8487f697 100 UInt_t ndigits=0;
101 AliL3DigitRowData *digits=0;
102#ifdef use_aliroot
103 fMemHandler = new AliL3FileHandler();
2357bb38 104 fMemHandler->Init(slice,patch);
8487f697 105 if(binary == kFALSE)
106 {
2357bb38 107 sprintf(fname,"%s/digitfile.root",fPath);
8487f697 108 fMemHandler->SetAliInput(fname);
109 digits = fMemHandler->AliDigits2Memory(ndigits);
110 }
2357bb38 111 else
112 {
113 sprintf(fname,"%sdigits_%d_%d.raw",fPath,fSlice,fPatch);
114 if(!fMemHandler->SetBinaryInput(fname))
115 {
116 cerr<<"AliL3Modeller::Init : Error opening file "<<fname<<endl;
117 return;
118 }
119 digits=(AliL3DigitRowData*)fMemHandler->CompBinary2Memory(ndigits);
120 }
8487f697 121#else
735e167e 122 fMemHandler = new AliL3MemHandler();
2357bb38 123 fMemHandler->Init(slice,patch);
8487f697 124 if(binary == kFALSE)
735e167e 125 {
8487f697 126 cerr<<"AliL3Modeller::Init : Compile with AliROOT if you want rootfile as input"<<endl;
735e167e 127 return;
128 }
8487f697 129 else
130 {
131 sprintf(fname,"%sdigits_%d_%d.raw",fPath,fSlice,fPatch);
132 if(!fMemHandler->SetBinaryInput(fname))
133 {
134 cerr<<"AliL3Modeller::Init : Error opening file "<<fname<<endl;
135 return;
136 }
137 }
138 digits=(AliL3DigitRowData*)fMemHandler->CompBinary2Memory(ndigits);
139#endif
735e167e 140
141 SetInputData(digits);
142}
143
029912b7 144void AliL3Modeller::FindClusters()
735e167e 145{
146
147 if(!fTracks)
148 {
149 cerr<<"AliL3Modeller::Process : No tracks"<<endl;
150 return;
151 }
152 if(!fRowData)
153 {
154 cerr<<"AliL3Modeller::Process : No data "<<endl;
155 return;
156 }
157
158 AliL3DigitRowData *rowPt = fRowData;
159 AliL3DigitData *digPt=0;
160
4a838220 161 Int_t ntimes = AliL3Transform::GetNTimeBins()+1;
162 Int_t npads = AliL3Transform::GetNPads(AliL3Transform::GetLastRow(fPatch))+1;//Max num of pads.
029912b7 163 Int_t bounds = ntimes*npads;
164 Digit *row = new Digit[bounds];
735e167e 165
166 Int_t seq_charge;
029912b7 167 Int_t pad,time,index;
735e167e 168 Short_t charge;
169 Cluster cluster;
170
4a838220 171 for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++)
735e167e 172 {
8487f697 173 fCurrentPadRow = i;
735e167e 174 memset((void*)row,0,ntimes*npads*sizeof(Digit));
175 digPt = (AliL3DigitData*)rowPt->fDigitData;
2357bb38 176 //cout<<"Loading row "<<i<<" with "<<(Int_t)rowPt->fNDigit<<" digits"<<endl;
735e167e 177 for(UInt_t j=0; j<rowPt->fNDigit; j++)
178 {
179 pad = digPt[j].fPad;
180 time = digPt[j].fTime;
181 charge = digPt[j].fCharge;
182 row[ntimes*pad+time].fCharge = charge;
183 row[ntimes*pad+time].fUsed = kFALSE;
029912b7 184 //cout<<"Row "<<i<<" pad "<<pad<<" time "<<time<<" charge "<<charge<<endl;
735e167e 185 }
186
187 for(Int_t k=0; k<fTracks->GetNTracks(); k++)
188 {
189 AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(k);
190 if(!track) continue;
029912b7 191
192 if(track->GetPadHit(i)<0 || track->GetTimeHit(i)<0 || track->GetOverlap(i)>=0)
95a00d93 193 {
4a838220 194 track->SetCluster(i,0,0,0,0,0,0); //The track has left the patch.
95a00d93 195 continue;
196 }
735e167e 197
198 Int_t hitpad = (Int_t)rint(track->GetPadHit(i));
199 Int_t hittime = (Int_t)rint(track->GetTimeHit(i));
029912b7 200 //cout<<"Checking track on row "<<i<<" with pad "<<hitpad<<" time "<<hittime<<endl;
735e167e 201 pad = hitpad;
202 time = hittime;
203 Int_t padsign=-1;
204 Int_t timesign=-1;
205
206 memset(&cluster,0,sizeof(Cluster));
207
4a838220 208 Int_t npads=0;
735e167e 209 while(1)//Process this padrow
210 {
4a838220 211 if(pad < 0 || pad >= AliL3Transform::GetNPads(i))
029912b7 212 {
92a876e2 213 //cout<<"Pad = "<<pad<<" on row "<<i<<endl;
4a838220 214 FillCluster(track,&cluster,i,npads);
029912b7 215 break;
216 }
735e167e 217 seq_charge=0;
218 timesign=-1;
219 time = hittime;
4a838220 220
735e167e 221 while(1) //Process sequence on this pad:
222 {
029912b7 223 if(time < 0) break;
224 index = ntimes*pad + time;
225 if(index < 0 || index >= bounds)
226 {
227 cerr<<"AliL3Modeller::FindClusters : Index out of range : "<<index
228 <<" on row "<<i<<" pad "<<pad<<" time "<<time<<endl;
229 break;
230 }
4a838220 231
029912b7 232 charge = row[index].fCharge;
2357bb38 233 if(charge==0 && timesign==-1) //zero charge on this timebin, perform checks:
029912b7 234 {
2357bb38 235 if(seq_charge==0 && abs(time-hittime) <= fTimeOverlap) //No charge found on this pad, look further.
029912b7 236 {
237 time--;
238 continue;
239 }
2357bb38 240 else //Boundary reached, or we have found one end of the sequence,->start looking in the other time direction
029912b7 241 {
242 time = hittime+1;
243 timesign=1;
244 continue;
245 }
246 }
2357bb38 247 else if(charge==0 && timesign==1)//zero charge on this timebin, perform checks:
029912b7 248 {
2357bb38 249 if(seq_charge==0 && abs(time-hittime) <= fTimeOverlap)//No charge found on this pad, look further
029912b7 250 {
251 time++;
252 continue;
253 }
2357bb38 254 else //Boundary reached, or we have found the other end of the sequence, stop looking on this pad.
92a876e2 255 {
2357bb38 256 //if(fCurrentPadRow==31)
257 // cerr<<"Breaking off at pad "<<pad<<" and time "<<time<<endl;
92a876e2 258 break;
259 }
029912b7 260 }
261
2357bb38 262 if(row[ntimes*pad+time].fUsed==kTRUE) //Don't use digits several times. This leads to mult. rec.tracks.
263 {
264 time += timesign;
265 continue;
266 }
735e167e 267
268 seq_charge += charge;
2357bb38 269
270 //Update the cluster parameters with this timebin
735e167e 271 cluster.fTime += time*charge;
272 cluster.fPad += pad*charge;
273 cluster.fCharge += charge;
274 cluster.fSigmaY2 += pad*pad*charge;
275 cluster.fSigmaZ2 += time*time*charge;
276
277 row[ntimes*pad+time].fUsed = kTRUE;
278 time += timesign;
279 }
4a838220 280
281
2357bb38 282 if(seq_charge)//There was something on this pad, so keep looking on the neighbouring pad
4a838220 283 {
284 pad += padsign;
285 npads++;
286 }
735e167e 287 else //Nothing more on this pad, goto next pad
288 {
289 if(padsign==-1)
95a00d93 290 {
2357bb38 291 if(cluster.fCharge==0 && abs(pad-hitpad) <= fPadOverlap && pad > 0)
95a00d93 292 {
293 pad--; //In this case, we haven't found anything yet,
294 } //so we will try to expand our search within the natural boundaries.
295 else
296 {
297 pad=hitpad+1;
298 padsign=1;
299 }
300 continue;
301 }
029912b7 302
303 else if(padsign==1)
95a00d93 304 {
2357bb38 305 if(cluster.fCharge==0 && abs(pad-hitpad) <= fPadOverlap && pad < AliL3Transform::GetNPads(i)-2)
029912b7 306 {
307 pad++; //In this case, we haven't found anything yet,
308 continue; //so we will try to expand our search within the natural boundaries.
309 }
310 else //We are out of range, or cluster if finished.
311 {
2357bb38 312 //if(fCurrentPadRow==31)
313 //cout<<"Out of range; charge "<<cluster.fCharge<<" paddiff "<<abs(pad-hitpad)<<endl;
4a838220 314 FillCluster(track,&cluster,i,npads);
029912b7 315 break;
316 }
95a00d93 317 }
735e167e 318 else //Nothing more in this cluster
319 {
2357bb38 320 //if(fCurrentPadRow==31)
029912b7 321 //cout<<"Filling final cluster"<<endl;
4a838220 322 FillCluster(track,&cluster,i,npads);
735e167e 323 break;
324 }
325 }
735e167e 326 }
029912b7 327 //cout<<"done"<<endl;
735e167e 328 }
95a00d93 329 FillZeros(rowPt,row);
735e167e 330 fMemHandler->UpdateRowPointer(rowPt);
735e167e 331 }
332 delete [] row;
029912b7 333 cout<<"done processing"<<endl;
334
335
336 //Debug:
337 for(Int_t i=0; i<fTracks->GetNTracks(); i++)
338 {
339 AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i);
340 if(!track) continue;
4a838220 341 if(track->GetNClusters() != AliL3Transform::GetNRows(fPatch))
342 cerr<<endl<<"Mismatching hitcounts; nclusters: "<<track->GetNClusters()<<" nrows "<<AliL3Transform::GetNRows(fPatch)<<endl<<endl;
029912b7 343 }
344
345}
346
4a838220 347void AliL3Modeller::FillCluster(AliL3ModelTrack *track,Cluster *cluster,Int_t row,Int_t npads)
029912b7 348{
349 if(cluster->fCharge==0)
350 {
4a838220 351 track->SetCluster(row,0,0,0,0,0,0);
029912b7 352 return;
353 }
354 Float_t fcharge = (Float_t)cluster->fCharge;
355 Float_t fpad = ((Float_t)cluster->fPad/fcharge);
356 Float_t ftime = ((Float_t)cluster->fTime/fcharge);
357 Float_t sigmaY2,sigmaZ2;
358 CalcClusterWidth(cluster,sigmaY2,sigmaZ2);
4a838220 359 track->SetCluster(row,fpad,ftime,fcharge,sigmaY2,sigmaZ2,npads);
735e167e 360}
361
95a00d93 362void AliL3Modeller::FillZeros(AliL3DigitRowData *rowPt,Digit *row)
363{
364 //Fill zero where data has been used.
365
4a838220 366 Int_t ntimes = AliL3Transform::GetNTimeBins()+1;
95a00d93 367 AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData;
368 for(UInt_t j=0; j<rowPt->fNDigit; j++)
369 {
370 Int_t pad = digPt[j].fPad;
371 Int_t time = digPt[j].fTime;
372 if(row[ntimes*pad+time].fUsed==kTRUE)
373 digPt[j].fCharge = 0;
374 }
375}
376
be6ddb10 377void AliL3Modeller::WriteRemaining()
95a00d93 378{
379 //Write remaining (nonzero) digits to file.
380
95a00d93 381 AliL3DigitRowData *rowPt;
382 rowPt = (AliL3DigitRowData*)fRowData;
383 Int_t digitcount=0;
4a838220 384 Int_t ndigits[(AliL3Transform::GetNRows(fPatch))];
385 for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++)
95a00d93 386 {
387 AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData;
4a838220 388 ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]=0;
95a00d93 389 for(UInt_t j=0; j<rowPt->fNDigit; j++)
390 {
391 if(digPt[j].fCharge==0) continue;
392 digitcount++;
4a838220 393 ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]++;
95a00d93 394 }
4a838220 395 //cout<<"Difference "<<(int)ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]<<" "<<(int)rowPt->fNDigit<<endl;
95a00d93 396 fMemHandler->UpdateRowPointer(rowPt);
397 }
398
4a838220 399 Int_t size = digitcount*sizeof(AliL3DigitData) + AliL3Transform::GetNRows(fPatch)*sizeof(AliL3DigitRowData);
95a00d93 400 Byte_t *data = new Byte_t[size];
401 memset(data,0,size);
402 AliL3DigitRowData *tempPt = (AliL3DigitRowData*)data;
403 rowPt = (AliL3DigitRowData*)fRowData;
404
4a838220 405 for(Int_t i=AliL3Transform::GetFirstRow(fPatch); i<=AliL3Transform::GetLastRow(fPatch); i++)
95a00d93 406 {
407 Int_t localcount=0;
408 tempPt->fRow = i;
4a838220 409 tempPt->fNDigit = ndigits[(i-AliL3Transform::GetFirstRow(fPatch))];
95a00d93 410 AliL3DigitData *digPt = (AliL3DigitData*)rowPt->fDigitData;
411 for(UInt_t j=0; j<rowPt->fNDigit; j++)
412 {
413 if(digPt[j].fCharge==0) continue;
4a838220 414 if(localcount >= ndigits[(i-AliL3Transform::GetFirstRow(fPatch))])
95a00d93 415 {
416 cerr<<"AliL3Modeller::WriteRemaining : Digitarray out of range!!"<<endl;
417 return;
418 }
419 tempPt->fDigitData[localcount].fCharge = digPt[j].fCharge;
420 tempPt->fDigitData[localcount].fPad = digPt[j].fPad;
421 tempPt->fDigitData[localcount].fTime = digPt[j].fTime;
4a838220 422
95a00d93 423 localcount++;
424 }
4a838220 425 if(ndigits[(i-AliL3Transform::GetFirstRow(fPatch))] != localcount)
95a00d93 426 {
427 cerr<<"AliL3Modeller::WriteRemaining : Mismatch in digitcount"<<endl;
428 return;
429 }
430 fMemHandler->UpdateRowPointer(rowPt);
431 Byte_t *tmp = (Byte_t*)tempPt;
4a838220 432 Int_t size = sizeof(AliL3DigitRowData) + ndigits[(i-AliL3Transform::GetFirstRow(fPatch))]*sizeof(AliL3DigitData);
95a00d93 433 tmp += size;
434 tempPt = (AliL3DigitRowData*)tmp;
435 }
436
be6ddb10 437 Char_t fname[100];
95a00d93 438 AliL3MemHandler *mem = new AliL3MemHandler();
2357bb38 439 sprintf(fname,"%s/comp/remains_%d_%d.raw",fPath,fSlice,fPatch);
be6ddb10 440 mem->SetBinaryOutput(fname);
4a838220 441 mem->Memory2CompBinary((UInt_t)AliL3Transform::GetNRows(fPatch),(AliL3DigitRowData*)data);
95a00d93 442 mem->CloseBinaryOutput();
443 delete mem;
029912b7 444 delete [] data;
95a00d93 445}
446
447
735e167e 448void AliL3Modeller::CalculateCrossingPoints()
449{
450 cout<<"Calculating crossing points on "<<fTracks->GetNTracks()<<" tracks"<<endl;
451 if(!fTracks)
452 {
453 cerr<<"AliL3Modeller::CalculateCrossingPoints(): No tracks"<<endl;
454 return;
455 }
456 Float_t hit[3];
029912b7 457
029912b7 458 Int_t sector,row;
4a838220 459 for(Int_t i=AliL3Transform::GetLastRow(fPatch); i>=AliL3Transform::GetFirstRow(fPatch); i--)
735e167e 460 {
461 for(Int_t j=0; j<fTracks->GetNTracks(); j++)
462 {
463 AliL3ModelTrack *track = (AliL3ModelTrack*)fTracks->GetCheckedTrack(j);
464 if(!track) continue;
4a838220 465
95a00d93 466 if(!track->GetCrossingPoint(i,hit))
735e167e 467 {
029912b7 468 cerr<<"AliL3Modeller::CalculateCrossingPoints : Track "<<j<<" does not intersect row "<<i<<" :"<<endl<<
469 "First point "<<track->GetFirstPointX()<<
470 " nhits "<<track->GetNHits()<<endl;//" tgl "<<track->GetTgl()<<" psi "<<track->GetPsi()<<" charge "<<track->GetCharge()<<endl;
471 //"Center "<<track->GetCenterX()<<" "<<track->GetCenterY()<<endl<<endl<<
472 //"--------"<<endl;
473 fTracks->Remove(j);
735e167e 474 continue;
475 }
2357bb38 476 //cout<<"X "<<hit[0]<<" Y "<<hit[1]<<" Z "<<hit[2]<<" tgl "<<track->GetTgl()<<endl;
95a00d93 477
4a838220 478 AliL3Transform::Slice2Sector(fSlice,i,sector,row);
479 AliL3Transform::Local2Raw(hit,sector,row);
2357bb38 480 //cout<<"Pad "<<hit[1]<<" time "<<hit[2]<<" in sector "<<sector<<" row "<<row<<endl;
4a838220 481 if(hit[1]<0 || hit[1]>AliL3Transform::GetNPads(i) ||
482 hit[2]<0 || hit[2]>AliL3Transform::GetNTimeBins())
95a00d93 483 {//Track is leaving the patch, so flag the track hits (<0)
484 track->SetPadHit(i,-1);
485 track->SetTimeHit(i,-1);
486 continue;
487 }
4a838220 488
489
735e167e 490 track->SetPadHit(i,hit[1]);
491 track->SetTimeHit(i,hit[2]);
95a00d93 492
735e167e 493 //if(hit[1]<0 || hit[2]>445)
95a00d93 494 //if(hit[2]<0 || hit[2]>445)
495 //cout<<"pad "<<hit[1]<<" time "<<hit[2]<<" pt "<<track->GetPt()<<" psi "<<track->GetPsi()<<" tgl "<<track->GetTgl()<<" firstpoint "<<track->GetFirstPointX()<<" "<<track->GetFirstPointY()<<" "<<track->GetFirstPointZ()<<endl;
735e167e 496 //cout<<"Crossing pad "<<hit[1]<<" time "<<hit[2]<<endl;
497 }
498 }
499 fTracks->Compress();
029912b7 500 cout<<"And there are "<<fTracks->GetNTracks()<<" tracks remaining"<<endl;
735e167e 501}
502
503void AliL3Modeller::CheckForOverlaps()
504{
505 //Flag the tracks that overlap
506
029912b7 507 cout<<"Checking for overlaps...";
2357bb38 508 Int_t counter=0;
735e167e 509 for(Int_t i=0; i<fTracks->GetNTracks(); i++)
510 {
511 AliL3ModelTrack *track1 = (AliL3ModelTrack*)fTracks->GetCheckedTrack(i);
512 if(!track1) continue;
513 for(Int_t j=i+1; j<fTracks->GetNTracks(); j++)
514 {
515 AliL3ModelTrack *track2 = (AliL3ModelTrack*)fTracks->GetCheckedTrack(j);
516 if(!track2) continue;
4a838220 517 for(Int_t k=AliL3Transform::GetFirstRow(fPatch); k<=AliL3Transform::GetLastRow(fPatch); k++)
735e167e 518 {
95a00d93 519 if(track1->GetPadHit(k)<0 || track1->GetTimeHit(k)<0 ||
520 track2->GetPadHit(k)<0 || track2->GetTimeHit(k)<0)
521 continue;
2357bb38 522
523 if(track1->GetOverlap(k)>=0 || track2->GetOverlap(k)>=0) continue;
524
525 if(abs((Int_t)rint(track1->GetPadHit(k))-(Int_t)rint(track2->GetPadHit(k))) <= fPadOverlap &&
526 abs((Int_t)rint(track1->GetTimeHit(k))-(Int_t)rint(track2->GetTimeHit(k))) <= fTimeOverlap)
735e167e 527 {
029912b7 528 track2->SetOverlap(k,i);
2357bb38 529 //track1->SetOverlap(k,j);
530 counter++;
735e167e 531 }
532 }
533 }
534 }
2357bb38 535 cout<<"found "<<counter<<" done"<<endl;
735e167e 536}
537
538
539void AliL3Modeller::CalcClusterWidth(Cluster *cl,Float_t &sigmaY2,Float_t &sigmaZ2)
540{
541
542 Float_t padw,timew;
4a838220 543
544 padw = AliL3Transform::GetPadPitchWidth(fPatch);
545
735e167e 546 Float_t charge = (Float_t)cl->fCharge;
547 Float_t pad = (Float_t)cl->fPad/charge;
548 Float_t time = (Float_t)cl->fTime/charge;
549 Float_t s2 = (Float_t)cl->fSigmaY2/charge - pad*pad;
4a838220 550
551 //Save the sigmas in pad and time:
552
553 sigmaY2 = (s2);// + 1./12);//*padw*padw;
554
555 /*Constants added by offline
556 if(s2 != 0)
735e167e 557 {
4a838220 558 sigmaY2 = sigmaY2*0.108;
559 if(fPatch<3)
560 sigmaY2 = sigmaY2*2.07;
735e167e 561 }
4a838220 562 */
563
735e167e 564 s2 = (Float_t)cl->fSigmaZ2/charge - time*time;
4a838220 565 timew = AliL3Transform::GetZWidth();
566 sigmaZ2 = (s2);// +1./12);//*timew*timew;
567
568
569
570 /*Constants added by offline
571 if(s2 != 0)
735e167e 572 {
4a838220 573 sigmaZ2 = sigmaZ2*0.169;
574 if(fPatch < 3)
575 sigmaZ2 = sigmaZ2*1.77;
735e167e 576 }
4a838220 577 */
735e167e 578}