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.25 2002/10/24 17:12:58 barbera
19 ITS tracking V1 integrated with the last version of ITS PID
21 Revision 1.24 2002/10/23 14:28:38 barbera
22 Fixes added to get into account the new magnetic field conversion factor automatically
24 Revision 1.23 2002/10/22 18:29:34 barbera
25 Tracking V1 ported to the HEAD
27 Revision 1.22 2002/10/22 14:45:36 alibrary
28 Introducing Riostream.h
30 Revision 1.21 2002/02/05 09:12:26 hristov
31 Small mods for gcc 3.02
33 Revision 1.20 2001/11/21 14:47:45 barbera
34 Some unuseful print-out commented out
36 Revision 1.19 2001/11/21 10:49:07 barbera
37 Bug correction suggested by Rene done
39 Revision 1.18 2001/11/20 15:46:17 barbera
40 Point coordinated are calculated in cylindrical reference frame once and for all at the beginning of tracking V1
42 Revision 1.10.2.1 2001/10/24 07:26:04 hristov
43 All the changes from the head are merged with the release
45 Revision 1.14 2001/10/24 07:19:57 hristov
46 Some pointer correctly initialised in one of the constructors
48 Revision 1.13 2001/10/21 19:17:12 hristov
49 Several pointers were set to zero in the default constructors to avoid memory management problems
51 Revision 1.12 2001/10/19 21:32:35 nilsen
52 Minor changes to remove compliation warning on gcc 2.92.2 compiler, and
53 cleanded up a little bit of code.
56 // The purpose of this class is to permorm the ITS tracking. The
57 // constructor has the task to inizialize some private members. The method
58 // DoTracking is written to be called by a macro. It gets the event number,
59 // the minimum and maximum order number of TPC tracks that are to be tracked
60 // trough the ITS, and the file where the recpoints are registered. The
61 // method Recursivetracking is a recursive function that performs the
62 // tracking trough the ITS The method Intersection found the layer, ladder
63 // and detector whre the intersection take place and caluclate the
64 // cohordinates of this intersection. It returns an integer that is 0 if the
65 // intersection has been found successfully. The two mwthods Kalmanfilter
66 // and kalmanfiltervert operate the kalmanfilter without and with the vertex
67 // imposition respectively. The authors thank Mariana Bondila to have help
68 // them to resolve some problems. July-2000
70 #include <Riostream.h>
71 #include <Riostream.h>
77 #include <TStopwatch.h>
79 #include "TParticle.h"
82 #include "AliITSsegmentationSSD.h"
83 #include "AliITSgeomSPD.h"
84 #include "AliITSgeomSDD.h"
85 #include "AliITSgeomSSD.h"
86 #include "AliITSgeom.h"
87 #include "AliITSRecPoint.h"
89 #include "AliKalmanTrack.h"
91 #include "AliITSTrackV1.h"
92 #include "AliITSIOTrack.h"
93 #include "AliITSRad.h"
94 #include "../TPC/AliTPCtracker.h"
95 #include "AliITSTrackerV1.h"
96 #include "AliITSVertex.h"
97 #include "AliITSPid.h"
99 ClassImp(AliITSTrackerV1)
100 //______________________________________________________________________
101 AliITSTrackerV1::AliITSTrackerV1() {
102 //Default constructor
114 for(ia=0; ia<6; ia++) {
132 //______________________________________________________________________
133 AliITSTrackerV1::AliITSTrackerV1(AliITS* IITTSS, Int_t evnumber, Bool_t flag) {
134 //Origin A. Badala' and G.S. Pappalardo:
135 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
136 // Class constructor. It does some initializations.
138 //PH Initialisation taken from the default constructor
153 Int_t imax = 200,jmax = 450;
154 frl = new AliITSRad(imax,jmax);
156 ////////// gets information on geometry /////////////////////////////
157 AliITSgeom *g1 = fITS->GetITSgeom();
162 for(ia=0; ia<6; ia++) {
163 fNlad[ia]=g1->GetNladders(ia+1);
164 fNdet[ia]=g1->GetNdetectors(ia+1);
165 //cout<<fNlad[i]<<" "<<fNdet[i]<<"\n";
168 //cout<<" mean radius = ";
170 for(ib=0; ib<6; ib++) {
171 g1->GetCenterThetaPhi(ib+1,ll,dd,det);
172 Double_t r1=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
173 g1->GetCenterThetaPhi(ib+1,ll,dd+1,det);
174 Double_t r2=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
175 fAvrad[ib]=(r1+r2)/2.;
176 //cout<<fAvrad[ib]<<" ";
178 //cout<<"\n"; getchar();
180 fDetx[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDx();
181 fDetz[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDz();
183 fDetx[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDx();
184 fDetz[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDz();
186 fDetx[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDx();
187 fDetz[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDz();
189 fDetx[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDx();
190 fDetz[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDz();
192 fDetx[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDx();
193 fDetz[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDz();
195 fDetx[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDx();
196 fDetz[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDz();
197 //cout<<" Detx Detz\n";
198 //for(Int_t la=0; la<6; la++) cout<<" "<<fDetx[la]<<" "<<
202 // allocate memory and define matrices fzmin, fzmax, fphimin and fphimax //
204 Double_t epszdrift=0.05;
206 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
207 Int_t im1, im2, im2max;
208 for(im1=0; im1<6; im1++) {
210 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
213 for(im1=0; im1<6; im1++) {
215 for(im2=0; im2<im2max; im2++) {
216 g1->GetCenterThetaPhi(im1+1,1,im2+1,det);
217 if(im2!=0) fzmin[im1][im2]=det(2)-fDetz[im1];
219 fzmin[im1][im2]=det(2)-(fDetz[im1])*epsz;
220 if(im2!=(im2max-1)) fzmax[im1][im2]=det(2)+fDetz[im1];
222 fzmax[im1][im2]=det(2)+fDetz[im1]*epsz;
223 if(im1==2 || im1==3) {
224 fzmin[im1][im2]-=epszdrift;
225 fzmax[im1][im2]+=epszdrift;
226 } // end if im1==2 || im1==3
230 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
231 for(im1=0;im1<6;im1++) {
233 fphimin[im1] = new Double_t[im2max];
234 fphimax[im1] = new Double_t[im2max];
237 fphidet = new Double_t*[6];
238 for(im1=0; im1<6; im1++) {
240 fphidet[im1] = new Double_t[im2max];
243 Double_t global[3],local[3];
244 Double_t pigre=TMath::Pi();
245 Double_t xmin,ymin,xmax,ymax;
247 for(im1=0; im1<6; im1++) {
249 for(im2=0; im2<im2max; im2++) {
251 g1->GetCenterThetaPhi(im1+1,im2+1,idet,det);
252 fphidet[im1][im2] = TMath::ATan2(Double_t(det(1)),
254 if(fphidet[im1][im2]<0.) fphidet[im1][im2]+=2.*pigre;
255 local[1]=local[2]=0.;
256 local[0]= -(fDetx[im1]);
257 if(im1==0) local[0]= (fDetx[im1]); //to take into account
258 // different reference system
259 g1->LtoG(im1+1,im2+1,idet,local,global);
260 xmax=global[0]; ymax=global[1];
261 local[0]= (fDetx[im1]);
262 if(im1==0) local[0]= -(fDetx[im1]);//take into account different
264 g1->LtoG(im1+1,im2+1,idet,local,global);
265 xmin=global[0]; ymin=global[1];
266 fphimin[im1][im2]= TMath::ATan2(ymin,xmin);
267 if(fphimin[im1][im2]<0.) fphimin[im1][im2]+=2.*pigre;
268 fphimax[im1][im2]= TMath::ATan2(ymax,xmax);
269 if(fphimax[im1][im2]<0.) fphimax[im1][im2]+=2.*pigre;
272 //////////////////////////////////////////////////////////////////////////////////////////////////////////
273 /////////////// allocate memory and define vector fNRecPoints and matrices fRecCylR, fRecCylPhi, fRecCylZ /////////////
274 gAlice->GetEvent(evnumber);
275 Int_t NumOfModules = g1->GetIndexMax();
276 fRecCylR = new Double_t *[NumOfModules];
277 fRecCylPhi = new Double_t *[NumOfModules];
278 fRecCylZ = new Double_t *[NumOfModules];
279 AliITSRecPoint *recp;
280 fNRecPoints = new Int_t[NumOfModules];
282 for(Int_t module=0; module<NumOfModules; module++) {
283 fITS->ResetRecPoints();
284 gAlice->TreeR()->GetEvent(module);
285 frecPoints=fITS->RecPoints();
286 Int_t nRecPoints=fNRecPoints[module]=frecPoints->GetEntries();
287 fRecCylR[module] = new Double_t[nRecPoints];
288 fRecCylPhi[module] = new Double_t[nRecPoints];
289 fRecCylZ[module] = new Double_t[nRecPoints];
291 for(ind=0; ind<fNRecPoints[module]; ind++) {
292 recp=(AliITSRecPoint*)frecPoints->UncheckedAt(ind);
293 // Float_t global[3], local[3];
294 Double_t global[3], local[3];
295 local[0]=recp->GetX();
297 local[2]= recp->GetZ();
298 g1->LtoG(module,local,global);
300 Double_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]); // r hit
301 Double_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi(); // phi hit
302 Double_t z = global[2]; // z hit
304 fRecCylR[module][ind]=r;
305 fRecCylPhi[module][ind]=phi;
306 fRecCylZ[module][ind]=z;
311 ////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
314 ////////// gets magnetic field factor //////////////////////////////
316 AliMagF * fieldPointer = gAlice->Field();
317 // fFieldFactor = (Double_t)fieldPointer->Factor();
318 fFieldFactor =(Double_t)fieldPointer-> SolenoidField()/10/.2;
319 // cout<< " field factor = "<<fFieldFactor<<"\n"; getchar();
321 //______________________________________________________________________
322 AliITSTrackerV1::AliITSTrackerV1(const AliITSTrackerV1 &cobj) {
323 // Origin A. Badala' and G.S. Pappalardo:
324 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
328 *fresult = *cobj.fresult;
329 fPtref = cobj.fPtref;
330 fChi2max = cobj.fChi2max;
331 **fvettid = **cobj.fvettid;
332 fflagvert = cobj.fflagvert;
333 Int_t imax=200,jmax=450;
334 frl = new AliITSRad(imax,jmax);
336 fFieldFactor = cobj.fFieldFactor;
337 Int_t i,im1,im2,im2max;
339 fNlad[i] = cobj.fNlad[i];
340 fNdet[i] = cobj.fNdet[i];
341 fAvrad[i] = cobj.fAvrad[i];
342 fDetx[i] = cobj.fDetx[i];
343 fDetz[i] = cobj.fDetz[i];
345 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
346 for(im1=0; im1<6; im1++) {
348 fzmin[im1] = new Double_t[im2max];
349 fzmax[im1] = new Double_t[im2max];
351 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
352 for(im1=0;im1<6;im1++) {
354 fphimin[im1] = new Double_t[im2max];
355 fphimax[im1] = new Double_t[im2max];
358 fphidet = new Double_t*[6];
359 for(im1=0; im1<6; im1++) {
361 fphidet[im1] = new Double_t[im2max];
363 for(im1=0; im1<6; im1++) {
365 for(im2=0; im2<im2max; im2++) {
366 fzmin[im1][im2]=cobj.fzmin[im1][im2];
367 fzmax[im1][im2]=cobj.fzmax[im1][im2];
370 for(im1=0; im1<6; im1++) {
372 for(im2=0; im2<im2max; im2++) {
373 fphimin[im1][im2]=cobj.fphimin[im1][im2];
374 fphimax[im1][im2]=cobj.fphimax[im1][im2];
375 fphidet[im1][im2]=cobj.fphidet[im1][im2];
380 AliITSgeom *g1 = fITS->GetITSgeom();
381 Int_t NumOfModules = g1->GetIndexMax();
383 fRecCylR = new Float_t *[NumOfModules];
384 fRecCylPhi = new Float_t *[NumOfModules];
385 fRecCylZ = new Float_t *[NumOfModules];
387 fRecCylR = new Double_t *[NumOfModules];
388 fRecCylPhi = new Double_t *[NumOfModules];
389 fRecCylZ = new Double_t *[NumOfModules];
390 fNRecPoints = new Int_t[NumOfModules];
391 for(Int_t module=0; module<NumOfModules; module++) {
392 Int_t nRecPoints=fNRecPoints[module]=cobj.fNRecPoints[module];
394 fRecCylR[module] = new Float_t[nRecPoints];
395 fRecCylPhi[module] = new Float_t[nRecPoints];
396 fRecCylZ[module] = new Float_t[nRecPoints];
398 fRecCylR[module] = new Double_t[nRecPoints];
399 fRecCylPhi[module] = new Double_t[nRecPoints];
400 fRecCylZ[module] = new Double_t[nRecPoints];
402 for(ind=0; ind<nRecPoints; ind++) {
403 fRecCylR[module][ind]=cobj.fRecCylR[module][ind];
404 fRecCylPhi[module][ind]=cobj.fRecCylPhi[module][ind];
405 fRecCylZ[module][ind]=cobj.fRecCylZ[module][ind];
410 void AliITSTrackerV1::DelMatrix(Int_t NumOfModules) {
411 for(Int_t mod=0; mod<NumOfModules; mod++) {
412 delete fRecCylR[mod];
413 delete fRecCylPhi[mod];
414 delete fRecCylZ[mod];
420 //______________________________________________________________________
421 AliITSTrackerV1::~AliITSTrackerV1(){
422 // Origin A. Badala' and G.S. Pappalardo:
423 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
427 for(Int_t i=0; i<6; i++) {
442 //______________________________________________________________________
443 AliITSTrackerV1 &AliITSTrackerV1::operator=(AliITSTrackerV1 obj) {
444 // Origin A. Badala' and G.S. Pappalardo:
445 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
446 // assignement operator
449 *fresult = *obj.fresult;
451 fChi2max = obj.fChi2max;
452 **fvettid = **obj.fvettid;
453 fflagvert = obj.fflagvert;
454 Int_t imax=200,jmax=450;
455 frl = new AliITSRad(imax,jmax);
457 fFieldFactor = obj.fFieldFactor;
460 fNlad[i] = obj.fNlad[i];
461 fNdet[i] = obj.fNdet[i];
462 fAvrad[i] = obj.fAvrad[i];
463 fDetx[i] = obj.fDetx[i];
464 fDetz[i] = obj.fDetz[i];
466 fzmin = new Double_t*[6];
467 fzmax = new Double_t*[6];
468 Int_t im1, im2, im2max;
469 for(im1=0; im1<6; im1++) {
471 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
473 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
474 for(im1=0;im1<6;im1++) {
476 fphimin[im1] = new Double_t[im2max];
477 fphimax[im1] = new Double_t[im2max];
480 fphidet = new Double_t*[6];
481 for(im1=0; im1<6; im1++) {
483 fphidet[im1] = new Double_t[im2max];
485 for(im1=0; im1<6; im1++) {
487 for(im2=0; im2<im2max; im2++) {
488 fzmin[im1][im2]=obj.fzmin[im1][im2];
489 fzmax[im1][im2]=obj.fzmax[im1][im2];
492 for(im1=0; im1<6; im1++) {
494 for(im2=0; im2<im2max; im2++) {
495 fphimin[im1][im2]=obj.fphimin[im1][im2];
496 fphimax[im1][im2]=obj.fphimax[im1][im2];
497 fphidet[im1][im2]=obj.fphidet[im1][im2];
501 AliITSgeom *g1 = fITS->GetITSgeom();
502 Int_t NumOfModules = g1->GetIndexMax();
503 fRecCylR = new Double_t *[NumOfModules];
504 fRecCylPhi = new Double_t *[NumOfModules];
505 fRecCylZ = new Double_t *[NumOfModules];
506 fNRecPoints = new Int_t[NumOfModules];
507 for(Int_t module=0; module<NumOfModules; module++) {
508 Int_t nRecPoints=fNRecPoints[module]=obj.fNRecPoints[module];
509 fRecCylR[module] = new Double_t[nRecPoints];
510 fRecCylPhi[module] = new Double_t[nRecPoints];
511 fRecCylZ[module] = new Double_t[nRecPoints];
513 for(ind=0; ind<nRecPoints; ind++) {
514 fRecCylR[module][ind]=obj.fRecCylR[module][ind];
515 fRecCylPhi[module][ind]=obj.fRecCylPhi[module][ind];
516 fRecCylZ[module][ind]=obj.fRecCylZ[module][ind];
523 //______________________________________________________________________
524 void AliITSTrackerV1::DoTracking(Int_t evNumber,Int_t minTr,Int_t maxTr,
525 TFile *file, Bool_t realmass) {
526 // Origin A. Badala' and G.S. Pappalardo:
527 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
528 // The method needs the event number, the minimum and maximum order
529 // number of TPC tracks that
530 // are to be tracked trough the ITS, and the file where the recpoints
532 // The method can be called by a macro. It preforms the tracking for
533 // all good TPC tracks
535 printf("begin DoTracking - file %p\n",file);
537 gAlice->GetEvent(evNumber); //modificato per gestire hbt
539 AliKalmanTrack *kkprov;
540 //kkprov->SetConvConst(100/0.299792458/0.2/fFieldFactor);
541 kkprov->SetConvConst(1000/0.299792458/gAlice->Field()->SolenoidField());
542 // cout<<" field = "<<gAlice->Field()->SolenoidField()<<endl;
545 TFile *cf=TFile::Open("AliTPCclusters.root");
546 AliTPCParam *digp= (AliTPCParam*)cf->Get("75x40_100x60_150x60");
547 if (!digp) { cerr<<"TPC parameters have not been found !\n"; getchar();}
550 AliTPCtracker *tracker = new AliTPCtracker(digp,evNumber);
553 tracker->LoadInnerSectors();
554 tracker->LoadOuterSectors();
557 TFile *tf=TFile::Open("AliTPCtracksSorted.root");
559 cerr<<"Can't open AliTPCtracksSorted.root !\n";
562 TObjArray tracks(200000);
564 sprintf(tname,"TreeT_TPC_%d",evNumber);
566 TTree *tracktree=(TTree*)tf->Get(tname);
567 if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n";}
568 TBranch *tbranch=tracktree->GetBranch("tracks");
569 Int_t nentr=(Int_t)tracktree->GetEntries();
572 AliITSRecPoint *recp; // oggi
573 AliTPCtrack *ioTrackTPC=0;
574 for (kk=0; kk<nentr; kk++) {
575 ioTrackTPC=new AliTPCtrack;
576 tbranch->SetAddress(&ioTrackTPC);
577 tracktree->GetEvent(kk);
578 tracker->CookLabel(ioTrackTPC,0.1);
579 tracks.AddLast(ioTrackTPC);
584 Int_t nt = tracks.GetEntriesFast();
585 cerr<<"Number of found tracks "<<nt<<endl;
588 TTree *tr=gAlice->TreeR();
589 Int_t nent=(Int_t)tr->GetEntries();
590 frecPoints = fITS->RecPoints();
594 Int_t *np = new Int_t[nent];
595 fvettid = new Int_t* [nent];
598 for (mod=0; mod<nent; mod++) {
600 fITS->ResetRecPoints();
601 gAlice->TreeR()->GetEvent(mod);
602 numbpoints = frecPoints->GetEntries();
603 totalpoints+=numbpoints;
604 np[mod] = numbpoints;
605 //cout<<" mod = "<<mod<<" numbpoints = "<<numbpoints<<"\n";getchar();
606 fvettid[mod] = new Int_t[numbpoints];
608 for (ii=0;ii<numbpoints; ii++) *(fvettid[mod]+ii)=0;
611 AliTPCtrack *track=0;
613 if(minTr < 0) {minTr = 0; maxTr = nt-1;}
617 TTree tracktree1("TreeT","Tree with ITS tracks");
618 AliITSIOTrack *ioTrack=0;
619 AliITSPid *pid=new AliITSPid(1000); // oggi
621 tracktree1.Branch("ITStracks","AliITSIOTrack",&ioTrack,32000,0);
623 TDatabasePDG * db = new TDatabasePDG;
626 for (j=minTr; j<=maxTr; j++) {
627 track=(AliTPCtrack*)tracks.UncheckedAt(j);
628 if (!track) continue;
630 /// mass definition ////////////////////////
631 Double_t mass=0.13956995;
632 Int_t pcode=211; // a pion by default
635 if(TMath::Abs(pcode)<20443) mass=db->GetParticle(pcode)->Mass();
638 mass = track->GetMass();
639 // cout << "Mass = " << mass << endl;
644 // new propagation to the end of TPC
646 // track->PropagateTo(xk,0.,0.); //Ne if it's still there //attenzione funziona solo se modifica in TPC
647 // Double_t xk=77.415;
648 track->PropagateTo(xk, 28.94, 1.204e-3);
650 track->PropagateTo(xk, 44.77,1.71); //Tedlar
652 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
654 track->PropagateTo(xk, 41.28, 0.029);//Nomex
656 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
658 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
661 // track->PropagateTo(xk,0.,0.); //C02
662 track->PropagateTo(xk,36.2,1.98e-3); //C02 //attenzione funziona solo se modifica in TPC
665 track->PropagateTo(xk, 24.01, 2.7); //Al
667 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
669 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
671 track->PropagateTo(xk, 41.28, 0.029); //Nomex
673 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
675 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
677 track->PropagateTo(xk, 24.01, 2.7); //Al
679 ////////////////////////////////////////////////////////////////////////////////////////////////////////
680 //AliITSTrackV1 trackITS(*track);
681 AliITSTrackV1 trackITS(*track, fFieldFactor);
682 //cout<<" fFieldFactor = "<<fFieldFactor<<"\n";
683 trackITS.PutMass(mass); //new to add mass to track
684 if(fresult){ delete fresult; fresult=0;}
685 fresult = new AliITSTrackV1(trackITS);
687 AliITSTrackV1 primaryTrack(trackITS);
688 vgeant=(*fresult).GetVertex();
690 // Definition of dv and zv for vertex constraint
691 Double_t sigmaDv=0.0050; Double_t sigmaZv=0.010;
692 //Double_t sigmaDv=0.0015; Double_t sigmaZv=0.0015;
693 Double_t uniform= gRandom->Uniform();
695 if(uniform<=0.5) signdv=-1.;
699 Double_t vr=TMath::Sqrt(vgeant(0)*vgeant(0)+ vgeant(1)*vgeant(1));
700 Double_t dv=gRandom->Gaus(signdv*vr,(Float_t)sigmaDv);
701 Double_t zv=gRandom->Gaus(vgeant(2),(Float_t)sigmaZv);
702 //cout<<" Dv e Zv = "<<dv<<" "<<zv<<"\n";
705 trackITS.SetsigmaDv(sigmaDv);
706 trackITS.SetsigmaZv(sigmaZv);
707 (*fresult).SetDv(dv);
708 (*fresult).SetZv(zv);
709 (*fresult).SetsigmaDv(sigmaDv);
710 (*fresult).SetsigmaZv(sigmaZv);
711 primaryTrack.SetDv(dv);
712 primaryTrack.SetZv(zv);
713 primaryTrack.SetsigmaDv(sigmaDv);
714 primaryTrack.SetsigmaZv(sigmaZv);
715 primaryTrack.PrimaryTrack(frl);
716 TVector d2=primaryTrack.Getd2();
717 TVector tgl2=primaryTrack.Gettgl2();
718 TVector dtgl=primaryTrack.Getdtgl();
719 trackITS.Setd2(d2); trackITS.Settgl2(tgl2);
720 trackITS.Setdtgl(dtgl);
721 (*fresult).Setd2(d2); (*fresult).Settgl2(tgl2);
722 (*fresult).Setdtgl(dtgl);
724 trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
725 (*result).SetVertex(vertex); (*result).SetErrorVertex(ervertex);
727 TList *list= new TList();
729 list->AddLast(&trackITS);
731 fPtref=TMath::Abs( (trackITS).GetPt() );
732 //cout<<" fPtref = " <<fPtref<<"\n";
733 if(fPtref>1.0) fChi2max=40.;
734 if(fPtref<=1.0) fChi2max=20.;
735 if(fPtref<0.4 ) fChi2max=100.;
736 if(fPtref<0.2 ) fChi2max=40.;
737 // if(fPtref<0.4 ) fChi2max=30.;
738 // if(fPtref<0.2 ) fChi2max=20.;
739 //if(fPtref<0.2 ) fChi2max=10.;
740 //if(fPtref<0.1 ) fChi2max=5.;
741 //cout << "\n Pt = " << fPtref <<"\n"; //stampa
742 RecursiveTracking(list);
747 TVector vecTotLabRef(18);
749 for(lay=5; lay>=0; lay--) {
750 TVector vecLabRef(3);
751 vecLabRef=(*fresult).GetLabTrack(lay);
752 Float_t clustZ=(*fresult).GetZclusterTrack( lay);
754 Int_t lpp=(Int_t)vecLabRef(k);
756 TParticle *p=(TParticle*) gAlice->Particle(lpp);
757 Int_t pcode=p->GetPdgCode();
758 if(pcode==11) vecLabRef(k)=p->GetFirstMother();
760 itot++; vecTotLabRef(itot)=vecLabRef(k);
761 if(vecLabRef(k)==0. && clustZ == -1.) vecTotLabRef(itot) =-3.;
766 (*fresult).Search(vecTotLabRef, labref, freq);
768 //if(freq < 6) labref=-labref; // cinque - sei
769 if(freq < 5) labref=-labref; // cinque - sei
770 (*fresult).SetLabel(labref);
772 // cout<<" progressive track number = "<<j<<"\r";
774 Int_t numOfCluster=(*fresult).GetNumClust();
775 //cout<<" progressive track number = "<<j<<"\n"; // stampa
776 Long_t labITS=(*fresult).GetLabel();
777 //cout << " ITS track label = " << labITS << "\n"; // stampa
778 Int_t lab=track->GetLabel();
779 //cout << " TPC track label = " << lab <<"\n"; // stampa
780 //propagation to vertex
783 if((*fresult).DoNotCross(rbeam)) continue; //no intersection with beampipe
784 (*fresult).Propagation(rbeam);
785 Double_t c00,c10,c11,c20,c21,c22,c30,c31,c32,c33,c40,c41,c42,c43,c44;
786 (*fresult).GetCElements(c00,
790 c40,c41,c42,c43,c44);
792 Double_t pt=TMath::Abs((*fresult).GetPt());
793 Double_t dr=(*fresult).GetD();
794 Double_t z=(*fresult).GetZ();
795 Double_t tgl=(*fresult).GetTgl();
796 Double_t c=(*fresult).GetC();
798 Double_t dz=z-(tgl/cy)*TMath::ASin((*fresult).Arga(rbeam));
800 // cout<<" dr e dz alla fine = "<<dr<<" "<<dz<<"\n"; getchar();
801 Double_t phi=(*fresult).Getphi();
802 Double_t phivertex = phi - TMath::ASin((*fresult).ArgA(rbeam));
803 Double_t duepi=2.*TMath::Pi();
804 if(phivertex>duepi) phivertex-=duepi;
805 if(phivertex<0.) phivertex+=duepi;
806 /////////////////////////////////////////////////////////////
807 Int_t idmodule,idpoint;
808 if(numOfCluster >=5) { // cinque - sei
809 //if(numOfCluster ==6) { // cinque - sei
810 AliITSIOTrack outTrack;
812 ioTrack->SetStatePhi(phi);
813 ioTrack->SetStateZ(z);
814 ioTrack->SetStateD(dr);
815 ioTrack->SetStateTgl(tgl);
816 ioTrack->SetStateC(c);
817 Double_t radius=(*fresult).Getrtrack();
818 ioTrack->SetRadius(radius);
820 if(c>0.) charge=-1; else charge=1;
821 ioTrack->SetCharge(charge);
822 Double_t trackmass=(*fresult).GetMass(); // oggi
823 ioTrack->SetMass(trackmass); // oggi
824 ioTrack->SetCovMatrix(c00,
828 c40,c41,c42,c43,c44);
829 Double_t px=pt*TMath::Cos(phivertex);
830 Double_t py=pt*TMath::Sin(phivertex);
832 Double_t xtrack=dr*TMath::Sin(phivertex);
833 Double_t ytrack=dr*TMath::Cos(phivertex);
834 Double_t ztrack=dz+vgeant(2);
838 ioTrack->SetX(xtrack);
839 ioTrack->SetY(ytrack);
840 ioTrack->SetZ(ztrack);
841 ioTrack->SetLabel(labITS);
842 ioTrack->SetTPCLabel(lab);
846 for(il=0;il<6; il++){
847 ioTrack->SetIdPoint(il,(*fresult).GetIdPoint(il));
848 ioTrack->SetIdModule(il,(*fresult).GetIdModule(il));
852 Float_t q[4]={-1.,-1.,-1.,-1.};
853 Float_t globaldedx=0.;
854 for (il=0;il<6;il++) {
855 idpoint=(*fresult).GetIdPoint(il);
856 idmodule=(*fresult).GetIdModule(il);
857 if(idmodule>0.) *(fvettid[idmodule]+idpoint)=1;
859 ioTrack->SetIdPoint(il,idpoint);
860 ioTrack->SetIdModule(il,idmodule);
861 //// for q definition
864 fITS->ResetRecPoints();
865 gAlice->TreeR()->GetEvent(idmodule);
866 recp=(AliITSRecPoint*)frecPoints->UncheckedAt(idpoint);
867 q[il-2]=recp->GetQ()*(*fresult).Getfcor(il-2);
871 q[0]/=280.; q[1]/=280.;
872 q[2]/=38.; q[3]/=38.;
874 // cout<<" q prima = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar();
879 for (il=0; il<3; il++) {
880 if (q[il]<=q[il+1]) continue;
882 q[il]=q[il+1]; q[il+1]=tmp;
888 // cout<<" q dopo = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar();
897 // cout<<" q dopo if = "<<q[0]<<" "<<q[1]<<" "<<q[2]<<" "<<q[3]<<"\n"; getchar();
899 globaldedx=(q[0]+q[1])/2.;
901 // if(q[3]> 0.) globaldedx=(q[0]+q[1]+q[2]+q[3])/4.;
902 // else globaldedx=(q[0]+q[1]+q[2])/3.;
904 ioTrack->SetdEdx(globaldedx);
905 ioTrack->SetPid(pid->GetPcode(ioTrack));
908 } // end if on numOfCluster
909 //gObjectTable->Print(); // stampa memoria
910 } // end for (int j=minTr; j<=maxTr; j++)
912 static Bool_t first=kTRUE;
915 tfile=new TFile("itstracks.root","RECREATE");
916 //cout<<"I have opened itstracks.root file "<<endl;
922 sprintf(hname,"TreeT%d",evNumber);
923 cout << "Number of saved ITS tracks " << tracktree1.GetEntries() << endl;
924 tracktree1.Write(hname);
926 TTree *fAli=gAlice->TreeK();
928 if (fAli) fileAli =fAli->GetCurrentFile();
930 ////////////////////////////////////////////////////////////////////
932 printf("delete vectors\n");
934 if(fvettid) delete [] fvettid;
935 if(fresult) {delete fresult; fresult=0;}
937 //______________________________________________________________________
938 void AliITSTrackerV1::RecursiveTracking(TList *trackITSlist) {
939 // This function perform the recursive tracking in ITS detectors
940 // reference is a pointer to the final best track
941 // Origin A. Badala' and G.S. Pappalardo:
942 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
943 // The authors thank Mariana Bondila to have help them to resolve some
944 // problems. July-2000
946 //Rlayer[0]=4.; Rlayer[1]=7.; Rlayer[2]=14.9;
947 // Rlayer[3]=23.8; Rlayer[4]=39.1; Rlayer[5]=43.6; //vecchio
949 //////////////////////
950 Float_t sigmaphil[6], sigmazl[6];
951 sigmaphil[0]=1.44e-6/(fAvrad[0]*fAvrad[0]);
952 sigmaphil[1]=1.44e-6/(fAvrad[1]*fAvrad[1]);
953 sigmaphil[2]=1.444e-5/(fAvrad[2]*fAvrad[2]);
954 sigmaphil[3]=1.444e-5/(fAvrad[3]*fAvrad[3]);
955 sigmaphil[4]=4e-6/(fAvrad[4]*fAvrad[4]);
956 sigmaphil[5]=4e-6/(fAvrad[5]*fAvrad[5]);
963 ///////////////////////////////////////////////////////////
965 AliITSgeom *g1 = fITS->GetITSgeom();
966 AliITSRecPoint *recp;
967 for(index =0; index<trackITSlist->GetSize(); index++) {
968 AliITSTrackV1 *trackITS = (AliITSTrackV1 *) trackITSlist->At(index);
969 if((*trackITS).GetLayer()==7) fresult->SetChi2(10.223e140);
970 // cout <<" Layer inizio = "<<(*trackITS).GetLayer()<<"\n";
971 // cout<<"fvtrack =" <<"\n";
972 // cout << (*trackITS)(0) << " "<<(*trackITS)(1)<<" "
973 // <<(*trackITS)(2)<<" "<<(*trackITS)(3)<<" "
974 // <<(*trackITS)(4)<<"\n";
975 // cout<< " rtrack = "<<(*trackITS).Getrtrack()<<"\n";
976 // cout<< " Pt = "<<(*trackITS).GetPt()<<"\n";
978 Double_t chi2Now, chi2Ref;
979 Float_t numClustRef = fresult->GetNumClust();
980 if((*trackITS).GetLayer()==1 ) {
981 chi2Now = trackITS->GetChi2();
982 Float_t numClustNow = trackITS->GetNumClust();
983 if(trackITS->GetNumClust())
984 chi2Now /= (Double_t)trackITS->GetNumClust();
985 chi2Ref = fresult->GetChi2();
986 if(fresult->GetNumClust())
987 chi2Ref /= (Double_t)fresult->GetNumClust();
988 //cout<<" chi2Now and chi2Ref = "<<chi2Now<<" "<<chi2Ref<<"\n";
989 if( numClustNow > numClustRef ) {*fresult = *trackITS;}
990 if((numClustNow == numClustRef )&&
991 (chi2Now < chi2Ref)) {
992 *fresult = *trackITS;
997 if(trackITS->Getfnoclust()>=2) continue;
998 Float_t numClustNow = trackITS->GetNumClust();
1000 chi2Now = trackITS->GetChi2();
1002 if(numClustNow<numClustRef && chi2Now>fresult->GetChi2()) continue;
1003 //cout<<" chi2Now = "<<chi2Now<<"\n";
1005 chi2Now/=numClustNow;
1006 if(fPtref > 1.0 && chi2Now > 30.) continue;
1007 if((fPtref >= 0.6 && fPtref<=1.0) && chi2Now > 40.) continue;
1008 // if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 40.) continue;
1009 // if(fPtref <= 0.2 && chi2Now > 8.) continue;
1010 if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 30.) continue;
1011 if(fPtref <= 0.2 && chi2Now > 7.) continue;
1012 /////////////////////////////
1015 Int_t layerInit = (*trackITS).GetLayer();
1016 Int_t layernew = layerInit - 2;// -1 for new layer, -1 for matrix index
1018 Int_t ladp, ladm, detp,detm,ladinters,detinters;
1019 Int_t layerfin=layerInit-1;
1020 // cout<<"Prima di intersection \n";
1021 Int_t outinters=Intersection(*trackITS,layerfin,ladinters,detinters);
1022 // cout<<" outinters = "<<outinters<<"\n";
1023 // cout<<" Layer ladder detector intersection ="
1024 // <<layerfin<<" "<<ladinters<<" "<<detinters<<"\n";
1025 // cout << " phiinters zinters = "<<(*trackITS)(0)
1026 // << " "<<(*trackITS)(1)<<"\n"; getchar();
1027 if(outinters==-1) continue;
1029 (*trackITS).SetLayer(layerfin); // oggi
1030 (*trackITS).Setfcor(); // oggi
1032 TVector toucLad(9), toucDet(9);
1033 Int_t lycur=layerfin;
1036 if(ladm <= 0) ladm=fNlad[layerfin-1];
1037 if(ladp > fNlad[layerfin-1]) ladp=1;
1042 toucLad(0)=ladinters; toucLad(1)=ladm; toucLad(2)=ladp;
1043 toucLad(3)=ladinters; toucLad(4)=ladm; toucLad(5)=ladp;
1044 toucLad(6)=ladinters; toucLad(7)=ladm; toucLad(8)=ladp;
1045 toucDet(0)=detinters; toucDet(1)=detinters; toucDet(2)=detinters;
1046 if(detm > 0 && detp <= fNdet[layerfin-1]) {
1048 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1049 toucDet(6)=detp; toucDet(7)=detp; toucDet(8)=detp;
1051 if(detm > 0 && detp > fNdet[layerfin-1]) {
1053 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1055 if(detm <= 0 && detp <= fNdet[layerfin-1]) {
1057 toucDet(3)=detp; toucDet(4)=detp; toucDet(5)=detp;
1060 Float_t epsphi=5.0, epsz=5.0;
1061 if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1062 // new definition of idetot e toucLad e toucDet to be
1063 // transformed in a method
1064 // these values could be modified
1065 Float_t pigre=TMath::Pi();
1066 Float_t rangephi=5., rangez=5.;
1067 if(layerfin==1 || layerfin ==2){
1068 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1069 (*trackITS).GetSigmaphi());
1070 rangez = 40.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1071 (*trackITS).GetSigmaZ());
1073 if(layerfin==3 || layerfin ==4){
1074 //rangephi=30.*fepsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1075 // (*trackITS).GetSigmaphi());
1076 //rangez = 40.*fepsz*TMath::Sqrt(sigmazl[layerfin-1]+
1077 // (*trackITS).GetSigmaZ());
1078 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1079 (*trackITS).GetSigmaphi());
1080 rangez = 50.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1081 (*trackITS).GetSigmaZ());
1083 if(layerfin==5 || layerfin ==6){
1084 rangephi=20.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1085 (*trackITS).GetSigmaphi());
1086 rangez =5.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1087 (*trackITS).GetSigmaZ());
1089 Float_t phinters, zinters;
1090 phinters=(*trackITS).Getphi();
1091 zinters=(*trackITS).GetZ();
1092 Float_t distz = 0.0;
1093 Float_t phicm, phicp, distphim, distphip;
1095 if(phinters>fphimax[layerfin-1][ladm-1]) phicm=phinters-2*pigre; //corretto il 20-11-2001
1096 distphim=TMath::Abs(phicm-fphimax[layerfin-1][ladm-1]); //corretto il 20-11-2001
1098 //cout<<" fNlad[layerfin-1] e ladp = "<<fNlad[layerfin-1]<<" "<<ladp<<endl;
1099 if(phinters>fphimin[layerfin-1][ladp-1]) phicp=phinters-2.*pigre; //corretto il 20-11-2001
1100 distphip=TMath::Abs(phicp-fphimin[layerfin-1][ladp-1]); //corretto il 20-11-2001
1104 toucLad(0)=ladinters; toucDet(0)=detinters;
1105 if(detm>0) distz=TMath::Abs(zinters-fzmax[layerfin-1][detm-1]);
1106 if(detm>0 && rangez>=distz){
1108 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detm;
1109 if(rangephi>=distphim){
1111 toucLad(idetot-1)=ladm;
1112 toucDet(idetot-1)=detinters;
1114 toucLad(idetot-1)=ladm;
1115 toucDet(idetot-1)=detm;
1117 if(rangephi>=distphip){
1119 toucLad(idetot-1)=ladp;
1120 toucDet(idetot-1)=detinters;
1122 toucLad(idetot-1)=ladp;
1123 toucDet(idetot-1)=detm;
1126 if(detp<=fNdet[layerfin-1])
1127 distz=TMath::Abs(zinters-fzmin[layerfin-1][detp-1]);
1128 if(detp<=fNdet[layerfin-1] && rangez>=distz){
1130 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detp;
1131 if(rangephi>=distphim){
1132 idetot++; toucLad(idetot-1)=ladm; toucDet(idetot-1)=detp;
1135 toucLad(idetot-1)=ladm;
1136 toucDet(idetot-1)=detinters;
1139 if(rangephi>=distphip){
1141 toucLad(idetot-1)=ladp;
1142 toucDet(idetot-1)=detp;
1145 toucLad(idetot-1)=ladp;
1146 toucDet(idetot-1)=detinters;
1149 } //end detm<fNdet[.......
1151 if(flagzmin == 0 && flagzmax==0){
1152 if(rangephi>=distphim){
1154 toucLad(idetot-1)=ladm;
1155 toucDet(idetot-1)=detinters;
1157 if(rangephi>=distphip){
1159 toucLad(idetot-1)=ladp;
1160 toucDet(idetot-1)=detinters;
1163 ////////////////////////////////////////////////////////////
1165 for (iriv=0; iriv<idetot; iriv++) { //for on detectors
1166 ///////////////////////////////////////////////////////
1167 /*** Rec points sorted by module *****/
1168 /**************************************/
1170 indexmod = g1->GetModuleIndex(lycur,(Int_t)toucLad(iriv),
1171 (Int_t)toucDet(iriv));
1172 fITS->ResetRecPoints();
1173 gAlice->TreeR()->GetEvent(indexmod);
1174 Int_t npoints=frecPoints->GetEntries();
1177 for(indnew=0; indnew<npoints; indnew++){
1178 if (*(fvettid[indexmod]+indnew)==0)
1179 recp =(AliITSRecPoint*)frecPoints->UncheckedAt(indnew);
1182 TVector cluster(3),vecclust(9);
1183 //vecclust(6)=vecclust(7)=vecclust(8)=-1.;
1185 // now vecclust is with cylindrical cohordinates
1186 vecclust(0)=(Float_t)fRecCylR[indexmod][indnew];
1187 vecclust(1)=(Float_t)fRecCylPhi[indexmod][indnew];
1188 vecclust(2)=(Float_t)fRecCylZ[indexmod][indnew];
1189 vecclust(3) = (Float_t)recp->fTracks[0];
1190 vecclust(4) = (Float_t)indnew;
1191 vecclust(5) = (Float_t)indexmod;
1192 vecclust(6) = (Float_t)recp->fTracks[0];
1193 vecclust(7) = (Float_t)recp->fTracks[1];
1194 vecclust(8) = (Float_t)recp->fTracks[2];
1195 sigma[0] = (Double_t) recp->GetSigmaX2();
1196 sigma[1] = (Double_t) recp->GetSigmaZ2();
1198 cluster(0)=fRecCylR[indexmod][indnew];
1199 cluster(1)=fRecCylPhi[indexmod][indnew];
1200 cluster(2)=fRecCylZ[indexmod][indnew];
1202 // cout<<" layer = "<<play<<"\n";
1203 // cout<<" cluster prima = "<<vecclust(0)<<" "
1204 // <<vecclust(1)<<" "
1205 // <<vecclust(2)<<"\n"; getchar();
1207 Float_t sigmatotphi, sigmatotz;
1208 // Float_t epsphi=5.0, epsz=5.0;
1209 //if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1210 Double_t rTrack=(*trackITS).Getrtrack();
1211 Double_t sigmaphi=sigma[0]/(rTrack*rTrack);
1212 sigmatotphi=epsphi*TMath::Sqrt(sigmaphi +
1213 (*trackITS).GetSigmaphi());
1214 sigmatotz=epsz*TMath::Sqrt(sigma[1] +
1215 (*trackITS).GetSigmaZ());
1216 //cout<<"cluster e sigmatotphi e track = "<<cluster(0)
1217 // <<" "<<cluster(1)<<" "<<sigmatotphi<<" "
1218 // <<vecclust(3)<<"\n";
1219 //if(vecclust(3)==481) getchar();
1220 if(cluster(1)<6. && (*trackITS).Getphi()>6.)
1221 cluster(1)=cluster(1)+(2.*TMath::Pi());
1222 if(cluster(1)>6. && (*trackITS).Getphi()<6.)
1223 cluster(1)=cluster(1)-(2.*TMath::Pi());
1224 if(TMath::Abs(cluster(1)-(*trackITS).Getphi())>sigmatotphi)
1226 // cout<<" supero sigmaphi \n";
1227 AliITSTrackV1 *newTrack = new AliITSTrackV1((*trackITS));
1228 //(*newTrack).SetLayer((*trackITS).GetLayer()-1);
1229 if (TMath::Abs(rTrack-cluster(0))/rTrack>1e-6)
1230 (*newTrack).Correct(Double_t(cluster(0)));
1231 //cout<<" cluster(2) e(*newTrack).GetZ()="<<cluster(2)<<" "
1232 // << (*newTrack).GetZ()<<"\n";
1233 if(TMath::Abs(cluster(2)-(*newTrack).GetZ()) > sigmatotz){
1237 Double_t sigmanew[2];
1238 sigmanew[0]= sigmaphi;
1239 sigmanew[1]=sigma[1];
1243 // Double_t chi2pred=newTrack->GetPredChi2(m,sigmanew);
1244 // cout<<" chi2pred = "<<chi2pred<<"\n";
1245 // if(chi2pred>fChi2max) continue; //aggiunto il 30-7-2001
1246 if(iriv == 0) flaghit=1;
1247 (*newTrack).AddMS(frl); // add the multiple scattering
1248 //matrix to the covariance matrix
1249 (*newTrack).AddEL(frl,1.,0);
1252 KalmanFilterVert(newTrack,cluster,sigmanew);
1253 //KalmanFilterVert(newTrack,cluster,sigmanew,chi2pred);
1255 KalmanFilter(newTrack,cluster,sigmanew);
1257 (*newTrack).PutCluster(layernew, vecclust);
1258 newTrack->AddClustInTrack();
1259 listoftrack.AddLast(newTrack);
1261 } // end of for on detectors (iriv)
1262 }//end if(outinters==0)
1264 if(flaghit==0 || outinters==-2) {
1265 AliITSTrackV1 *newTrack = new AliITSTrackV1(*trackITS);
1266 (*newTrack).Setfnoclust();
1267 //(*newTrack).SetLayer((*trackITS).GetLayer()-1);
1268 (*newTrack).AddMS(frl); // add the multiple scattering matrix
1269 // to the covariance matrix
1270 (*newTrack).AddEL(frl,1.,0);
1271 listoftrack.AddLast(newTrack);
1274 //gObjectTable->Print(); // stampa memoria
1276 RecursiveTracking(&listoftrack);
1277 listoftrack.Delete();
1278 } // end of for on tracks (index)
1280 //gObjectTable->Print(); // stampa memoria
1282 //______________________________________________________________________
1283 Int_t AliITSTrackerV1::Intersection(AliITSTrackV1 &track,Int_t layer,
1284 Int_t &ladder,Int_t &detector) {
1285 // Origin A. Badala' and G.S. Pappalardo
1286 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1287 // Found the intersection and the detector
1289 Double_t rk=fAvrad[layer-1];
1290 if(track.DoNotCross(rk)){ /*cout<< " Do not cross \n";*/ return -1;}
1291 track.Propagation(rk);
1292 Double_t zinters=track.GetZ();
1293 Double_t phinters=track.Getphi();
1294 //cout<<"zinters = "<<zinters<<" phinters = "<<phinters<<"\n";
1298 TVector distZCenter(2);
1302 for(iD = 1; iD<= fNdet[layer-1]; iD++) {
1303 if(zinters > fzmin[layer-1][iD-1] && zinters <= fzmax[layer-1][iD-1]) {
1305 cout<< " Errore su iz in Intersection \n";
1308 listDet(iz)= iD; distZCenter(iz)=TMath::Abs(zinters-det(2));
1314 if(iz==0) {/* cout<< " No detector along Z \n";*/ return -2;}
1315 detector=Int_t (listDet(0));
1316 if(iz>1 && (distZCenter(0)>distZCenter(1))) detector=Int_t (listDet(1));
1319 TVector distPhiCenter(2);
1321 Double_t pigre=TMath::Pi();
1323 for(iLd = 1; iLd<= fNlad[layer-1]; iLd++) {
1324 Double_t phimin=fphimin[layer-1][iLd-1];
1325 Double_t phimax=fphimax[layer-1][iLd-1];
1326 Double_t phidet=fphidet[layer-1][iLd-1];
1327 Double_t phiconfr=phinters;
1329 //if(phimin <5.5) {cout<<" Error in Intersection for phi \n";
1332 if(phinters>(1.5*pigre)) phiconfr=phinters-(2.*pigre);
1333 if(phidet>(1.5*pigre)) phidet-=(2.*pigre);
1335 if(phiconfr>phimin && phiconfr<= phimax) {
1337 cout<< " Errore su ip in Intersection \n"; getchar();
1340 distPhiCenter(ip)=TMath::Abs(phiconfr-phidet); ip++;
1344 if(ip==0) { cout<< " No detector along phi \n"; getchar();}
1345 ladder=Int_t (listLad(0));
1346 if(ip>1 && (distPhiCenter(0)>distPhiCenter(1))) ladder=Int_t (listLad(1));
1349 //______________________________________________________________________
1350 void AliITSTrackerV1::KalmanFilter(AliITSTrackV1 *newTrack,TVector &cluster,
1352 //Origin A. Badala' and G.S. Pappalardo:
1353 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1354 // Kalman filter without vertex constraint
1355 ////// Evaluation of the measurement vector ////////////////////////
1357 Double_t rk,phik,zk;
1358 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1360 //////////////////////// Evaluation of the error matrix V /////////
1361 Double_t v00=sigma[0];
1362 Double_t v11=sigma[1];
1363 ////////////////////////////////////////////////////////////////////
1364 Double_t cin00,cin10,cin20,cin30,cin40,cin11,cin21,cin31,cin41,cin22,
1365 cin32,cin42,cin33,cin43,cin44;
1367 newTrack->GetCElements(cin00,
1370 cin30,cin31,cin32,cin33,
1371 cin40,cin41,cin42,cin43,cin44); //get C matrix
1372 Double_t rold00=cin00+v00;
1373 Double_t rold10=cin10;
1374 Double_t rold11=cin11+v11;
1375 ////////////////////// R matrix inversion /////////////////////////
1376 Double_t det=rold00*rold11-rold10*rold10;
1377 Double_t r00=rold11/det;
1378 Double_t r10=-rold10/det;
1379 Double_t r11=rold00/det;
1380 ////////////////////////////////////////////////////////////////////
1381 Double_t k00=cin00*r00+cin10*r10;
1382 Double_t k01=cin00*r10+cin10*r11;
1383 Double_t k10=cin10*r00+cin11*r10;
1384 Double_t k11=cin10*r10+cin11*r11;
1385 Double_t k20=cin20*r00+cin21*r10;
1386 Double_t k21=cin20*r10+cin21*r11;
1387 Double_t k30=cin30*r00+cin31*r10;
1388 Double_t k31=cin30*r10+cin31*r11;
1389 Double_t k40=cin40*r00+cin41*r10;
1390 Double_t k41=cin40*r10+cin41*r11;
1391 Double_t x0,x1,x2,x3,x4;
1392 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1393 Double_t savex0=x0, savex1=x1;
1394 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1);
1395 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1);
1396 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1);
1397 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1);
1398 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1);
1399 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1400 c00=cin00-k00*cin00-k01*cin10;
1401 c10=cin10-k00*cin10-k01*cin11;
1402 c20=cin20-k00*cin20-k01*cin21;
1403 c30=cin30-k00*cin30-k01*cin31;
1404 c40=cin40-k00*cin40-k01*cin41;
1405 c11=cin11-k10*cin10-k11*cin11;
1406 c21=cin21-k10*cin20-k11*cin21;
1407 c31=cin31-k10*cin30-k11*cin31;
1408 c41=cin41-k10*cin40-k11*cin41;
1409 c22=cin22-k20*cin20-k21*cin21;
1410 c32=cin32-k20*cin30-k21*cin31;
1411 c42=cin42-k20*cin40-k21*cin41;
1412 c33=cin33-k30*cin30-k31*cin31;
1413 c43=cin43-k30*cin40-k31*cin41;
1414 c44=cin44-k40*cin40-k41*cin41;
1415 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1416 newTrack->PutCElements(c00,
1420 c40,c41,c42,c43,c44); // put in track the
1422 Double_t vmcold00=v00-c00;
1423 Double_t vmcold10=-c10;
1424 Double_t vmcold11=v11-c11;
1425 ////////////////////// Matrix vmc inversion ///////////////////////
1426 det=vmcold00*vmcold11-vmcold10*vmcold10;
1427 Double_t vmc00=vmcold11/det;
1428 Double_t vmc10=-vmcold10/det;
1429 Double_t vmc11=vmcold00/det;
1430 ////////////////////////////////////////////////////////////////////
1431 Double_t chi2=(m[0]-x0)*( vmc00*(m[0]-x0) + 2.*vmc10*(m[1]-x1) ) +
1432 (m[1]-x1)*vmc11*(m[1]-x1);
1433 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1435 //----------------------------------------------------------------------
1436 //void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1437 // TVector &cluster,Double_t sigma[2]){
1438 void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1439 TVector &cluster,Double_t sigma[2]
1440 /*, Double_t chi2pred*/){
1441 //Origin A. Badala' and G.S. Pappalardo:
1442 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1443 // Kalman filter with vertex constraint
1444 ///////////////////// Evaluation of the measurement vector m ////////
1446 Double_t rk,phik,zk;
1447 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1449 Double_t cc=(*newTrack).GetC();
1450 Double_t zv=(*newTrack).GetZv();
1451 Double_t dv=(*newTrack).GetDv();
1453 Double_t tgl= (zk-zv)*cy/TMath::ASin(cy*rk);
1455 /////////////////////// Evaluation of the error matrix V //////////
1456 Int_t layer=newTrack->GetLayer();
1457 Double_t v00=sigma[0];
1458 Double_t v11=sigma[1];
1459 Double_t v31=sigma[1]/rk;
1460 Double_t sigmaDv=newTrack->GetsigmaDv();
1461 Double_t v22=sigmaDv*sigmaDv + newTrack->Getd2(layer-1);
1462 Double_t v32=newTrack->Getdtgl(layer-1);
1463 Double_t sigmaZv=newTrack->GetsigmaZv();
1464 Double_t v33=(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+newTrack->Gettgl2(layer-1);
1465 //////////////////////////////////////////////////////////////////
1466 Double_t cin00,cin10,cin11,cin20,cin21,cin22,
1467 cin30,cin31,cin32,cin33,cin40,cin41,cin42,cin43,cin44;
1468 newTrack->GetCElements(cin00,
1471 cin30,cin31,cin32,cin33,
1472 cin40,cin41,cin42,cin43,cin44); //get C matrix
1480 r[3][1]=cin31+sigma[1]/rk;
1481 r[2][2]=cin22+sigmaDv*sigmaDv+newTrack->Getd2(layer-1);
1482 r[3][2]=cin32+newTrack->Getdtgl(layer-1);
1483 r[3][3]=cin33+(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+
1484 newTrack->Gettgl2(layer-1);
1485 r[0][1]=r[1][0]; r[0][2]=r[2][0]; r[0][3]=r[3][0];
1486 r[1][2]=r[2][1]; r[1][3]=r[3][1]; r[2][3]=r[3][2];
1487 ///////////////////// Matrix R inversion //////////////////////////
1491 Int_t ll[kn],mm[kn];
1494 for(k=0; k<kn; k++) {
1498 for(j=k; j<kn ; j++) {
1499 for (i=j; i<kn; i++) {
1500 if(TMath::Abs(big) < TMath::Abs(r[i][j]) ) {
1510 for(i=0; i<kn; i++) {
1519 for(j=0; j<kn; j++) {
1528 cout << "Singular matrix\n";
1530 for(i=0; i<kn; i++) {
1531 if(i == k) { continue; }
1532 r[i][k]=r[i][k]/(-big);
1535 for(i=0; i<kn; i++) {
1537 for(j=0; j<kn; j++) {
1538 if(i == k || j == k) continue;
1539 r[i][j]=hold*r[k][j]+r[i][j];
1543 for(j=0; j<kn; j++) {
1544 if(j == k) continue;
1545 r[k][j]=r[k][j]/big;
1553 for(k=kn-1; k>=0; k--) {
1556 for (j=0; j<kn; j++) {
1564 for (i=0; i<kn; i++) {
1571 ////////////////////////////////////////////////////////////////////
1572 Double_t k00=cin00*r[0][0]+cin10*r[1][0]+cin20*r[2][0]+cin30*r[3][0];
1573 Double_t k01=cin00*r[1][0]+cin10*r[1][1]+cin20*r[2][1]+cin30*r[3][1];
1574 Double_t k02=cin00*r[2][0]+cin10*r[2][1]+cin20*r[2][2]+cin30*r[3][2];
1575 Double_t k03=cin00*r[3][0]+cin10*r[3][1]+cin20*r[3][2]+cin30*r[3][3];
1576 Double_t k10=cin10*r[0][0]+cin11*r[1][0]+cin21*r[2][0]+cin31*r[3][0];
1577 Double_t k11=cin10*r[1][0]+cin11*r[1][1]+cin21*r[2][1]+cin31*r[3][1];
1578 Double_t k12=cin10*r[2][0]+cin11*r[2][1]+cin21*r[2][2]+cin31*r[3][2];
1579 Double_t k13=cin10*r[3][0]+cin11*r[3][1]+cin21*r[3][2]+cin31*r[3][3];
1580 Double_t k20=cin20*r[0][0]+cin21*r[1][0]+cin22*r[2][0]+cin32*r[3][0];
1581 Double_t k21=cin20*r[1][0]+cin21*r[1][1]+cin22*r[2][1]+cin32*r[3][1];
1582 Double_t k22=cin20*r[2][0]+cin21*r[2][1]+cin22*r[2][2]+cin32*r[3][2];
1583 Double_t k23=cin20*r[3][0]+cin21*r[3][1]+cin22*r[3][2]+cin32*r[3][3];
1584 Double_t k30=cin30*r[0][0]+cin31*r[1][0]+cin32*r[2][0]+cin33*r[3][0];
1585 Double_t k31=cin30*r[1][0]+cin31*r[1][1]+cin32*r[2][1]+cin33*r[3][1];
1586 Double_t k32=cin30*r[2][0]+cin31*r[2][1]+cin32*r[2][2]+cin33*r[3][2];
1587 Double_t k33=cin30*r[3][0]+cin31*r[3][1]+cin32*r[3][2]+cin33*r[3][3];
1588 Double_t k40=cin40*r[0][0]+cin41*r[1][0]+cin42*r[2][0]+cin43*r[3][0];
1589 Double_t k41=cin40*r[1][0]+cin41*r[1][1]+cin42*r[2][1]+cin43*r[3][1];
1590 Double_t k42=cin40*r[2][0]+cin41*r[2][1]+cin42*r[2][2]+cin43*r[3][2];
1591 Double_t k43=cin40*r[3][0]+cin41*r[3][1]+cin42*r[3][2]+cin43*r[3][3];
1593 Double_t x0,x1,x2,x3,x4;
1594 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1595 Double_t savex0=x0, savex1=x1, savex2=x2, savex3=x3;
1596 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1)+k02*(m[2]-savex2)+
1598 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1)+k12*(m[2]-savex2)+
1600 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1)+k22*(m[2]-savex2)+
1602 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1)+k32*(m[2]-savex2)+
1604 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1)+k42*(m[2]-savex2)+
1606 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1607 c00=cin00-k00*cin00-k01*cin10-k02*cin20-k03*cin30;
1608 c10=cin10-k00*cin10-k01*cin11-k02*cin21-k03*cin31;
1609 c20=cin20-k00*cin20-k01*cin21-k02*cin22-k03*cin32;
1610 c30=cin30-k00*cin30-k01*cin31-k02*cin32-k03*cin33;
1611 c40=cin40-k00*cin40-k01*cin41-k02*cin42-k03*cin43;
1612 c11=cin11-k10*cin10-k11*cin11-k12*cin21-k13*cin31;
1613 c21=cin21-k10*cin20-k11*cin21-k12*cin22-k13*cin32;
1614 c31=cin31-k10*cin30-k11*cin31-k12*cin32-k13*cin33;
1615 c41=cin41-k10*cin40-k11*cin41-k12*cin42-k13*cin43;
1616 c22=cin22-k20*cin20-k21*cin21-k22*cin22-k23*cin32;
1617 c32=cin32-k20*cin30-k21*cin31-k22*cin32-k23*cin33;
1618 c42=cin42-k20*cin40-k21*cin41-k22*cin42-k23*cin43;
1619 c33=cin33-k30*cin30-k31*cin31-k32*cin32-k33*cin33;
1620 c43=cin43-k30*cin40-k31*cin41-k32*cin42-k33*cin43;
1621 c44=cin44-k40*cin40-k41*cin41-k42*cin42-k43*cin43;
1623 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1624 newTrack->PutCElements(c00,
1628 c40,c41,c42,c43,c44); // put in track the
1631 vmc[0][0]=v00-c00; vmc[1][0]=-c10; vmc[2][0]=-c20; vmc[3][0]=-c30;
1632 vmc[1][1]=v11-c11; vmc[2][1]=-c21; vmc[3][1]=v31-c31;
1633 vmc[2][2]=v22-c22; vmc[3][2]=v32-c32;
1635 vmc[0][1]=vmc[1][0]; vmc[0][2]=vmc[2][0]; vmc[0][3]=vmc[3][0];
1636 vmc[1][2]=vmc[2][1]; vmc[1][3]=vmc[3][1];
1637 vmc[2][3]=vmc[3][2];
1638 /////////////////////// vmc matrix inversion ///////////////////////
1640 for(k=0; k<kn; k++) {
1644 for(j=k; j<kn ; j++) {
1645 for (i=j; i<kn; i++) {
1646 if(TMath::Abs(big) < TMath::Abs(vmc[i][j]) ) {
1656 for(i=0; i<kn; i++) {
1658 vmc[k][i]=vmc[j][i];
1665 for(j=0; j<kn; j++) {
1667 vmc[j][k]=vmc[j][i];
1674 cout << "Singular matrix\n";
1676 for(i=0; i<kn; i++) {
1677 if(i == k) continue;
1678 vmc[i][k]=vmc[i][k]/(-big);
1681 for(i=0; i<kn; i++) {
1683 for(j=0; j<kn; j++) {
1684 if(i == k || j == k) continue;
1685 vmc[i][j]=hold*vmc[k][j]+vmc[i][j];
1689 for(j=0; j<kn; j++) {
1690 if(j == k) continue;
1691 vmc[k][j]=vmc[k][j]/big;
1699 for(k=kn-1; k>=0; k--) {
1702 for (j=0; j<kn; j++) {
1704 vmc[j][k]=-vmc[j][i];
1710 for (i=0; i<kn; i++) {
1712 vmc[k][i]=-vmc[j][i];
1717 ////////////////////////////////////////////////////////////////////
1718 Double_t chi2=(m[0]-x0)*( vmc[0][0]*(m[0]-x0) + 2.*vmc[1][0]*(m[1]-x1) +
1719 2.*vmc[2][0]*(m[2]-x2)+ 2.*vmc[3][0]*(m[3]-x3) )+
1720 (m[1]-x1)* ( vmc[1][1]*(m[1]-x1) + 2.*vmc[2][1]*(m[2]-x2)+
1721 2.*vmc[3][1]*(m[3]-x3) ) +
1722 (m[2]-x2)* ( vmc[2][2]*(m[2]-x2)+ 2.*vmc[3][2]*(m[3]-x3) ) +
1723 (m[3]-x3)*vmc[3][3]*(m[3]-x3);
1724 //cout<<" chi2 kalman = "<<chi2<<"\n"; getchar();
1725 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1726 // newTrack->SetChi2(newTrack->GetChi2()+chi2pred);