]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCTrackHits.cxx
Bugfix in AliPoints2Memory
[u/mrichter/AliRoot.git] / TPC / AliTPCTrackHits.cxx
CommitLineData
b6895dd8 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/*
2289aadd 17$Log$
8217f168 18Revision 1.4 2001/01/10 09:35:27 kowal2
19Changes to obey the coding rules.
20
633d3715 21Revision 1.3 2000/11/02 10:22:50 kowal2
22Logs added
23
b6895dd8 24*/
cf98c13f 25///////////////////////////////////////////////////////////////////////////////
26// //
27// Time Projection Chamber track hits object //
28//
29// Origin: Marian Ivanov , GSI Darmstadt
30//
31// Class for storing simulated AliTPCHits for given track //
32// -average compression comparing to classical ClonesArray is //
33// around 5-7 (depending on the required hit precision) //
34// //
35//Begin_Html
36/*
37<img src="gif/AliTPCTrackHits.gif">
38*/
39//End_Html
40// //
41// //
42///////////////////////////////////////////////////////////////////////////////
43
44#include "TVector3.h"
45#include "TClonesArray.h"
46#include "AliTPCTrackHits.h"
47#include "AliTPC.h"
48
49#include <iostream.h>
50
51
52
53ClassImp(AliTPCTrackHits)
54LClassImp(AliTrackHitsInfo)
55LClassImp(AliTrackHitsParam)
56LClassImp(AliHitInfo)
57
58Int_t AliTrackHitsInfo::fgCounter1 =0;
59Int_t AliTrackHitsInfo::fgCounter2 =0;
60Int_t AliTrackHitsParam::fgCounter1 =0;
61Int_t AliTrackHitsParam::fgCounter2 =0;
62Int_t AliHitInfo::fgCounter1 =0;
63Int_t AliHitInfo::fgCounter2 =0;
64Int_t AliTPCTrackHits::fgCounter1 =0;
65Int_t AliTPCTrackHits::fgCounter2 =0;
66const Double_t AliTPCTrackHits::fgkPrecision=1e-6; //precision
67const Double_t AliTPCTrackHits::fgkPrecision2=1e-20; //precision
68
69
70/************************************************************/
71// Interface classes //
72#include "AliTPCTrackHitsInterfaces.h"
73
74
75
76
77struct AliTPCCurrentHit {
78 AliTPChit fHit;
79 UInt_t fInfoIndex;// - current info pointer
80 UInt_t fParamIndex;// - current param pointer
81 UInt_t fStackIndex; // - current hit stack index
82 Double_t fR; //current Radius
83 Bool_t fStatus; //current status
84};
85
86
87struct AliTPCTempHitInfo {
88 enum { fkStackSize = 10000};
89 AliTPCTempHitInfo();
90 void NewParam(Double_t r, Double_t z, Double_t fi, Int_t q);
91 void SetHit(Double_t r, Double_t z, Double_t fi, Int_t q);
92 Double_t * GetPosition(Int_t index){return &fPositionStack[index*3];}
93 void UpdateParam(Double_t maxdelta); //recal
94 void Fit2(Double_t fSumY, Double_t fSumYX, Double_t fSumYX2,
95 Double_t fSumX, Double_t fSumX2, Double_t fSumX3,
96 Double_t fSumX4, Int_t n,
97 Double_t &a, Double_t &b, Double_t &c);
98 void Fit(AliTrackHitsParam * param);
99 Double_t fSumDr; //
100 Double_t fSumDr2; //
101 Double_t fSumDr3; //
102 Double_t fSumDr4; //
103 Double_t fSumDFi; //
104 Double_t fSumDFiDr; //
105 Double_t fSumDFiDr2;//
106 Double_t fSumDZ; //
107 Double_t fSumDZDr; //
108 Double_t fSumDZDr2; //
109 Double_t fOldR; //previos r
110 Double_t fPositionStack[3*fkStackSize]; //position stack
111 UInt_t fQStack[fkStackSize]; //Q stack
112 UInt_t fStackIndex; //current stack index
113 UInt_t fInfoIndex; //current track info index
114 UInt_t fParamIndex; //current track parameters index
115 AliTrackHitsInfo * fInfo; //current track info
116 AliTrackHitsParam * fParam; //current track param
117};
118
119
120AliTPCTempHitInfo::AliTPCTempHitInfo()
121{
122 //
123 //set to default value
124 fSumDr=fSumDr2=fSumDr3=fSumDr4=
125 fSumDFi=fSumDFiDr=fSumDFiDr2=
126 fSumDZ=fSumDZDr=fSumDZDr2=0;
127 fStackIndex = 0;
128 fInfoIndex = 0;
129 fParamIndex = 0;
130}
131
132
133void AliTPCTempHitInfo::NewParam(Double_t r, Double_t z, Double_t fi, Int_t q)
134{
135 //
136 //reset stack and sum parameters
137 //store line initial point
138 fSumDr=fSumDr2=fSumDr3=fSumDr4=
139 fSumDFi=fSumDFiDr=fSumDFiDr2=
140 fSumDZ=fSumDZDr=fSumDZDr2=0;
141 fStackIndex=0;
142 fParam->fR = r;
143 fOldR = r;
144 fParam->fZ = z;
145 fParam->fFi = fi;
146 fParam->fAn = 0.;
147 fParam->fAd = 0.;
148 fParam->fTheta =0.;
149 fParam->fThetaD =0.;
150 SetHit(r,z,fi,q);
151}
152
153void AliTPCTempHitInfo::SetHit(Double_t r, Double_t z, Double_t fi, Int_t q)
154{
155 //
156 //add hit to the stack
157 //recalculate new estimete of line parameters
158 Double_t *f = GetPosition(fStackIndex);
159 f[0] = r;
160 f[1] = z;
161 f[2] = fi;
162 fQStack[fStackIndex]=q;
163 if (fStackIndex==0) return;
164 Double_t dr = (r-fParam->fR);
165 if (TMath::Abs(dr)<AliTPCTrackHits::fgkPrecision) dr =AliTPCTrackHits::fgkPrecision;
166 Double_t dfi = fi-fParam->fFi;
167 Double_t dz = z -fParam->fZ;
168 Double_t dr2 =dr*dr;
169 Double_t dr3 =dr2*dr;
170 Double_t dr4 =dr3*dr;
171 fSumDr +=dr;
172 fSumDr2+=dr2;
173 fSumDr3+=dr3;
174 fSumDr4+=dr4;
175 fSumDFi +=dfi;
176 fSumDFiDr+=dfi*dr;
177 fSumDFiDr2+=dfi*dr2;
178 fSumDZ +=dz;
179 fSumDZDr+=dz*dr;
180 fSumDZDr2+=dz*dr2;
181
182 //update fit parameters
183 //
184 Double_t det = fSumDr2*fSumDr4-fSumDr3*fSumDr3;
185 if (TMath::Abs(det)<AliTPCTrackHits::fgkPrecision2) return;
186 if ( ( fStackIndex>1 ) ){
187 fParam->fAn = (fSumDr4*fSumDFiDr-fSumDr3*fSumDFiDr2)/det;
188 fParam->fAd = (fSumDr2*fSumDFiDr2-fSumDr3*fSumDFiDr)/det;
189 }
190 else
191 fParam->fAn = fSumDFiDr/fSumDr2;
192 if ( ( fStackIndex>1 ) ){
193 fParam->fTheta = (fSumDr4*fSumDZDr-fSumDr3*fSumDZDr2)/det;
194 fParam->fThetaD= (fSumDr2*fSumDZDr2-fSumDr3*fSumDZDr)/det;
195 }
196 else
197 fParam->fTheta = fSumDZDr/fSumDr2;
198}
199
200
201void AliTPCTempHitInfo::UpdateParam(Double_t maxdelta)
202{
203 //recalc parameters not fixing origin point
204 if (fStackIndex>5){
205 Double_t a,b,c;
206 a=b=c=0;
207 Fit2(fSumDFi, fSumDFiDr, fSumDFiDr2, fSumDr,fSumDr2,fSumDr3,fSumDr4,
208 fStackIndex, a,b,c);
209 if (TMath::Abs(a)<maxdelta){
210 fParam->fFi +=a/fParam->fR;
211 fParam->fAn = b;
212 fParam->fAd = c;
213 }
214 Fit2(fSumDZ, fSumDZDr, fSumDZDr2, fSumDr,fSumDr2,fSumDr3,fSumDr4,
215 fStackIndex, a,b,c) ;
216 if (TMath::Abs(a)<maxdelta){
217 fParam->fZ +=a;
218 fParam->fTheta = b;
219 fParam->fThetaD = c;
220 }
221 }
222
223}
224void AliTPCTempHitInfo::Fit2(Double_t fSumY, Double_t fSumYX, Double_t fSumYX2,
225 Double_t fSumX, Double_t fSumX2, Double_t fSumX3,
226 Double_t fSumX4, Int_t n,
227 Double_t &a, Double_t &b, Double_t &c)
228{
229 //fit of second order
230 Double_t det =
231 n* (fSumX2*fSumX4-fSumX3*fSumX3) -
232 fSumX* (fSumX*fSumX4-fSumX3*fSumX2)+
233 fSumX2* (fSumX*fSumX3-fSumX2*fSumX2);
234
235 if (TMath::Abs(det)> AliTPCTrackHits::fgkPrecision) {
236 a =
237 (fSumY * (fSumX2*fSumX4-fSumX3*fSumX3)-
238 fSumX *(fSumYX*fSumX4-fSumYX2*fSumX3)+
239 fSumX2*(fSumYX*fSumX3-fSumYX2*fSumX2))/det;
240 b=
241 (n*(fSumYX*fSumX4-fSumX3*fSumYX2)-
242 fSumY*(fSumX*fSumX4-fSumX3*fSumX2)+
243 fSumX2*(fSumX*fSumYX2-fSumYX*fSumX2))/det;
244 c=
245 (n*(fSumX2*fSumYX2-fSumYX*fSumX3)-
246 fSumX*(fSumX*fSumYX2-fSumYX*fSumX2)+
247 fSumY*(fSumX*fSumX3-fSumX2*fSumX2))/det;
248 }
249}
250
251void AliTPCTempHitInfo::Fit(AliTrackHitsParam * param)
252{
253 // fit fixing first and the last point
254 //result stored in new param
255 Double_t dx2 = (GetPosition(fStackIndex))[0]-fParam->fR;
256 Double_t det = fSumDr4+dx2*fSumDr2-2*dx2*fSumDr3;
257 if ( (TMath::Abs(det)> AliTPCTrackHits::fgkPrecision) &&
258 ((TMath::Abs(dx2)> AliTPCTrackHits::fgkPrecision))){
259 Double_t dfi2 = (GetPosition(fStackIndex))[1]-fParam->fFi;
260 param->fAd = (fSumDFiDr2+dfi2*fSumDr-dx2*fSumDFiDr-dfi2*fSumDr3/dx2)/det;
261 param->fAn = (dfi2-param->fAd*dx2*dx2)/dx2;
262
263 Double_t dz2 = (GetPosition(fStackIndex))[1]-fParam->fZ;
264 param->fTheta = (fSumDZDr2+dz2*fSumDr-dx2*fSumDZDr-dz2*fSumDr3/dx2)/det;
265 param->fTheta = (dz2-param->fAd*dx2*dx2)/dx2;
266 }
267
268}
269
270
271
272
273AliTPCTrackHits::AliTPCTrackHits()
274{
275 //
276 //default constructor
277 //
278 const Float_t kHitPrecision=0.002; //default precision for hit position in cm
279 const Float_t kStep =0.003; //30 mum step
280 const UShort_t kMaxDistance =100; //maximum distance 100
281
282 fPrecision=kHitPrecision; //precision in cm
283 fStep = kStep; //step size
284 fMaxDistance = kMaxDistance; //maximum distance
285 fTempInfo =0;
286 fTrackHitsInfo = new AliObjectArray("AliTrackHitsInfo");
287 fTrackHitsParam = new AliObjectArray("AliTrackHitsParam");
633d3715 288 fHitsPosAndQ = new TArrayOfArrayVStack("AliHitInfo");
cf98c13f 289 fCurrentHit = new AliTPCCurrentHit;
290 fgCounter1++;
291 fgCounter2++;
292
293}
294
295AliTPCTrackHits::~AliTPCTrackHits()
296{
297 //
298 //default destructor
299 //
300 if (fTrackHitsInfo) delete fTrackHitsInfo;
301 if (fTrackHitsParam) delete fTrackHitsParam;
302 if (fHitsPosAndQ) delete fHitsPosAndQ;
303 if (fCurrentHit) delete fCurrentHit;
304 if (fTempInfo) delete fTempInfo;
305 fgCounter1--;
306}
307
308void AliTPCTrackHits::Clear()
309{
310 //
311 //clear object
633d3715 312 fTrackHitsInfo->Clear();
313 fTrackHitsParam->Clear();
314 //fTrackHitsInfo->Resize(0);
315 //fTrackHitsParam->Resize(0);
cf98c13f 316 fHitsPosAndQ->Clear();
317
318 if (fTempInfo){
319 delete fTempInfo;
320 fTempInfo =0;
321 }
322}
323
324
325void AliTPCTrackHits::AddHitKartez(Int_t volumeID, Int_t trackID, Double_t x,
326 Double_t y, Double_t z,Int_t q)
327{
328 Double_t r = TMath::Sqrt(x*x+y*y);
329 Double_t fi = TMath::ACos(x/r);
330 if (y<0) fi*=-1.;
331 AddHit(volumeID,trackID,r,z,fi,q);
332}
333
334void AliTPCTrackHits::AddHit(Int_t volumeID, Int_t trackID,
335 Double_t r, Double_t z, Double_t fi, Int_t q)
336{
337 //
338 Bool_t diff=kFALSE;
339 if (!fTempInfo) { //initialsation of track
340 fTempInfo = new AliTPCTempHitInfo;
341 //
342 if (fTrackHitsInfo->GetCapacity()<10) fTrackHitsInfo->Reserve(10);
343 fTrackHitsInfo->Resize(1);
344 fTempInfo->fInfoIndex =0;
345 if (fTrackHitsParam->GetCapacity()<100) fTrackHitsParam->Reserve(100);
346 fTrackHitsParam->Resize(1);
347 //
348 fTempInfo->fInfo =
349 (AliTrackHitsInfo*) (fTrackHitsInfo->At(0));
350 fTempInfo->fInfo->fVolumeID = volumeID;
351 fTempInfo->fInfo->fTrackID = trackID;
352 fTempInfo->fInfo->fHitParamIndex =0;
353 fTempInfo->fInfoIndex = 0;
354 //
355 fTempInfo->fParam =
356 (AliTrackHitsParam*) (fTrackHitsParam->At(0));
357 fTempInfo->fParamIndex = 0;
358 fTempInfo->NewParam(r,z,fi,q);
359 return;
360 }
361
362 Int_t size = fHitsPosAndQ->GetSize();
363 if (size>(Int_t)fTempInfo->fParamIndex) {
364 fTempInfo->fParamIndex++;
365 if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize())
366 fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);
367 fTempInfo->fParam =
368 (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));
369 fTempInfo->NewParam(r,z,fi,q);
370 return;
371 }
372
373
374 // if new volume or new trackID
375 if ( (volumeID!=fTempInfo->fInfo->fVolumeID) ||
376 (trackID!=fTempInfo->fInfo->fTrackID)){
377 diff=kTRUE;
378
379 FlushHitStack(kTRUE);
380
381 fTempInfo->fInfoIndex++;
382 if (fTempInfo->fInfoIndex+1>fTrackHitsInfo->GetSize())
383 fTrackHitsInfo->Resize(fTempInfo->fInfoIndex+1);
384 fTempInfo->fInfo =
385 (AliTrackHitsInfo*) (fTrackHitsInfo->At(fTempInfo->fInfoIndex));
386 fTempInfo->fInfo->fVolumeID = volumeID;
387 fTempInfo->fInfo->fTrackID = trackID;
388 fTempInfo->fInfo->fHitParamIndex =fTempInfo->fParamIndex+1;
389 // FlushHitStack(kTRUE);
390
391 fTempInfo->fParamIndex++;
392 if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize())
393 fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);
394 fTempInfo->fParam =
395 (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));
396 fTempInfo->NewParam(r,z,fi,q);
397 return;
398 }
399
400 //calculate current fit precission to next point
401 AliTrackHitsParam &param = *(fTempInfo->fParam);
402 Double_t dd=0;
403 Double_t dl=0;
404 Double_t ratio=0;
405 Double_t dr,dz,dfi,ddz,ddfi;
406 Double_t drhit,ddl;
407 dr=dz=dfi=ddz=ddfi=0;
408 drhit = r-fTempInfo->fOldR;
409 {
410 //Double_t dfi2 = param.fAn+2*param.fAd*(r-param.fR);
411 Double_t dfi2 = param.fAn;
412 dfi2*=dfi2*fTempInfo->fOldR*fTempInfo->fOldR;
413 //Double_t ddz2 = param.fTheta+2*param.fThetaD*(r-param.fR);
414 Double_t ddz2 = param.fTheta;
415 ddz2*=ddz2;
416 ratio = TMath::Sqrt(1.+ dfi2+ ddz2);
417 }
8217f168 418 //
419
420
421 dl = (TMath::Abs(drhit*ratio/fStep)<32000) ? fStep * Short_t(TMath::Nint(drhit*ratio/fStep)):0;
cf98c13f 422 ddl = dl - drhit*ratio;
423 fTempInfo->fOldR += dl/ratio;
424
425 if (fTempInfo->fStackIndex>2){
426 dr = r-param.fR;
427 dz = z-param.fZ;
428 dfi = fi-param.fFi;
429 ddz = dr*param.fTheta+dr*dr*param.fThetaD-dz;
430 ddfi= dr*param.fAn+dr*dr*param.fAd-dfi;
431 dd = TMath::Sqrt(ddz*ddz+r*r*ddfi*ddfi+ddl*ddl);
432 //
433 }
434 //safety factor 1.25
435 if ( ( (dd*1.25>fPrecision) ) ||
436 (fTempInfo->fStackIndex+4>fTempInfo->fkStackSize) ||
437 (TMath::Abs(dl/fStep)>fMaxDistance) )
438 diff=kTRUE;
439 else{
440 fTempInfo->fStackIndex++;
441 fTempInfo->SetHit(r,z,fi,q);
442 return;
443 }
444 //if parameter changed
445 if (FlushHitStack(kFALSE)){ //if full buffer flushed
446 fTempInfo->fParamIndex++;
447 if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize())
448 fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);
449 fTempInfo->fParam =
450 (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));
451 fTempInfo->NewParam(r,z,fi,q);
452 }
453 else{
454 fTempInfo->fStackIndex++;
455 fTempInfo->SetHit(r,z,fi,q);
456 }
457}
458
459Bool_t AliTPCTrackHits::FlushHitStack(Bool_t force)
460{
461 //
462 //write fHitsPosAndQ information from the stack to te arrays
463 if (!fTempInfo) return kFALSE;
464 Int_t size = fHitsPosAndQ->GetSize();
465
466 if ( (size>0)&&(size!=(Int_t)fTempInfo->fParamIndex)) return kFALSE;
467
468 if (fHitsPosAndQ->Push(fTempInfo->fStackIndex+1)!=fTempInfo->fParamIndex){
469 cout<<"internal error - contact MI\n";
470 return kFALSE;
471 }
472 AliHitInfo * info;
473
474 AliTrackHitsParam & param = *(fTempInfo->fParam);
475 //recalculate track parameter not fixing first point
476 fTempInfo->UpdateParam(fStep/4.);
477 //fTempInfo->Fit(fTempInfo->fParam); //- fixing the first and the last point
478
479 Double_t oldr = param.fR;
480 //cout<<"C3"<<fTempInfo->fStackIndex<<"\n"<<flush;
481 UInt_t i;
482 Double_t dd;
483 for (i=0; i <= fTempInfo->fStackIndex; i++){
484 Double_t * position = fTempInfo->GetPosition(i);
485 Double_t dr = position[0]-oldr;
486 Double_t ratio;
487 {
488 //Double_t dfi2 = param.fAn+2*param.fAd*(position[0]-param.fR);
489 Double_t dfi2 = param.fAn;
490 dfi2*=dfi2*oldr*oldr;
491 //Double_t ddz2 = param.fTheta+2*param.fThetaD*(position[0]-param.fR);
492 Double_t ddz2 = param.fTheta;
493 ddz2*=ddz2;
494 ratio = TMath::Sqrt(1.+ dfi2+ ddz2);
495 }
496
8217f168 497 Double_t dl = (TMath::Abs(dr*ratio/fStep)<32000) ? fStep*(Short_t)TMath::Nint(dr*ratio/fStep):0;
cf98c13f 498 dr = dl/ratio;
499 oldr+=dr;
500 //calculate precission
501 AliTrackHitsParam &param = *(fTempInfo->fParam);
502 //real deltas
503 Double_t dr1= position[0]-param.fR;
504 Double_t dz = position[1]-param.fZ;
505
506 Double_t dfi = position[2]-param.fFi;
507 //extrapolated deltas
508 Double_t dr2 = oldr-param.fR;
509 Double_t ddr = dr2-dr1;
510 Double_t ddz = dr2*param.fTheta+dr2*dr2*param.fThetaD-dz;
511 Double_t ddfi= dr2*param.fAn+dr2*dr2*param.fAd-dfi;
512 dd = TMath::Sqrt(ddz*ddz+oldr*oldr*ddfi*ddfi+ddr*ddr);
513
514 if ( (dd>fPrecision) ){
515 if (i==0){
516 param.fAn = 0;
517 param.fAd = 0;
518 param.fTheta =0;
519 param.fThetaD =0;
520 Double_t ddz = dr2*param.fTheta+dr2*dr2*param.fThetaD-dz;
521 Double_t ddfi= dr2*param.fAn+dr2*dr2*param.fAd-dfi;
522 dl = 0;
523 dd = TMath::Sqrt(ddz*ddz+oldr*oldr*ddfi*ddfi+ddr*ddr);
524 }
525 else
526 break;
527 }
528
529 info = (AliHitInfo*)(fHitsPosAndQ->At(fTempInfo->fParamIndex,i));
8217f168 530 info->fHitDistance = (TMath::Abs(dl/fStep)<32000) ?Short_t(TMath::Nint(dl/fStep)):0;
cf98c13f 531 info->fCharge = Short_t(fTempInfo->fQStack[i]);
532 /*
533 cout<<"C2";
534 cout<<" "<<fTempInfo->fStackIndex<<" \t";
535 cout<<" "<<i<<" \t";
536 cout<<position[0]<<"\t";
537 cout<<position[1]<<"\t";
538 cout<<position[2]<<"\t";
539 cout<<param.fAn<<"\t";
540 cout<<param.fTheta<<"\t";
541 cout<<dr1<<"\t"<<ddr<<"\t"<<ddz<<"\t"<<ddfi<<"\t"<<dd<<"\n"<<flush;
542 */
543 }
544
545 if (i<=fTempInfo->fStackIndex){ //if previous iteration not succesfull
546 fHitsPosAndQ->Resize(fTempInfo->fParamIndex,i);
547 //
548 fTempInfo->fParamIndex++;
549 if (fTempInfo->fParamIndex+1>fTrackHitsParam->GetSize())
550 fTrackHitsParam->Resize(fTempInfo->fParamIndex+1);
551 fTempInfo->fParam =
552 (AliTrackHitsParam*) (fTrackHitsParam->At(fTempInfo->fParamIndex));
553 Double_t * p = fTempInfo->GetPosition(i);
554 UInt_t index2 = fTempInfo->fStackIndex;
555 fTempInfo->NewParam(p[0],p[1],p[2],fTempInfo->fQStack[i]);
556 if (i+1<=index2) FlushHitStack2(i+1,index2);
557
558 if (force) return FlushHitStack(kTRUE);
559 return kFALSE;
560 }
561 return kTRUE;
562}
563
564
565void AliTPCTrackHits::FlushHitStack2(Int_t index1, Int_t index2)
566{
567 //
568 // second iteration flush stack
569 // call only for hits where first iteration were not succesfully interpolated
570 Double_t * positionstack = new Double_t[3*(index2-index1+1)];
571 UInt_t * qstack = new UInt_t[index2-index1+1];
572 memcpy(positionstack, &fTempInfo->fPositionStack[3*index1],
573 (3*(index2-index1+1))*sizeof(Double_t));
574 memcpy(qstack, &fTempInfo->fQStack[index1],(index2-index1+1)*sizeof(UInt_t));
575 Double_t *p = positionstack;
576 for (Int_t j=0; j<=index2-index1;j++){
577 fTempInfo->fStackIndex++;
578 fTempInfo->SetHit(p[3*j+0],p[3*j+1],p[3*j+2],qstack[j]);
579 }
580 delete []positionstack;
581 delete []qstack;
582}
583
584
585
586
587
588
589
590
591
592Bool_t AliTPCTrackHits::First()
593{
594 //
595 //set Current hit for the first hit
596 //
597 AliTrackHitsInfo *info = (AliTrackHitsInfo *)fTrackHitsInfo->At(0);
598 AliTrackHitsParam *param = (AliTrackHitsParam *)fTrackHitsParam->At(0);
599 AliHitInfo * hinfo = (AliHitInfo *)fHitsPosAndQ->At(0,0);
600
601 if (!(info) || !(param) || !(hinfo) ) {
602 fCurrentHit->fStatus = kFALSE;
603 return kFALSE;
604 }
605
606 fCurrentHit->fInfoIndex = 0;
607 fCurrentHit->fParamIndex = 0;
608 fCurrentHit->fStackIndex = 0;
609
610 fCurrentHit->fHit.fSector = info->fVolumeID;
611 fCurrentHit->fHit.SetTrack(info->fTrackID);
612 fCurrentHit->fHit.SetX(param->fR*TMath::Cos(param->fFi));
613 fCurrentHit->fHit.SetY(param->fR*TMath::Sin(param->fFi));
614 fCurrentHit->fHit.SetZ(param->fZ);
615 fCurrentHit->fHit.fQ = hinfo->fCharge;
616
617 fCurrentHit->fR = param->fR;
618
619 return fCurrentHit->fStatus = kTRUE;
620}
621
622Bool_t AliTPCTrackHits::Next()
623{
624 //
625 //
626 if (!(fCurrentHit->fStatus))
627 return kFALSE;
628
629 fCurrentHit->fStackIndex++;
630 AliHitInfo * hinfo = (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex,
631 fCurrentHit->fStackIndex);
632 AliTrackHitsInfo *info = (AliTrackHitsInfo *)fTrackHitsInfo->At(fCurrentHit->fInfoIndex);
633 AliTrackHitsParam *param = (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex);
634
635 if (!hinfo) {
636 hinfo = (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex+1, 0);
637 if (!hinfo)
638 return fCurrentHit->fStatus = kFALSE;
639 if (hinfo){
640 fCurrentHit->fParamIndex++;
641 fCurrentHit->fStackIndex = 0;
642 param = (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex);
643 if (!param)
644 return fCurrentHit->fStatus = kFALSE;
645 fCurrentHit->fR = param->fR;
646
647 if ((fCurrentHit->fInfoIndex+1<fTrackHitsInfo->GetSize())
648 &&((info+1)->fHitParamIndex<=fCurrentHit->fParamIndex)){
649 fCurrentHit->fInfoIndex++;
650 info = (AliTrackHitsInfo *)fTrackHitsInfo->At(fCurrentHit->fInfoIndex);
651 if (!info)
652 return fCurrentHit->fStatus = kFALSE;
653 fCurrentHit->fHit.fSector = info->fVolumeID;
654 fCurrentHit->fHit.SetTrack(info->fTrackID);
655 }
656 }
657 }
658 Double_t ratio;
659 {
660 // Double_t dfi2 = param->fAn+2*param->fAd*(fCurrentHit->fR-param->fR);
661 Double_t dfi2 = param->fAn;
662 dfi2*=dfi2*fCurrentHit->fR*fCurrentHit->fR;
663 // Double_t ddz2 = param->fTheta+2*param->fThetaD*(fCurrentHit->fR-param->fR);
664 Double_t ddz2 = param->fTheta;
665 ddz2*=ddz2;
666 ratio = TMath::Sqrt(1.+ dfi2+ ddz2);
667 }
668
669 fCurrentHit->fHit.fQ = hinfo->fCharge;
670 fCurrentHit->fR += fStep*hinfo->fHitDistance/ratio;
671 Double_t dR = fCurrentHit->fR - param->fR;
672 //Double_t dR =0;
673 Double_t fi = param->fFi + (param->fAn*dR+param->fAd*dR*dR);
674 Double_t z = param->fZ + (param->fTheta*dR+param->fThetaD*dR*dR);
675
676 fCurrentHit->fHit.SetX(fCurrentHit->fR*TMath::Cos(fi));
677 fCurrentHit->fHit.SetY(fCurrentHit->fR*TMath::Sin(fi));
678 fCurrentHit->fHit.SetZ(z);
679 return kTRUE;
680}
681
682AliTPChit * AliTPCTrackHits::GetHit()
683{
684 //
685 return (fCurrentHit->fStatus)? &fCurrentHit->fHit:0;
686 //return &fCurrentHit->fHit;
687
688}
689
690AliTrackHitsParam * AliTPCTrackHits::GetParam()
691{
692 //
693 return (fCurrentHit->fStatus)? (AliTrackHitsParam *)fTrackHitsParam->At(fCurrentHit->fParamIndex) :0;
694}
695
696AliHitInfo * AliTPCTrackHits::GetHitInfo()
697{
698 //
699 return (fCurrentHit->fStatus)?
700 (AliHitInfo *)fHitsPosAndQ->At(fCurrentHit->fParamIndex,fCurrentHit->fStackIndex) :0;
701}
702
703
704