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 Revision 1.22 2002/10/22 14:45:36 alibrary
19 Introducing Riostream.h
21 Revision 1.21 2002/02/05 09:12:26 hristov
22 Small mods for gcc 3.02
24 Revision 1.20 2001/11/21 14:47:45 barbera
25 Some unuseful print-out commented out
27 Revision 1.19 2001/11/21 10:49:07 barbera
28 Bug correction suggested by Rene done
30 Revision 1.18 2001/11/20 15:46:17 barbera
31 Point coordinated are calculated in cylindrical reference frame once and for all at the beginning of tracking V1
33 Revision 1.10.2.1 2001/10/24 07:26:04 hristov
34 All the changes from the head are merged with the release
36 Revision 1.14 2001/10/24 07:19:57 hristov
37 Some pointer correctly initialised in one of the constructors
39 Revision 1.13 2001/10/21 19:17:12 hristov
40 Several pointers were set to zero in the default constructors to avoid memory management problems
42 Revision 1.12 2001/10/19 21:32:35 nilsen
43 Minor changes to remove compliation warning on gcc 2.92.2 compiler, and
44 cleanded up a little bit of code.
47 // The purpose of this class is to permorm the ITS tracking. The
48 // constructor has the task to inizialize some private members. The method
49 // DoTracking is written to be called by a macro. It gets the event number,
50 // the minimum and maximum order number of TPC tracks that are to be tracked
51 // trough the ITS, and the file where the recpoints are registered. The
52 // method Recursivetracking is a recursive function that performs the
53 // tracking trough the ITS The method Intersection found the layer, ladder
54 // and detector whre the intersection take place and caluclate the
55 // cohordinates of this intersection. It returns an integer that is 0 if the
56 // intersection has been found successfully. The two mwthods Kalmanfilter
57 // and kalmanfiltervert operate the kalmanfilter without and with the vertex
58 // imposition respectively. The authors thank Mariana Bondila to have help
59 // them to resolve some problems. July-2000
61 #include <Riostream.h>
62 #include <Riostream.h>
68 #include <TStopwatch.h>
70 #include "TParticle.h"
73 #include "AliITSsegmentationSSD.h"
74 #include "AliITSgeomSPD.h"
75 #include "AliITSgeomSDD.h"
76 #include "AliITSgeomSSD.h"
77 #include "AliITSgeom.h"
78 #include "AliITSRecPoint.h"
80 #include "AliKalmanTrack.h"
82 #include "AliITSTrackV1.h"
83 #include "AliITSIOTrack.h"
84 #include "AliITSRad.h"
85 #include "../TPC/AliTPCtracker.h"
86 #include "AliITSTrackerV1.h"
87 #include "AliITSVertex.h"
89 ClassImp(AliITSTrackerV1)
90 //______________________________________________________________________
91 AliITSTrackerV1::AliITSTrackerV1() {
104 for(ia=0; ia<6; ia++) {
122 //______________________________________________________________________
123 AliITSTrackerV1::AliITSTrackerV1(AliITS* IITTSS, Int_t evnumber, Bool_t flag) {
124 //Origin A. Badala' and G.S. Pappalardo:
125 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
126 // Class constructor. It does some initializations.
128 //PH Initialisation taken from the default constructor
143 Int_t imax = 200,jmax = 450;
144 frl = new AliITSRad(imax,jmax);
146 ////////// gets information on geometry /////////////////////////////
147 AliITSgeom *g1 = fITS->GetITSgeom();
152 for(ia=0; ia<6; ia++) {
153 fNlad[ia]=g1->GetNladders(ia+1);
154 fNdet[ia]=g1->GetNdetectors(ia+1);
155 //cout<<fNlad[i]<<" "<<fNdet[i]<<"\n";
158 //cout<<" mean radius = ";
160 for(ib=0; ib<6; ib++) {
161 g1->GetCenterThetaPhi(ib+1,ll,dd,det);
162 Double_t r1=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
163 g1->GetCenterThetaPhi(ib+1,ll,dd+1,det);
164 Double_t r2=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
165 fAvrad[ib]=(r1+r2)/2.;
166 //cout<<fAvrad[ib]<<" ";
168 //cout<<"\n"; getchar();
170 fDetx[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDx();
171 fDetz[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDz();
173 fDetx[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDx();
174 fDetz[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDz();
176 fDetx[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDx();
177 fDetz[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDz();
179 fDetx[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDx();
180 fDetz[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDz();
182 fDetx[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDx();
183 fDetz[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDz();
185 fDetx[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDx();
186 fDetz[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDz();
187 //cout<<" Detx Detz\n";
188 //for(Int_t la=0; la<6; la++) cout<<" "<<fDetx[la]<<" "<<
192 // allocate memory and define matrices fzmin, fzmax, fphimin and fphimax //
194 Double_t epszdrift=0.05;
196 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
197 Int_t im1, im2, im2max;
198 for(im1=0; im1<6; im1++) {
200 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
203 for(im1=0; im1<6; im1++) {
205 for(im2=0; im2<im2max; im2++) {
206 g1->GetCenterThetaPhi(im1+1,1,im2+1,det);
207 if(im2!=0) fzmin[im1][im2]=det(2)-fDetz[im1];
209 fzmin[im1][im2]=det(2)-(fDetz[im1])*epsz;
210 if(im2!=(im2max-1)) fzmax[im1][im2]=det(2)+fDetz[im1];
212 fzmax[im1][im2]=det(2)+fDetz[im1]*epsz;
213 if(im1==2 || im1==3) {
214 fzmin[im1][im2]-=epszdrift;
215 fzmax[im1][im2]+=epszdrift;
216 } // end if im1==2 || im1==3
220 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
221 for(im1=0;im1<6;im1++) {
223 fphimin[im1] = new Double_t[im2max];
224 fphimax[im1] = new Double_t[im2max];
227 fphidet = new Double_t*[6];
228 for(im1=0; im1<6; im1++) {
230 fphidet[im1] = new Double_t[im2max];
233 //Float_t global[3],local[3];
234 Double_t global[3],local[3];
235 Double_t pigre=TMath::Pi();
236 Double_t xmin,ymin,xmax,ymax;
238 for(im1=0; im1<6; im1++) {
240 for(im2=0; im2<im2max; im2++) {
242 g1->GetCenterThetaPhi(im1+1,im2+1,idet,det);
243 fphidet[im1][im2] = TMath::ATan2(Double_t(det(1)),
245 if(fphidet[im1][im2]<0.) fphidet[im1][im2]+=2.*pigre;
246 local[1]=local[2]=0.;
247 local[0]= -(fDetx[im1]);
248 if(im1==0) local[0]= (fDetx[im1]); //to take into account
249 // different reference system
250 g1->LtoG(im1+1,im2+1,idet,local,global);
251 xmax=global[0]; ymax=global[1];
252 local[0]= (fDetx[im1]);
253 if(im1==0) local[0]= -(fDetx[im1]);//take into account different
255 g1->LtoG(im1+1,im2+1,idet,local,global);
256 xmin=global[0]; ymin=global[1];
257 fphimin[im1][im2]= TMath::ATan2(ymin,xmin);
258 if(fphimin[im1][im2]<0.) fphimin[im1][im2]+=2.*pigre;
259 fphimax[im1][im2]= TMath::ATan2(ymax,xmax);
260 if(fphimax[im1][im2]<0.) fphimax[im1][im2]+=2.*pigre;
263 //////////////////////////////////////////////////////////////////////////////////////////////////////////
264 /////////////// allocate memory and define vector fNRecPoints and matrices fRecCylR, fRecCylPhi, fRecCylZ /////////////
265 gAlice->GetEvent(evnumber);
266 Int_t NumOfModules = g1->GetIndexMax();
267 //fRecCylR = new Float_t *[NumOfModules];
268 fRecCylR = new Double_t *[NumOfModules];
269 //fRecCylPhi = new Float_t *[NumOfModules];
270 fRecCylPhi = new Double_t *[NumOfModules];
271 //fRecCylZ = new Float_t *[NumOfModules];
272 fRecCylZ = new Double_t *[NumOfModules];
273 AliITSRecPoint *recp;
274 fNRecPoints = new Int_t[NumOfModules];
276 for(Int_t module=0; module<NumOfModules; module++) {
277 fITS->ResetRecPoints();
278 gAlice->TreeR()->GetEvent(module);
279 frecPoints=fITS->RecPoints();
280 Int_t nRecPoints=fNRecPoints[module]=frecPoints->GetEntries();
282 fRecCylR[module] = new Float_t[nRecPoints];
283 fRecCylPhi[module] = new Float_t[nRecPoints];
284 fRecCylZ[module] = new Float_t[nRecPoints];
286 fRecCylR[module] = new Double_t[nRecPoints];
287 fRecCylPhi[module] = new Double_t[nRecPoints];
288 fRecCylZ[module] = new Double_t[nRecPoints];
290 for(ind=0; ind<fNRecPoints[module]; ind++) {
291 recp=(AliITSRecPoint*)frecPoints->UncheckedAt(ind);
292 // Float_t global[3], local[3];
293 Double_t global[3], local[3];
294 local[0]=recp->GetX();
296 local[2]= recp->GetZ();
297 g1->LtoG(module,local,global);
299 Float_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]); // r hit
300 Float_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi(); // phi hit
301 Float_t z = global[2]; // z hit
304 Double_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]); // r hit
305 Double_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi(); // phi hit
306 Double_t z = global[2]; // z hit
308 fRecCylR[module][ind]=r;
309 fRecCylPhi[module][ind]=phi;
310 fRecCylZ[module][ind]=z;
315 ////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
318 ////////// gets magnetic field factor //////////////////////////////
320 AliMagF * fieldPointer = gAlice->Field();
321 fFieldFactor = (Double_t)fieldPointer->Factor();
322 //cout<< " field factor = "<<fFieldFactor<<"\n"; getchar();
324 //______________________________________________________________________
325 AliITSTrackerV1::AliITSTrackerV1(const AliITSTrackerV1 &cobj) {
326 // Origin A. Badala' and G.S. Pappalardo:
327 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
331 *fresult = *cobj.fresult;
332 fPtref = cobj.fPtref;
333 fChi2max = cobj.fChi2max;
334 **fvettid = **cobj.fvettid;
335 fflagvert = cobj.fflagvert;
336 Int_t imax=200,jmax=450;
337 frl = new AliITSRad(imax,jmax);
339 fFieldFactor = cobj.fFieldFactor;
340 Int_t i,im1,im2,im2max;
342 fNlad[i] = cobj.fNlad[i];
343 fNdet[i] = cobj.fNdet[i];
344 fAvrad[i] = cobj.fAvrad[i];
345 fDetx[i] = cobj.fDetx[i];
346 fDetz[i] = cobj.fDetz[i];
348 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
349 for(im1=0; im1<6; im1++) {
351 fzmin[im1] = new Double_t[im2max];
352 fzmax[im1] = new Double_t[im2max];
354 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
355 for(im1=0;im1<6;im1++) {
357 fphimin[im1] = new Double_t[im2max];
358 fphimax[im1] = new Double_t[im2max];
361 fphidet = new Double_t*[6];
362 for(im1=0; im1<6; im1++) {
364 fphidet[im1] = new Double_t[im2max];
366 for(im1=0; im1<6; im1++) {
368 for(im2=0; im2<im2max; im2++) {
369 fzmin[im1][im2]=cobj.fzmin[im1][im2];
370 fzmax[im1][im2]=cobj.fzmax[im1][im2];
373 for(im1=0; im1<6; im1++) {
375 for(im2=0; im2<im2max; im2++) {
376 fphimin[im1][im2]=cobj.fphimin[im1][im2];
377 fphimax[im1][im2]=cobj.fphimax[im1][im2];
378 fphidet[im1][im2]=cobj.fphidet[im1][im2];
383 AliITSgeom *g1 = fITS->GetITSgeom();
384 Int_t NumOfModules = g1->GetIndexMax();
386 fRecCylR = new Float_t *[NumOfModules];
387 fRecCylPhi = new Float_t *[NumOfModules];
388 fRecCylZ = new Float_t *[NumOfModules];
390 fRecCylR = new Double_t *[NumOfModules];
391 fRecCylPhi = new Double_t *[NumOfModules];
392 fRecCylZ = new Double_t *[NumOfModules];
393 fNRecPoints = new Int_t[NumOfModules];
394 for(Int_t module=0; module<NumOfModules; module++) {
395 Int_t nRecPoints=fNRecPoints[module]=cobj.fNRecPoints[module];
397 fRecCylR[module] = new Float_t[nRecPoints];
398 fRecCylPhi[module] = new Float_t[nRecPoints];
399 fRecCylZ[module] = new Float_t[nRecPoints];
401 fRecCylR[module] = new Double_t[nRecPoints];
402 fRecCylPhi[module] = new Double_t[nRecPoints];
403 fRecCylZ[module] = new Double_t[nRecPoints];
405 for(ind=0; ind<nRecPoints; ind++) {
406 fRecCylR[module][ind]=cobj.fRecCylR[module][ind];
407 fRecCylPhi[module][ind]=cobj.fRecCylPhi[module][ind];
408 fRecCylZ[module][ind]=cobj.fRecCylZ[module][ind];
413 void AliITSTrackerV1::DelMatrix(Int_t NumOfModules) {
414 for(Int_t mod=0; mod<NumOfModules; mod++) {
415 delete fRecCylR[mod];
416 delete fRecCylPhi[mod];
417 delete fRecCylZ[mod];
423 //______________________________________________________________________
424 AliITSTrackerV1::~AliITSTrackerV1(){
425 // Origin A. Badala' and G.S. Pappalardo:
426 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
430 for(Int_t i=0; i<6; i++) {
445 //______________________________________________________________________
446 AliITSTrackerV1 &AliITSTrackerV1::operator=(AliITSTrackerV1 obj) {
447 // Origin A. Badala' and G.S. Pappalardo:
448 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
449 // assignement operator
452 *fresult = *obj.fresult;
454 fChi2max = obj.fChi2max;
455 **fvettid = **obj.fvettid;
456 fflagvert = obj.fflagvert;
457 Int_t imax=200,jmax=450;
458 frl = new AliITSRad(imax,jmax);
460 fFieldFactor = obj.fFieldFactor;
463 fNlad[i] = obj.fNlad[i];
464 fNdet[i] = obj.fNdet[i];
465 fAvrad[i] = obj.fAvrad[i];
466 fDetx[i] = obj.fDetx[i];
467 fDetz[i] = obj.fDetz[i];
469 fzmin = new Double_t*[6];
470 fzmax = new Double_t*[6];
471 Int_t im1, im2, im2max;
472 for(im1=0; im1<6; im1++) {
474 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
476 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
477 for(im1=0;im1<6;im1++) {
479 fphimin[im1] = new Double_t[im2max];
480 fphimax[im1] = new Double_t[im2max];
483 fphidet = new Double_t*[6];
484 for(im1=0; im1<6; im1++) {
486 fphidet[im1] = new Double_t[im2max];
488 for(im1=0; im1<6; im1++) {
490 for(im2=0; im2<im2max; im2++) {
491 fzmin[im1][im2]=obj.fzmin[im1][im2];
492 fzmax[im1][im2]=obj.fzmax[im1][im2];
495 for(im1=0; im1<6; im1++) {
497 for(im2=0; im2<im2max; im2++) {
498 fphimin[im1][im2]=obj.fphimin[im1][im2];
499 fphimax[im1][im2]=obj.fphimax[im1][im2];
500 fphidet[im1][im2]=obj.fphidet[im1][im2];
504 AliITSgeom *g1 = fITS->GetITSgeom();
505 Int_t NumOfModules = g1->GetIndexMax();
507 fRecCylR = new Float_t *[NumOfModules];
508 fRecCylPhi = new Float_t *[NumOfModules];
509 fRecCylZ = new Float_t *[NumOfModules];
511 fRecCylR = new Double_t *[NumOfModules];
512 fRecCylPhi = new Double_t *[NumOfModules];
513 fRecCylZ = new Double_t *[NumOfModules];
514 fNRecPoints = new Int_t[NumOfModules];
515 for(Int_t module=0; module<NumOfModules; module++) {
516 Int_t nRecPoints=fNRecPoints[module]=obj.fNRecPoints[module];
518 fRecCylR[module] = new Float_t[nRecPoints];
519 fRecCylPhi[module] = new Float_t[nRecPoints];
520 fRecCylZ[module] = new Float_t[nRecPoints];
522 fRecCylR[module] = new Double_t[nRecPoints];
523 fRecCylPhi[module] = new Double_t[nRecPoints];
524 fRecCylZ[module] = new Double_t[nRecPoints];
526 for(ind=0; ind<nRecPoints; ind++) {
527 fRecCylR[module][ind]=obj.fRecCylR[module][ind];
528 fRecCylPhi[module][ind]=obj.fRecCylPhi[module][ind];
529 fRecCylZ[module][ind]=obj.fRecCylZ[module][ind];
536 //______________________________________________________________________
537 void AliITSTrackerV1::DoTracking(Int_t evNumber,Int_t minTr,Int_t maxTr,
538 TFile *file, Bool_t realmass) {
539 // Origin A. Badala' and G.S. Pappalardo:
540 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
541 // The method needs the event number, the minimum and maximum order
542 // number of TPC tracks that
543 // are to be tracked trough the ITS, and the file where the recpoints
545 // The method can be called by a macro. It preforms the tracking for
546 // all good TPC tracks
548 printf("begin DoTracking - file %p\n",file);
550 gAlice->GetEvent(evNumber); //modificato per gestire hbt
552 AliKalmanTrack *kkprov;
553 kkprov->SetConvConst(100/0.299792458/0.2/fFieldFactor);
555 TFile *cf=TFile::Open("AliTPCclusters.root");
556 AliTPCParam *digp= (AliTPCParam*)cf->Get("75x40_100x60_150x60");
557 if (!digp) { cerr<<"TPC parameters have not been found !\n"; getchar();}
560 AliTPCtracker *tracker = new AliTPCtracker(digp,evNumber);
563 tracker->LoadInnerSectors();
564 tracker->LoadOuterSectors();
567 TFile *tf=TFile::Open("AliTPCtracksSorted.root");
569 cerr<<"Can't open AliTPCtracksSorted.root !\n";
572 TObjArray tracks(200000);
574 sprintf(tname,"TreeT_TPC_%d",evNumber);
576 TTree *tracktree=(TTree*)tf->Get(tname);
577 if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n";}
578 TBranch *tbranch=tracktree->GetBranch("tracks");
579 Int_t nentr=(Int_t)tracktree->GetEntries();
582 AliTPCtrack *ioTrackTPC=0;
583 for (kk=0; kk<nentr; kk++) {
584 ioTrackTPC=new AliTPCtrack;
585 tbranch->SetAddress(&ioTrackTPC);
586 tracktree->GetEvent(kk);
587 tracker->CookLabel(ioTrackTPC,0.1);
588 tracks.AddLast(ioTrackTPC);
593 Int_t nt = tracks.GetEntriesFast();
594 cerr<<"Number of found tracks "<<nt<<endl;
597 TTree *tr=gAlice->TreeR();
598 Int_t nent=(Int_t)tr->GetEntries();
599 frecPoints = fITS->RecPoints();
603 Int_t *np = new Int_t[nent];
604 fvettid = new Int_t* [nent];
607 for (mod=0; mod<nent; mod++) {
609 fITS->ResetRecPoints();
610 gAlice->TreeR()->GetEvent(mod);
611 numbpoints = frecPoints->GetEntries();
612 totalpoints+=numbpoints;
613 np[mod] = numbpoints;
614 //cout<<" mod = "<<mod<<" numbpoints = "<<numbpoints<<"\n";getchar();
615 fvettid[mod] = new Int_t[numbpoints];
617 for (ii=0;ii<numbpoints; ii++) *(fvettid[mod]+ii)=0;
620 AliTPCtrack *track=0;
622 if(minTr < 0) {minTr = 0; maxTr = nt-1;}
626 TTree tracktree1("TreeT","Tree with ITS tracks");
627 AliITSIOTrack *ioTrack=0;
628 tracktree1.Branch("ITStracks","AliITSIOTrack",&ioTrack,32000,0);
630 TDatabasePDG * db = new TDatabasePDG;
633 for (j=minTr; j<=maxTr; j++) {
634 track=(AliTPCtrack*)tracks.UncheckedAt(j);
635 if (!track) continue;
637 /// mass definition ////////////////////////
638 Double_t mass=0.13956995;
639 Int_t pcode=211; // a pion by default
641 Int_t TPClabel=TMath::Abs( track->GetLabel() );
642 TParticle *p = (TParticle*)gAlice->Particle(TPClabel);
643 pcode=p->GetPdgCode();
644 // Int_t mothercode=p->GetFirstMother();
645 //if(mothercode>0 ) numofsecondaries++; else numofprimaries++;
647 //if(!pcode) pcode=211;
648 if(TMath::Abs(pcode)<20443) mass=db->GetParticle(pcode)->Mass();
651 ///////////////////////////////////////////////
653 ////// propagation to the end of TPC //////////////
655 track->PropagateTo(xk, 28.94, 1.204e-3,mass); //Ne
657 track->PropagateTo(xk, 44.77, 1.71,mass); //Tedlar
659 track->PropagateTo(xk, 44.86, 1.45,mass); //kevlar
661 track->PropagateTo(xk, 41.28, 0.029,mass); //Nomex
663 track->PropagateTo(xk,36.2,1.98e-3,mass); //C02
665 track->PropagateTo(xk, 24.01, 2.7,mass); //Al
667 track->PropagateTo(xk, 44.77, 1.71,mass); //Tedlar
669 track->PropagateTo(xk, 44.86, 1.45,mass); //kevlar
671 track->PropagateTo(xk, 41.28, 0.029,mass); //Nomex
672 ////////////////////////////////////////////////////////////////////
674 // new propagation to the end of TPC
676 // track->PropagateTo(xk,0.,0.); //Ne if it's still there //attenzione funziona solo se modifica in TPC
677 // Double_t xk=77.415;
678 track->PropagateTo(xk, 28.94, 1.204e-3);
680 track->PropagateTo(xk, 44.77,1.71); //Tedlar
682 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
684 track->PropagateTo(xk, 41.28, 0.029);//Nomex
686 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
688 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
691 // track->PropagateTo(xk,0.,0.); //C02
692 track->PropagateTo(xk,36.2,1.98e-3); //C02 //attenzione funziona solo se modifica in TPC
695 track->PropagateTo(xk, 24.01, 2.7); //Al
697 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
699 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
701 track->PropagateTo(xk, 41.28, 0.029); //Nomex
703 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
705 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
707 track->PropagateTo(xk, 24.01, 2.7); //Al
708 ////////////////////////////////////////////////////////////////////////////////////////////////////////
709 AliITSTrackV1 trackITS(*track);
710 trackITS.PutMass(mass); //new to add mass to track
711 if(fresult){ delete fresult; fresult=0;}
712 fresult = new AliITSTrackV1(trackITS);
714 AliITSTrackV1 primaryTrack(trackITS);
715 vgeant=(*fresult).GetVertex();
717 // Definition of dv and zv for vertex constraint
718 Double_t sigmaDv=0.0050; Double_t sigmaZv=0.010;
719 //Double_t sigmaDv=0.0015; Double_t sigmaZv=0.0015;
720 Double_t uniform= gRandom->Uniform();
722 if(uniform<=0.5) signdv=-1.;
726 Double_t vr=TMath::Sqrt(vgeant(0)*vgeant(0)+ vgeant(1)*vgeant(1));
727 Double_t dv=gRandom->Gaus(signdv*vr,(Float_t)sigmaDv);
728 Double_t zv=gRandom->Gaus(vgeant(2),(Float_t)sigmaZv);
729 //cout<<" Dv e Zv = "<<dv<<" "<<zv<<"\n";
732 trackITS.SetsigmaDv(sigmaDv);
733 trackITS.SetsigmaZv(sigmaZv);
734 (*fresult).SetDv(dv);
735 (*fresult).SetZv(zv);
736 (*fresult).SetsigmaDv(sigmaDv);
737 (*fresult).SetsigmaZv(sigmaZv);
738 primaryTrack.SetDv(dv);
739 primaryTrack.SetZv(zv);
740 primaryTrack.SetsigmaDv(sigmaDv);
741 primaryTrack.SetsigmaZv(sigmaZv);
742 primaryTrack.PrimaryTrack(frl);
743 TVector d2=primaryTrack.Getd2();
744 TVector tgl2=primaryTrack.Gettgl2();
745 TVector dtgl=primaryTrack.Getdtgl();
746 trackITS.Setd2(d2); trackITS.Settgl2(tgl2);
747 trackITS.Setdtgl(dtgl);
748 (*fresult).Setd2(d2); (*fresult).Settgl2(tgl2);
749 (*fresult).Setdtgl(dtgl);
751 trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
752 (*result).SetVertex(vertex); (*result).SetErrorVertex(ervertex);
754 TList *list= new TList();
756 list->AddLast(&trackITS);
758 fPtref=TMath::Abs( (trackITS).GetPt() );
759 if(fPtref>1.0) fChi2max=40.;
760 if(fPtref<=1.0) fChi2max=20.;
761 if(fPtref<0.4 ) fChi2max=100.;
762 if(fPtref<0.2 ) fChi2max=40.;
763 // if(fPtref<0.4 ) fChi2max=30.;
764 // if(fPtref<0.2 ) fChi2max=20.;
765 //if(fPtref<0.2 ) fChi2max=10.;
766 //if(fPtref<0.1 ) fChi2max=5.;
767 //cout << "\n Pt = " << fPtref <<"\n"; //stampa
768 RecursiveTracking(list);
773 TVector vecTotLabRef(18);
775 for(lay=5; lay>=0; lay--) {
776 TVector vecLabRef(3);
777 vecLabRef=(*fresult).GetLabTrack(lay);
778 Float_t clustZ=(*fresult).GetZclusterTrack( lay);
780 Int_t lpp=(Int_t)vecLabRef(k);
782 TParticle *p=(TParticle*) gAlice->Particle(lpp);
783 Int_t pcode=p->GetPdgCode();
784 if(pcode==11) vecLabRef(k)=p->GetFirstMother();
786 itot++; vecTotLabRef(itot)=vecLabRef(k);
787 if(vecLabRef(k)==0. && clustZ == -1.) vecTotLabRef(itot) =-3.;
792 (*fresult).Search(vecTotLabRef, labref, freq);
794 //if(freq < 6) labref=-labref; // cinque - sei
795 if(freq < 5) labref=-labref; // cinque - sei
796 (*fresult).SetLabel(labref);
798 // cout<<" progressive track number = "<<j<<"\r";
800 Int_t numOfCluster=(*fresult).GetNumClust();
801 //cout<<" progressive track number = "<<j<<"\n"; // stampa
802 Long_t labITS=(*fresult).GetLabel();
803 //cout << " ITS track label = " << labITS << "\n"; // stampa
804 Int_t lab=track->GetLabel();
805 //cout << " TPC track label = " << lab <<"\n"; // stampa
806 //propagation to vertex
809 if((*fresult).DoNotCross(rbeam)) continue; //no intersection with beampipe
810 (*fresult).Propagation(rbeam);
811 Double_t c00,c10,c11,c20,c21,c22,c30,c31,c32,c33,c40,c41,c42,c43,c44;
812 (*fresult).GetCElements(c00,
816 c40,c41,c42,c43,c44);
818 Double_t pt=TMath::Abs((*fresult).GetPt());
819 Double_t dr=(*fresult).GetD();
820 Double_t z=(*fresult).GetZ();
821 Double_t tgl=(*fresult).GetTgl();
822 Double_t c=(*fresult).GetC();
824 Double_t dz=z-(tgl/cy)*TMath::ASin((*fresult).Arga(rbeam));
826 // cout<<" dr e dz alla fine = "<<dr<<" "<<dz<<"\n"; getchar();
827 Double_t phi=(*fresult).Getphi();
828 Double_t phivertex = phi - TMath::ASin((*fresult).ArgA(rbeam));
829 Double_t duepi=2.*TMath::Pi();
830 if(phivertex>duepi) phivertex-=duepi;
831 if(phivertex<0.) phivertex+=duepi;
832 /////////////////////////////////////////////////////////////
833 Int_t idmodule,idpoint;
834 if(numOfCluster >=5) { // cinque - sei
835 //if(numOfCluster ==6) { // cinque - sei
836 AliITSIOTrack outTrack;
838 ioTrack->SetStatePhi(phi);
839 ioTrack->SetStateZ(z);
840 ioTrack->SetStateD(dr);
841 ioTrack->SetStateTgl(tgl);
842 ioTrack->SetStateC(c);
843 Double_t radius=(*fresult).Getrtrack();
844 ioTrack->SetRadius(radius);
846 if(c>0.) charge=-1; else charge=1;
847 ioTrack->SetCharge(charge);
848 ioTrack->SetCovMatrix(c00,
852 c40,c41,c42,c43,c44);
853 Double_t px=pt*TMath::Cos(phivertex);
854 Double_t py=pt*TMath::Sin(phivertex);
856 Double_t xtrack=dr*TMath::Sin(phivertex);
857 Double_t ytrack=dr*TMath::Cos(phivertex);
858 Double_t ztrack=dz+vgeant(2);
862 ioTrack->SetX(xtrack);
863 ioTrack->SetY(ytrack);
864 ioTrack->SetZ(ztrack);
865 ioTrack->SetLabel(labITS);
866 ioTrack->SetTPCLabel(lab);
869 for(il=0;il<6; il++){
870 ioTrack->SetIdPoint(il,(*fresult).GetIdPoint(il));
871 ioTrack->SetIdModule(il,(*fresult).GetIdModule(il));
874 for (il=0;il<6;il++) {
875 idpoint=(*fresult).GetIdPoint(il);
876 idmodule=(*fresult).GetIdModule(il);
877 if(idmodule>0.) *(fvettid[idmodule]+idpoint)=1;
879 ioTrack->SetIdPoint(il,idpoint);
880 ioTrack->SetIdModule(il,idmodule);
882 } // end if on numOfCluster
883 //gObjectTable->Print(); // stampa memoria
884 } // end for (int j=minTr; j<=maxTr; j++)
886 static Bool_t first=kTRUE;
889 tfile=new TFile("itstracks.root","RECREATE");
890 //cout<<"I have opened itstracks.root file "<<endl;
896 sprintf(hname,"TreeT%d",evNumber);
897 tracktree1.Write(hname);
899 TTree *fAli=gAlice->TreeK();
901 if (fAli) fileAli =fAli->GetCurrentFile();
903 ////////////////////////////////////////////////////////////////////
905 printf("delete vectors\n");
907 if(fvettid) delete [] fvettid;
908 if(fresult) {delete fresult; fresult=0;}
910 //______________________________________________________________________
911 void AliITSTrackerV1::RecursiveTracking(TList *trackITSlist) {
912 // This function perform the recursive tracking in ITS detectors
913 // reference is a pointer to the final best track
914 // Origin A. Badala' and G.S. Pappalardo:
915 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
916 // The authors thank Mariana Bondila to have help them to resolve some
917 // problems. July-2000
919 //Rlayer[0]=4.; Rlayer[1]=7.; Rlayer[2]=14.9;
920 // Rlayer[3]=23.8; Rlayer[4]=39.1; Rlayer[5]=43.6; //vecchio
922 //////////////////////
923 Float_t sigmaphil[6], sigmazl[6];
924 sigmaphil[0]=1.44e-6/(fAvrad[0]*fAvrad[0]);
925 sigmaphil[1]=1.44e-6/(fAvrad[1]*fAvrad[1]);
926 sigmaphil[2]=1.444e-5/(fAvrad[2]*fAvrad[2]);
927 sigmaphil[3]=1.444e-5/(fAvrad[3]*fAvrad[3]);
928 sigmaphil[4]=4e-6/(fAvrad[4]*fAvrad[4]);
929 sigmaphil[5]=4e-6/(fAvrad[5]*fAvrad[5]);
936 ///////////////////////////////////////////////////////////
938 AliITSgeom *g1 = fITS->GetITSgeom();
939 AliITSRecPoint *recp;
940 for(index =0; index<trackITSlist->GetSize(); index++) {
941 AliITSTrackV1 *trackITS = (AliITSTrackV1 *) trackITSlist->At(index);
942 if((*trackITS).GetLayer()==7) fresult->SetChi2(10.223e140);
943 // cout <<" Layer inizio = "<<(*trackITS).GetLayer()<<"\n";
944 // cout<<"fvtrack =" <<"\n";
945 // cout << (*trackITS)(0) << " "<<(*trackITS)(1)<<" "
946 // <<(*trackITS)(2)<<" "<<(*trackITS)(3)<<" "
947 // <<(*trackITS)(4)<<"\n";
948 // cout<< " rtrack = "<<(*trackITS).Getrtrack()<<"\n";
949 // cout<< " Pt = "<<(*trackITS).GetPt()<<"\n";
951 Double_t chi2Now, chi2Ref;
952 Float_t numClustRef = fresult->GetNumClust();
953 if((*trackITS).GetLayer()==1 ) {
954 chi2Now = trackITS->GetChi2();
955 Float_t numClustNow = trackITS->GetNumClust();
956 if(trackITS->GetNumClust())
957 chi2Now /= (Double_t)trackITS->GetNumClust();
958 chi2Ref = fresult->GetChi2();
959 if(fresult->GetNumClust())
960 chi2Ref /= (Double_t)fresult->GetNumClust();
961 //cout<<" chi2Now and chi2Ref = "<<chi2Now<<" "<<chi2Ref<<"\n";
962 if( numClustNow > numClustRef ) {*fresult = *trackITS;}
963 if((numClustNow == numClustRef )&&
964 (chi2Now < chi2Ref)) {
965 *fresult = *trackITS;
970 if(trackITS->Getfnoclust()>=2) continue;
971 Float_t numClustNow = trackITS->GetNumClust();
973 chi2Now = trackITS->GetChi2();
975 if(numClustNow<numClustRef && chi2Now>fresult->GetChi2()) continue;
976 //cout<<" chi2Now = "<<chi2Now<<"\n";
978 chi2Now/=numClustNow;
979 if(fPtref > 1.0 && chi2Now > 30.) continue;
980 if((fPtref >= 0.6 && fPtref<=1.0) && chi2Now > 40.) continue;
981 // if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 40.) continue;
982 // if(fPtref <= 0.2 && chi2Now > 8.) continue;
983 if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 30.) continue;
984 if(fPtref <= 0.2 && chi2Now > 7.) continue;
985 /////////////////////////////
988 Int_t layerInit = (*trackITS).GetLayer();
989 Int_t layernew = layerInit - 2;// -1 for new layer, -1 for matrix index
991 Int_t ladp, ladm, detp,detm,ladinters,detinters;
992 Int_t layerfin=layerInit-1;
993 // cout<<"Prima di intersection \n";
994 Int_t outinters=Intersection(*trackITS,layerfin,ladinters,detinters);
995 // cout<<" outinters = "<<outinters<<"\n";
996 // cout<<" Layer ladder detector intersection ="
997 // <<layerfin<<" "<<ladinters<<" "<<detinters<<"\n";
998 // cout << " phiinters zinters = "<<(*trackITS)(0)
999 // << " "<<(*trackITS)(1)<<"\n"; getchar();
1000 if(outinters==-1) continue;
1003 TVector toucLad(9), toucDet(9);
1004 Int_t lycur=layerfin;
1007 if(ladm <= 0) ladm=fNlad[layerfin-1];
1008 if(ladp > fNlad[layerfin-1]) ladp=1;
1013 toucLad(0)=ladinters; toucLad(1)=ladm; toucLad(2)=ladp;
1014 toucLad(3)=ladinters; toucLad(4)=ladm; toucLad(5)=ladp;
1015 toucLad(6)=ladinters; toucLad(7)=ladm; toucLad(8)=ladp;
1016 toucDet(0)=detinters; toucDet(1)=detinters; toucDet(2)=detinters;
1017 if(detm > 0 && detp <= fNdet[layerfin-1]) {
1019 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1020 toucDet(6)=detp; toucDet(7)=detp; toucDet(8)=detp;
1022 if(detm > 0 && detp > fNdet[layerfin-1]) {
1024 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1026 if(detm <= 0 && detp <= fNdet[layerfin-1]) {
1028 toucDet(3)=detp; toucDet(4)=detp; toucDet(5)=detp;
1031 Float_t epsphi=5.0, epsz=5.0;
1032 if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1033 // new definition of idetot e toucLad e toucDet to be
1034 // transformed in a method
1035 // these values could be modified
1036 Float_t pigre=TMath::Pi();
1037 Float_t rangephi=5., rangez=5.;
1038 if(layerfin==1 || layerfin ==2){
1039 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1040 (*trackITS).GetSigmaphi());
1041 rangez = 40.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1042 (*trackITS).GetSigmaZ());
1044 if(layerfin==3 || layerfin ==4){
1045 //rangephi=30.*fepsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1046 // (*trackITS).GetSigmaphi());
1047 //rangez = 40.*fepsz*TMath::Sqrt(sigmazl[layerfin-1]+
1048 // (*trackITS).GetSigmaZ());
1049 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1050 (*trackITS).GetSigmaphi());
1051 rangez = 50.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1052 (*trackITS).GetSigmaZ());
1054 if(layerfin==5 || layerfin ==6){
1055 rangephi=20.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1056 (*trackITS).GetSigmaphi());
1057 rangez =5.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1058 (*trackITS).GetSigmaZ());
1060 Float_t phinters, zinters;
1061 phinters=(*trackITS).Getphi();
1062 zinters=(*trackITS).GetZ();
1063 Float_t distz = 0.0;
1064 Float_t phicm, phicp, distphim, distphip;
1066 if(phinters>fphimax[layerfin-1][ladm-1]) phicm=phinters-2*pigre; //corretto il 20-11-2001
1067 distphim=TMath::Abs(phicm-fphimax[layerfin-1][ladm-1]); //corretto il 20-11-2001
1069 //cout<<" fNlad[layerfin-1] e ladp = "<<fNlad[layerfin-1]<<" "<<ladp<<endl;
1070 if(phinters>fphimin[layerfin-1][ladp-1]) phicp=phinters-2.*pigre; //corretto il 20-11-2001
1071 distphip=TMath::Abs(phicp-fphimin[layerfin-1][ladp-1]); //corretto il 20-11-2001
1075 toucLad(0)=ladinters; toucDet(0)=detinters;
1076 if(detm>0) distz=TMath::Abs(zinters-fzmax[layerfin-1][detm-1]);
1077 if(detm>0 && rangez>=distz){
1079 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detm;
1080 if(rangephi>=distphim){
1082 toucLad(idetot-1)=ladm;
1083 toucDet(idetot-1)=detinters;
1085 toucLad(idetot-1)=ladm;
1086 toucDet(idetot-1)=detm;
1088 if(rangephi>=distphip){
1090 toucLad(idetot-1)=ladp;
1091 toucDet(idetot-1)=detinters;
1093 toucLad(idetot-1)=ladp;
1094 toucDet(idetot-1)=detm;
1097 if(detp<=fNdet[layerfin-1])
1098 distz=TMath::Abs(zinters-fzmin[layerfin-1][detp-1]);
1099 if(detp<=fNdet[layerfin-1] && rangez>=distz){
1101 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detp;
1102 if(rangephi>=distphim){
1103 idetot++; toucLad(idetot-1)=ladm; toucDet(idetot-1)=detp;
1106 toucLad(idetot-1)=ladm;
1107 toucDet(idetot-1)=detinters;
1110 if(rangephi>=distphip){
1112 toucLad(idetot-1)=ladp;
1113 toucDet(idetot-1)=detp;
1116 toucLad(idetot-1)=ladp;
1117 toucDet(idetot-1)=detinters;
1120 } //end detm<fNdet[.......
1122 if(flagzmin == 0 && flagzmax==0){
1123 if(rangephi>=distphim){
1125 toucLad(idetot-1)=ladm;
1126 toucDet(idetot-1)=detinters;
1128 if(rangephi>=distphip){
1130 toucLad(idetot-1)=ladp;
1131 toucDet(idetot-1)=detinters;
1134 ////////////////////////////////////////////////////////////
1136 for (iriv=0; iriv<idetot; iriv++) { //for on detectors
1137 ///////////////////////////////////////////////////////
1138 /*** Rec points sorted by module *****/
1139 /**************************************/
1141 indexmod = g1->GetModuleIndex(lycur,(Int_t)toucLad(iriv),
1142 (Int_t)toucDet(iriv));
1143 fITS->ResetRecPoints();
1144 gAlice->TreeR()->GetEvent(indexmod);
1145 Int_t npoints=frecPoints->GetEntries();
1148 for(indnew=0; indnew<npoints; indnew++){
1149 if (*(fvettid[indexmod]+indnew)==0)
1150 recp =(AliITSRecPoint*)frecPoints->UncheckedAt(indnew);
1153 TVector cluster(3),vecclust(9);
1154 //vecclust(6)=vecclust(7)=vecclust(8)=-1.;
1156 // now vecclust is with cylindrical cohordinates
1157 vecclust(0)=(Float_t)fRecCylR[indexmod][indnew];
1158 vecclust(1)=(Float_t)fRecCylPhi[indexmod][indnew];
1159 vecclust(2)=(Float_t)fRecCylZ[indexmod][indnew];
1160 vecclust(3) = (Float_t)recp->fTracks[0];
1161 vecclust(4) = (Float_t)indnew;
1162 vecclust(5) = (Float_t)indexmod;
1163 vecclust(6) = (Float_t)recp->fTracks[0];
1164 vecclust(7) = (Float_t)recp->fTracks[1];
1165 vecclust(8) = (Float_t)recp->fTracks[2];
1166 sigma[0] = (Double_t) recp->GetSigmaX2();
1167 sigma[1] = (Double_t) recp->GetSigmaZ2();
1169 cluster(0)=fRecCylR[indexmod][indnew];
1170 cluster(1)=fRecCylPhi[indexmod][indnew];
1171 cluster(2)=fRecCylZ[indexmod][indnew];
1173 // cout<<" layer = "<<play<<"\n";
1174 // cout<<" cluster prima = "<<vecclust(0)<<" "
1175 // <<vecclust(1)<<" "
1176 // <<vecclust(2)<<"\n"; getchar();
1178 Float_t sigmatotphi, sigmatotz;
1179 // Float_t epsphi=5.0, epsz=5.0;
1180 //if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1181 Double_t rTrack=(*trackITS).Getrtrack();
1182 Double_t sigmaphi=sigma[0]/(rTrack*rTrack);
1183 sigmatotphi=epsphi*TMath::Sqrt(sigmaphi +
1184 (*trackITS).GetSigmaphi());
1185 sigmatotz=epsz*TMath::Sqrt(sigma[1] +
1186 (*trackITS).GetSigmaZ());
1187 //cout<<"cluster e sigmatotphi e track = "<<cluster(0)
1188 // <<" "<<cluster(1)<<" "<<sigmatotphi<<" "
1189 // <<vecclust(3)<<"\n";
1190 //if(vecclust(3)==481) getchar();
1191 if(cluster(1)<6. && (*trackITS).Getphi()>6.)
1192 cluster(1)=cluster(1)+(2.*TMath::Pi());
1193 if(cluster(1)>6. && (*trackITS).Getphi()<6.)
1194 cluster(1)=cluster(1)-(2.*TMath::Pi());
1195 if(TMath::Abs(cluster(1)-(*trackITS).Getphi())>sigmatotphi)
1197 // cout<<" supero sigmaphi \n";
1198 AliITSTrackV1 *newTrack = new AliITSTrackV1((*trackITS));
1199 (*newTrack).SetLayer((*trackITS).GetLayer()-1);
1200 if (TMath::Abs(rTrack-cluster(0))/rTrack>1e-6)
1201 (*newTrack).Correct(Double_t(cluster(0)));
1202 //cout<<" cluster(2) e(*newTrack).GetZ()="<<cluster(2)<<" "
1203 // << (*newTrack).GetZ()<<"\n";
1204 if(TMath::Abs(cluster(2)-(*newTrack).GetZ()) > sigmatotz){
1208 Double_t sigmanew[2];
1209 sigmanew[0]= sigmaphi;
1210 sigmanew[1]=sigma[1];
1214 // Double_t chi2pred=newTrack->GetPredChi2(m,sigmanew);
1215 // cout<<" chi2pred = "<<chi2pred<<"\n";
1216 // if(chi2pred>fChi2max) continue; //aggiunto il 30-7-2001
1217 if(iriv == 0) flaghit=1;
1218 (*newTrack).AddMS(frl); // add the multiple scattering
1219 //matrix to the covariance matrix
1220 (*newTrack).AddEL(frl,1.,0);
1223 KalmanFilterVert(newTrack,cluster,sigmanew);
1224 //KalmanFilterVert(newTrack,cluster,sigmanew,chi2pred);
1226 KalmanFilter(newTrack,cluster,sigmanew);
1228 (*newTrack).PutCluster(layernew, vecclust);
1229 newTrack->AddClustInTrack();
1230 listoftrack.AddLast(newTrack);
1232 } // end of for on detectors (iriv)
1233 }//end if(outinters==0)
1235 if(flaghit==0 || outinters==-2) {
1236 AliITSTrackV1 *newTrack = new AliITSTrackV1(*trackITS);
1237 (*newTrack).Setfnoclust();
1238 (*newTrack).SetLayer((*trackITS).GetLayer()-1);
1239 (*newTrack).AddMS(frl); // add the multiple scattering matrix
1240 // to the covariance matrix
1241 (*newTrack).AddEL(frl,1.,0);
1242 listoftrack.AddLast(newTrack);
1245 //gObjectTable->Print(); // stampa memoria
1247 RecursiveTracking(&listoftrack);
1248 listoftrack.Delete();
1249 } // end of for on tracks (index)
1251 //gObjectTable->Print(); // stampa memoria
1253 //______________________________________________________________________
1254 Int_t AliITSTrackerV1::Intersection(AliITSTrackV1 &track,Int_t layer,
1255 Int_t &ladder,Int_t &detector) {
1256 // Origin A. Badala' and G.S. Pappalardo
1257 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1258 // Found the intersection and the detector
1260 Double_t rk=fAvrad[layer-1];
1261 if(track.DoNotCross(rk)){ /*cout<< " Do not cross \n";*/ return -1;}
1262 track.Propagation(rk);
1263 Double_t zinters=track.GetZ();
1264 Double_t phinters=track.Getphi();
1265 //cout<<"zinters = "<<zinters<<" phinters = "<<phinters<<"\n";
1269 TVector distZCenter(2);
1273 for(iD = 1; iD<= fNdet[layer-1]; iD++) {
1274 if(zinters > fzmin[layer-1][iD-1] && zinters <= fzmax[layer-1][iD-1]) {
1276 cout<< " Errore su iz in Intersection \n";
1279 listDet(iz)= iD; distZCenter(iz)=TMath::Abs(zinters-det(2));
1285 if(iz==0) {/* cout<< " No detector along Z \n";*/ return -2;}
1286 detector=Int_t (listDet(0));
1287 if(iz>1 && (distZCenter(0)>distZCenter(1))) detector=Int_t (listDet(1));
1290 TVector distPhiCenter(2);
1292 Double_t pigre=TMath::Pi();
1294 for(iLd = 1; iLd<= fNlad[layer-1]; iLd++) {
1295 Double_t phimin=fphimin[layer-1][iLd-1];
1296 Double_t phimax=fphimax[layer-1][iLd-1];
1297 Double_t phidet=fphidet[layer-1][iLd-1];
1298 Double_t phiconfr=phinters;
1300 //if(phimin <5.5) {cout<<" Error in Intersection for phi \n";
1303 if(phinters>(1.5*pigre)) phiconfr=phinters-(2.*pigre);
1304 if(phidet>(1.5*pigre)) phidet-=(2.*pigre);
1306 if(phiconfr>phimin && phiconfr<= phimax) {
1308 cout<< " Errore su ip in Intersection \n"; getchar();
1311 distPhiCenter(ip)=TMath::Abs(phiconfr-phidet); ip++;
1315 if(ip==0) { cout<< " No detector along phi \n"; getchar();}
1316 ladder=Int_t (listLad(0));
1317 if(ip>1 && (distPhiCenter(0)>distPhiCenter(1))) ladder=Int_t (listLad(1));
1320 //______________________________________________________________________
1321 void AliITSTrackerV1::KalmanFilter(AliITSTrackV1 *newTrack,TVector &cluster,
1323 //Origin A. Badala' and G.S. Pappalardo:
1324 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1325 // Kalman filter without vertex constraint
1326 ////// Evaluation of the measurement vector ////////////////////////
1328 Double_t rk,phik,zk;
1329 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1331 //////////////////////// Evaluation of the error matrix V /////////
1332 Double_t v00=sigma[0];
1333 Double_t v11=sigma[1];
1334 ////////////////////////////////////////////////////////////////////
1335 Double_t cin00,cin10,cin20,cin30,cin40,cin11,cin21,cin31,cin41,cin22,
1336 cin32,cin42,cin33,cin43,cin44;
1338 newTrack->GetCElements(cin00,
1341 cin30,cin31,cin32,cin33,
1342 cin40,cin41,cin42,cin43,cin44); //get C matrix
1343 Double_t rold00=cin00+v00;
1344 Double_t rold10=cin10;
1345 Double_t rold11=cin11+v11;
1346 ////////////////////// R matrix inversion /////////////////////////
1347 Double_t det=rold00*rold11-rold10*rold10;
1348 Double_t r00=rold11/det;
1349 Double_t r10=-rold10/det;
1350 Double_t r11=rold00/det;
1351 ////////////////////////////////////////////////////////////////////
1352 Double_t k00=cin00*r00+cin10*r10;
1353 Double_t k01=cin00*r10+cin10*r11;
1354 Double_t k10=cin10*r00+cin11*r10;
1355 Double_t k11=cin10*r10+cin11*r11;
1356 Double_t k20=cin20*r00+cin21*r10;
1357 Double_t k21=cin20*r10+cin21*r11;
1358 Double_t k30=cin30*r00+cin31*r10;
1359 Double_t k31=cin30*r10+cin31*r11;
1360 Double_t k40=cin40*r00+cin41*r10;
1361 Double_t k41=cin40*r10+cin41*r11;
1362 Double_t x0,x1,x2,x3,x4;
1363 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1364 Double_t savex0=x0, savex1=x1;
1365 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1);
1366 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1);
1367 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1);
1368 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1);
1369 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1);
1370 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1371 c00=cin00-k00*cin00-k01*cin10;
1372 c10=cin10-k00*cin10-k01*cin11;
1373 c20=cin20-k00*cin20-k01*cin21;
1374 c30=cin30-k00*cin30-k01*cin31;
1375 c40=cin40-k00*cin40-k01*cin41;
1376 c11=cin11-k10*cin10-k11*cin11;
1377 c21=cin21-k10*cin20-k11*cin21;
1378 c31=cin31-k10*cin30-k11*cin31;
1379 c41=cin41-k10*cin40-k11*cin41;
1380 c22=cin22-k20*cin20-k21*cin21;
1381 c32=cin32-k20*cin30-k21*cin31;
1382 c42=cin42-k20*cin40-k21*cin41;
1383 c33=cin33-k30*cin30-k31*cin31;
1384 c43=cin43-k30*cin40-k31*cin41;
1385 c44=cin44-k40*cin40-k41*cin41;
1386 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1387 newTrack->PutCElements(c00,
1391 c40,c41,c42,c43,c44); // put in track the
1393 Double_t vmcold00=v00-c00;
1394 Double_t vmcold10=-c10;
1395 Double_t vmcold11=v11-c11;
1396 ////////////////////// Matrix vmc inversion ///////////////////////
1397 det=vmcold00*vmcold11-vmcold10*vmcold10;
1398 Double_t vmc00=vmcold11/det;
1399 Double_t vmc10=-vmcold10/det;
1400 Double_t vmc11=vmcold00/det;
1401 ////////////////////////////////////////////////////////////////////
1402 Double_t chi2=(m[0]-x0)*( vmc00*(m[0]-x0) + 2.*vmc10*(m[1]-x1) ) +
1403 (m[1]-x1)*vmc11*(m[1]-x1);
1404 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1406 //----------------------------------------------------------------------
1407 //void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1408 // TVector &cluster,Double_t sigma[2]){
1409 void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1410 TVector &cluster,Double_t sigma[2]
1411 /*, Double_t chi2pred*/){
1412 //Origin A. Badala' and G.S. Pappalardo:
1413 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1414 // Kalman filter with vertex constraint
1415 ///////////////////// Evaluation of the measurement vector m ////////
1417 Double_t rk,phik,zk;
1418 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1420 Double_t cc=(*newTrack).GetC();
1421 Double_t zv=(*newTrack).GetZv();
1422 Double_t dv=(*newTrack).GetDv();
1424 Double_t tgl= (zk-zv)*cy/TMath::ASin(cy*rk);
1426 /////////////////////// Evaluation of the error matrix V //////////
1427 Int_t layer=newTrack->GetLayer();
1428 Double_t v00=sigma[0];
1429 Double_t v11=sigma[1];
1430 Double_t v31=sigma[1]/rk;
1431 Double_t sigmaDv=newTrack->GetsigmaDv();
1432 Double_t v22=sigmaDv*sigmaDv + newTrack->Getd2(layer-1);
1433 Double_t v32=newTrack->Getdtgl(layer-1);
1434 Double_t sigmaZv=newTrack->GetsigmaZv();
1435 Double_t v33=(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+newTrack->Gettgl2(layer-1);
1436 //////////////////////////////////////////////////////////////////
1437 Double_t cin00,cin10,cin11,cin20,cin21,cin22,
1438 cin30,cin31,cin32,cin33,cin40,cin41,cin42,cin43,cin44;
1439 newTrack->GetCElements(cin00,
1442 cin30,cin31,cin32,cin33,
1443 cin40,cin41,cin42,cin43,cin44); //get C matrix
1451 r[3][1]=cin31+sigma[1]/rk;
1452 r[2][2]=cin22+sigmaDv*sigmaDv+newTrack->Getd2(layer-1);
1453 r[3][2]=cin32+newTrack->Getdtgl(layer-1);
1454 r[3][3]=cin33+(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+
1455 newTrack->Gettgl2(layer-1);
1456 r[0][1]=r[1][0]; r[0][2]=r[2][0]; r[0][3]=r[3][0];
1457 r[1][2]=r[2][1]; r[1][3]=r[3][1]; r[2][3]=r[3][2];
1458 ///////////////////// Matrix R inversion //////////////////////////
1462 Int_t ll[kn],mm[kn];
1465 for(k=0; k<kn; k++) {
1469 for(j=k; j<kn ; j++) {
1470 for (i=j; i<kn; i++) {
1471 if(TMath::Abs(big) < TMath::Abs(r[i][j]) ) {
1481 for(i=0; i<kn; i++) {
1490 for(j=0; j<kn; j++) {
1499 cout << "Singular matrix\n";
1501 for(i=0; i<kn; i++) {
1502 if(i == k) { continue; }
1503 r[i][k]=r[i][k]/(-big);
1506 for(i=0; i<kn; i++) {
1508 for(j=0; j<kn; j++) {
1509 if(i == k || j == k) continue;
1510 r[i][j]=hold*r[k][j]+r[i][j];
1514 for(j=0; j<kn; j++) {
1515 if(j == k) continue;
1516 r[k][j]=r[k][j]/big;
1524 for(k=kn-1; k>=0; k--) {
1527 for (j=0; j<kn; j++) {
1535 for (i=0; i<kn; i++) {
1542 ////////////////////////////////////////////////////////////////////
1543 Double_t k00=cin00*r[0][0]+cin10*r[1][0]+cin20*r[2][0]+cin30*r[3][0];
1544 Double_t k01=cin00*r[1][0]+cin10*r[1][1]+cin20*r[2][1]+cin30*r[3][1];
1545 Double_t k02=cin00*r[2][0]+cin10*r[2][1]+cin20*r[2][2]+cin30*r[3][2];
1546 Double_t k03=cin00*r[3][0]+cin10*r[3][1]+cin20*r[3][2]+cin30*r[3][3];
1547 Double_t k10=cin10*r[0][0]+cin11*r[1][0]+cin21*r[2][0]+cin31*r[3][0];
1548 Double_t k11=cin10*r[1][0]+cin11*r[1][1]+cin21*r[2][1]+cin31*r[3][1];
1549 Double_t k12=cin10*r[2][0]+cin11*r[2][1]+cin21*r[2][2]+cin31*r[3][2];
1550 Double_t k13=cin10*r[3][0]+cin11*r[3][1]+cin21*r[3][2]+cin31*r[3][3];
1551 Double_t k20=cin20*r[0][0]+cin21*r[1][0]+cin22*r[2][0]+cin32*r[3][0];
1552 Double_t k21=cin20*r[1][0]+cin21*r[1][1]+cin22*r[2][1]+cin32*r[3][1];
1553 Double_t k22=cin20*r[2][0]+cin21*r[2][1]+cin22*r[2][2]+cin32*r[3][2];
1554 Double_t k23=cin20*r[3][0]+cin21*r[3][1]+cin22*r[3][2]+cin32*r[3][3];
1555 Double_t k30=cin30*r[0][0]+cin31*r[1][0]+cin32*r[2][0]+cin33*r[3][0];
1556 Double_t k31=cin30*r[1][0]+cin31*r[1][1]+cin32*r[2][1]+cin33*r[3][1];
1557 Double_t k32=cin30*r[2][0]+cin31*r[2][1]+cin32*r[2][2]+cin33*r[3][2];
1558 Double_t k33=cin30*r[3][0]+cin31*r[3][1]+cin32*r[3][2]+cin33*r[3][3];
1559 Double_t k40=cin40*r[0][0]+cin41*r[1][0]+cin42*r[2][0]+cin43*r[3][0];
1560 Double_t k41=cin40*r[1][0]+cin41*r[1][1]+cin42*r[2][1]+cin43*r[3][1];
1561 Double_t k42=cin40*r[2][0]+cin41*r[2][1]+cin42*r[2][2]+cin43*r[3][2];
1562 Double_t k43=cin40*r[3][0]+cin41*r[3][1]+cin42*r[3][2]+cin43*r[3][3];
1564 Double_t x0,x1,x2,x3,x4;
1565 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1566 Double_t savex0=x0, savex1=x1, savex2=x2, savex3=x3;
1567 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1)+k02*(m[2]-savex2)+
1569 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1)+k12*(m[2]-savex2)+
1571 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1)+k22*(m[2]-savex2)+
1573 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1)+k32*(m[2]-savex2)+
1575 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1)+k42*(m[2]-savex2)+
1577 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1578 c00=cin00-k00*cin00-k01*cin10-k02*cin20-k03*cin30;
1579 c10=cin10-k00*cin10-k01*cin11-k02*cin21-k03*cin31;
1580 c20=cin20-k00*cin20-k01*cin21-k02*cin22-k03*cin32;
1581 c30=cin30-k00*cin30-k01*cin31-k02*cin32-k03*cin33;
1582 c40=cin40-k00*cin40-k01*cin41-k02*cin42-k03*cin43;
1583 c11=cin11-k10*cin10-k11*cin11-k12*cin21-k13*cin31;
1584 c21=cin21-k10*cin20-k11*cin21-k12*cin22-k13*cin32;
1585 c31=cin31-k10*cin30-k11*cin31-k12*cin32-k13*cin33;
1586 c41=cin41-k10*cin40-k11*cin41-k12*cin42-k13*cin43;
1587 c22=cin22-k20*cin20-k21*cin21-k22*cin22-k23*cin32;
1588 c32=cin32-k20*cin30-k21*cin31-k22*cin32-k23*cin33;
1589 c42=cin42-k20*cin40-k21*cin41-k22*cin42-k23*cin43;
1590 c33=cin33-k30*cin30-k31*cin31-k32*cin32-k33*cin33;
1591 c43=cin43-k30*cin40-k31*cin41-k32*cin42-k33*cin43;
1592 c44=cin44-k40*cin40-k41*cin41-k42*cin42-k43*cin43;
1594 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1595 newTrack->PutCElements(c00,
1599 c40,c41,c42,c43,c44); // put in track the
1602 vmc[0][0]=v00-c00; vmc[1][0]=-c10; vmc[2][0]=-c20; vmc[3][0]=-c30;
1603 vmc[1][1]=v11-c11; vmc[2][1]=-c21; vmc[3][1]=v31-c31;
1604 vmc[2][2]=v22-c22; vmc[3][2]=v32-c32;
1606 vmc[0][1]=vmc[1][0]; vmc[0][2]=vmc[2][0]; vmc[0][3]=vmc[3][0];
1607 vmc[1][2]=vmc[2][1]; vmc[1][3]=vmc[3][1];
1608 vmc[2][3]=vmc[3][2];
1609 /////////////////////// vmc matrix inversion ///////////////////////
1611 for(k=0; k<kn; k++) {
1615 for(j=k; j<kn ; j++) {
1616 for (i=j; i<kn; i++) {
1617 if(TMath::Abs(big) < TMath::Abs(vmc[i][j]) ) {
1627 for(i=0; i<kn; i++) {
1629 vmc[k][i]=vmc[j][i];
1636 for(j=0; j<kn; j++) {
1638 vmc[j][k]=vmc[j][i];
1645 cout << "Singular matrix\n";
1647 for(i=0; i<kn; i++) {
1648 if(i == k) continue;
1649 vmc[i][k]=vmc[i][k]/(-big);
1652 for(i=0; i<kn; i++) {
1654 for(j=0; j<kn; j++) {
1655 if(i == k || j == k) continue;
1656 vmc[i][j]=hold*vmc[k][j]+vmc[i][j];
1660 for(j=0; j<kn; j++) {
1661 if(j == k) continue;
1662 vmc[k][j]=vmc[k][j]/big;
1670 for(k=kn-1; k>=0; k--) {
1673 for (j=0; j<kn; j++) {
1675 vmc[j][k]=-vmc[j][i];
1681 for (i=0; i<kn; i++) {
1683 vmc[k][i]=-vmc[j][i];
1688 ////////////////////////////////////////////////////////////////////
1689 Double_t chi2=(m[0]-x0)*( vmc[0][0]*(m[0]-x0) + 2.*vmc[1][0]*(m[1]-x1) +
1690 2.*vmc[2][0]*(m[2]-x2)+ 2.*vmc[3][0]*(m[3]-x3) )+
1691 (m[1]-x1)* ( vmc[1][1]*(m[1]-x1) + 2.*vmc[2][1]*(m[2]-x2)+
1692 2.*vmc[3][1]*(m[3]-x3) ) +
1693 (m[2]-x2)* ( vmc[2][2]*(m[2]-x2)+ 2.*vmc[3][2]*(m[3]-x3) ) +
1694 (m[3]-x3)*vmc[3][3]*(m[3]-x3);
1695 //cout<<" chi2 kalman = "<<chi2<<"\n"; getchar();
1696 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1697 // newTrack->SetChi2(newTrack->GetChi2()+chi2pred);