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.23 2002/10/22 18:29:34 barbera
19 Tracking V1 ported to the HEAD
21 Revision 1.22 2002/10/22 14:45:36 alibrary
22 Introducing Riostream.h
24 Revision 1.21 2002/02/05 09:12:26 hristov
25 Small mods for gcc 3.02
27 Revision 1.20 2001/11/21 14:47:45 barbera
28 Some unuseful print-out commented out
30 Revision 1.19 2001/11/21 10:49:07 barbera
31 Bug correction suggested by Rene done
33 Revision 1.18 2001/11/20 15:46:17 barbera
34 Point coordinated are calculated in cylindrical reference frame once and for all at the beginning of tracking V1
36 Revision 1.10.2.1 2001/10/24 07:26:04 hristov
37 All the changes from the head are merged with the release
39 Revision 1.14 2001/10/24 07:19:57 hristov
40 Some pointer correctly initialised in one of the constructors
42 Revision 1.13 2001/10/21 19:17:12 hristov
43 Several pointers were set to zero in the default constructors to avoid memory management problems
45 Revision 1.12 2001/10/19 21:32:35 nilsen
46 Minor changes to remove compliation warning on gcc 2.92.2 compiler, and
47 cleanded up a little bit of code.
50 // The purpose of this class is to permorm the ITS tracking. The
51 // constructor has the task to inizialize some private members. The method
52 // DoTracking is written to be called by a macro. It gets the event number,
53 // the minimum and maximum order number of TPC tracks that are to be tracked
54 // trough the ITS, and the file where the recpoints are registered. The
55 // method Recursivetracking is a recursive function that performs the
56 // tracking trough the ITS The method Intersection found the layer, ladder
57 // and detector whre the intersection take place and caluclate the
58 // cohordinates of this intersection. It returns an integer that is 0 if the
59 // intersection has been found successfully. The two mwthods Kalmanfilter
60 // and kalmanfiltervert operate the kalmanfilter without and with the vertex
61 // imposition respectively. The authors thank Mariana Bondila to have help
62 // them to resolve some problems. July-2000
64 #include <Riostream.h>
65 #include <Riostream.h>
71 #include <TStopwatch.h>
73 #include "TParticle.h"
76 #include "AliITSsegmentationSSD.h"
77 #include "AliITSgeomSPD.h"
78 #include "AliITSgeomSDD.h"
79 #include "AliITSgeomSSD.h"
80 #include "AliITSgeom.h"
81 #include "AliITSRecPoint.h"
83 #include "AliKalmanTrack.h"
85 #include "AliITSTrackV1.h"
86 #include "AliITSIOTrack.h"
87 #include "AliITSRad.h"
88 #include "../TPC/AliTPCtracker.h"
89 #include "AliITSTrackerV1.h"
90 #include "AliITSVertex.h"
92 ClassImp(AliITSTrackerV1)
93 //______________________________________________________________________
94 AliITSTrackerV1::AliITSTrackerV1() {
107 for(ia=0; ia<6; ia++) {
125 //______________________________________________________________________
126 AliITSTrackerV1::AliITSTrackerV1(AliITS* IITTSS, Int_t evnumber, Bool_t flag) {
127 //Origin A. Badala' and G.S. Pappalardo:
128 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
129 // Class constructor. It does some initializations.
131 //PH Initialisation taken from the default constructor
146 Int_t imax = 200,jmax = 450;
147 frl = new AliITSRad(imax,jmax);
149 ////////// gets information on geometry /////////////////////////////
150 AliITSgeom *g1 = fITS->GetITSgeom();
155 for(ia=0; ia<6; ia++) {
156 fNlad[ia]=g1->GetNladders(ia+1);
157 fNdet[ia]=g1->GetNdetectors(ia+1);
158 //cout<<fNlad[i]<<" "<<fNdet[i]<<"\n";
161 //cout<<" mean radius = ";
163 for(ib=0; ib<6; ib++) {
164 g1->GetCenterThetaPhi(ib+1,ll,dd,det);
165 Double_t r1=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
166 g1->GetCenterThetaPhi(ib+1,ll,dd+1,det);
167 Double_t r2=TMath::Sqrt(det(0)*det(0)+det(1)*det(1));
168 fAvrad[ib]=(r1+r2)/2.;
169 //cout<<fAvrad[ib]<<" ";
171 //cout<<"\n"; getchar();
173 fDetx[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDx();
174 fDetz[0] = ((AliITSgeomSPD*)(g1->GetShape(1, ll, dd)))->GetDz();
176 fDetx[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDx();
177 fDetz[1] = ((AliITSgeomSPD*)(g1->GetShape(2, ll, dd)))->GetDz();
179 fDetx[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDx();
180 fDetz[2] = ((AliITSgeomSDD*)(g1->GetShape(3, ll, dd)))->GetDz();
182 fDetx[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDx();
183 fDetz[3] = ((AliITSgeomSDD*)(g1->GetShape(4, ll, dd)))->GetDz();
185 fDetx[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDx();
186 fDetz[4] = ((AliITSgeomSSD*)(g1->GetShape(5, ll, dd)))->GetDz();
188 fDetx[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDx();
189 fDetz[5] = ((AliITSgeomSSD*)(g1->GetShape(6, ll, dd)))->GetDz();
190 //cout<<" Detx Detz\n";
191 //for(Int_t la=0; la<6; la++) cout<<" "<<fDetx[la]<<" "<<
195 // allocate memory and define matrices fzmin, fzmax, fphimin and fphimax //
197 Double_t epszdrift=0.05;
199 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
200 Int_t im1, im2, im2max;
201 for(im1=0; im1<6; im1++) {
203 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
206 for(im1=0; im1<6; im1++) {
208 for(im2=0; im2<im2max; im2++) {
209 g1->GetCenterThetaPhi(im1+1,1,im2+1,det);
210 if(im2!=0) fzmin[im1][im2]=det(2)-fDetz[im1];
212 fzmin[im1][im2]=det(2)-(fDetz[im1])*epsz;
213 if(im2!=(im2max-1)) fzmax[im1][im2]=det(2)+fDetz[im1];
215 fzmax[im1][im2]=det(2)+fDetz[im1]*epsz;
216 if(im1==2 || im1==3) {
217 fzmin[im1][im2]-=epszdrift;
218 fzmax[im1][im2]+=epszdrift;
219 } // end if im1==2 || im1==3
223 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
224 for(im1=0;im1<6;im1++) {
226 fphimin[im1] = new Double_t[im2max];
227 fphimax[im1] = new Double_t[im2max];
230 fphidet = new Double_t*[6];
231 for(im1=0; im1<6; im1++) {
233 fphidet[im1] = new Double_t[im2max];
236 //Float_t global[3],local[3];
237 Double_t global[3],local[3];
238 Double_t pigre=TMath::Pi();
239 Double_t xmin,ymin,xmax,ymax;
241 for(im1=0; im1<6; im1++) {
243 for(im2=0; im2<im2max; im2++) {
245 g1->GetCenterThetaPhi(im1+1,im2+1,idet,det);
246 fphidet[im1][im2] = TMath::ATan2(Double_t(det(1)),
248 if(fphidet[im1][im2]<0.) fphidet[im1][im2]+=2.*pigre;
249 local[1]=local[2]=0.;
250 local[0]= -(fDetx[im1]);
251 if(im1==0) local[0]= (fDetx[im1]); //to take into account
252 // different reference system
253 g1->LtoG(im1+1,im2+1,idet,local,global);
254 xmax=global[0]; ymax=global[1];
255 local[0]= (fDetx[im1]);
256 if(im1==0) local[0]= -(fDetx[im1]);//take into account different
258 g1->LtoG(im1+1,im2+1,idet,local,global);
259 xmin=global[0]; ymin=global[1];
260 fphimin[im1][im2]= TMath::ATan2(ymin,xmin);
261 if(fphimin[im1][im2]<0.) fphimin[im1][im2]+=2.*pigre;
262 fphimax[im1][im2]= TMath::ATan2(ymax,xmax);
263 if(fphimax[im1][im2]<0.) fphimax[im1][im2]+=2.*pigre;
266 //////////////////////////////////////////////////////////////////////////////////////////////////////////
267 /////////////// allocate memory and define vector fNRecPoints and matrices fRecCylR, fRecCylPhi, fRecCylZ /////////////
268 gAlice->GetEvent(evnumber);
269 Int_t NumOfModules = g1->GetIndexMax();
270 //fRecCylR = new Float_t *[NumOfModules];
271 fRecCylR = new Double_t *[NumOfModules];
272 //fRecCylPhi = new Float_t *[NumOfModules];
273 fRecCylPhi = new Double_t *[NumOfModules];
274 //fRecCylZ = new Float_t *[NumOfModules];
275 fRecCylZ = new Double_t *[NumOfModules];
276 AliITSRecPoint *recp;
277 fNRecPoints = new Int_t[NumOfModules];
279 for(Int_t module=0; module<NumOfModules; module++) {
280 fITS->ResetRecPoints();
281 gAlice->TreeR()->GetEvent(module);
282 frecPoints=fITS->RecPoints();
283 Int_t nRecPoints=fNRecPoints[module]=frecPoints->GetEntries();
285 fRecCylR[module] = new Float_t[nRecPoints];
286 fRecCylPhi[module] = new Float_t[nRecPoints];
287 fRecCylZ[module] = new Float_t[nRecPoints];
289 fRecCylR[module] = new Double_t[nRecPoints];
290 fRecCylPhi[module] = new Double_t[nRecPoints];
291 fRecCylZ[module] = new Double_t[nRecPoints];
293 for(ind=0; ind<fNRecPoints[module]; ind++) {
294 recp=(AliITSRecPoint*)frecPoints->UncheckedAt(ind);
295 // Float_t global[3], local[3];
296 Double_t global[3], local[3];
297 local[0]=recp->GetX();
299 local[2]= recp->GetZ();
300 g1->LtoG(module,local,global);
302 Float_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]); // r hit
303 Float_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi(); // phi hit
304 Float_t z = global[2]; // z hit
307 Double_t r = TMath::Sqrt(global[0]*global[0]+global[1]*global[1]); // r hit
308 Double_t phi = TMath::ATan2(global[1],global[0]); if(phi<0.) phi+=2.*TMath::Pi(); // phi hit
309 Double_t z = global[2]; // z hit
311 fRecCylR[module][ind]=r;
312 fRecCylPhi[module][ind]=phi;
313 fRecCylZ[module][ind]=z;
318 ////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
321 ////////// gets magnetic field factor //////////////////////////////
323 AliMagF * fieldPointer = gAlice->Field();
324 // fFieldFactor = (Double_t)fieldPointer->Factor();
325 fFieldFactor =(Double_t)fieldPointer-> SolenoidField()/10/.2;
326 // cout<< " field factor = "<<fFieldFactor<<"\n"; getchar();
328 //______________________________________________________________________
329 AliITSTrackerV1::AliITSTrackerV1(const AliITSTrackerV1 &cobj) {
330 // Origin A. Badala' and G.S. Pappalardo:
331 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
335 *fresult = *cobj.fresult;
336 fPtref = cobj.fPtref;
337 fChi2max = cobj.fChi2max;
338 **fvettid = **cobj.fvettid;
339 fflagvert = cobj.fflagvert;
340 Int_t imax=200,jmax=450;
341 frl = new AliITSRad(imax,jmax);
343 fFieldFactor = cobj.fFieldFactor;
344 Int_t i,im1,im2,im2max;
346 fNlad[i] = cobj.fNlad[i];
347 fNdet[i] = cobj.fNdet[i];
348 fAvrad[i] = cobj.fAvrad[i];
349 fDetx[i] = cobj.fDetx[i];
350 fDetz[i] = cobj.fDetz[i];
352 fzmin = new Double_t*[6]; fzmax = new Double_t*[6];
353 for(im1=0; im1<6; im1++) {
355 fzmin[im1] = new Double_t[im2max];
356 fzmax[im1] = new Double_t[im2max];
358 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
359 for(im1=0;im1<6;im1++) {
361 fphimin[im1] = new Double_t[im2max];
362 fphimax[im1] = new Double_t[im2max];
365 fphidet = new Double_t*[6];
366 for(im1=0; im1<6; im1++) {
368 fphidet[im1] = new Double_t[im2max];
370 for(im1=0; im1<6; im1++) {
372 for(im2=0; im2<im2max; im2++) {
373 fzmin[im1][im2]=cobj.fzmin[im1][im2];
374 fzmax[im1][im2]=cobj.fzmax[im1][im2];
377 for(im1=0; im1<6; im1++) {
379 for(im2=0; im2<im2max; im2++) {
380 fphimin[im1][im2]=cobj.fphimin[im1][im2];
381 fphimax[im1][im2]=cobj.fphimax[im1][im2];
382 fphidet[im1][im2]=cobj.fphidet[im1][im2];
387 AliITSgeom *g1 = fITS->GetITSgeom();
388 Int_t NumOfModules = g1->GetIndexMax();
390 fRecCylR = new Float_t *[NumOfModules];
391 fRecCylPhi = new Float_t *[NumOfModules];
392 fRecCylZ = new Float_t *[NumOfModules];
394 fRecCylR = new Double_t *[NumOfModules];
395 fRecCylPhi = new Double_t *[NumOfModules];
396 fRecCylZ = new Double_t *[NumOfModules];
397 fNRecPoints = new Int_t[NumOfModules];
398 for(Int_t module=0; module<NumOfModules; module++) {
399 Int_t nRecPoints=fNRecPoints[module]=cobj.fNRecPoints[module];
401 fRecCylR[module] = new Float_t[nRecPoints];
402 fRecCylPhi[module] = new Float_t[nRecPoints];
403 fRecCylZ[module] = new Float_t[nRecPoints];
405 fRecCylR[module] = new Double_t[nRecPoints];
406 fRecCylPhi[module] = new Double_t[nRecPoints];
407 fRecCylZ[module] = new Double_t[nRecPoints];
409 for(ind=0; ind<nRecPoints; ind++) {
410 fRecCylR[module][ind]=cobj.fRecCylR[module][ind];
411 fRecCylPhi[module][ind]=cobj.fRecCylPhi[module][ind];
412 fRecCylZ[module][ind]=cobj.fRecCylZ[module][ind];
417 void AliITSTrackerV1::DelMatrix(Int_t NumOfModules) {
418 for(Int_t mod=0; mod<NumOfModules; mod++) {
419 delete fRecCylR[mod];
420 delete fRecCylPhi[mod];
421 delete fRecCylZ[mod];
427 //______________________________________________________________________
428 AliITSTrackerV1::~AliITSTrackerV1(){
429 // Origin A. Badala' and G.S. Pappalardo:
430 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
434 for(Int_t i=0; i<6; i++) {
449 //______________________________________________________________________
450 AliITSTrackerV1 &AliITSTrackerV1::operator=(AliITSTrackerV1 obj) {
451 // Origin A. Badala' and G.S. Pappalardo:
452 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
453 // assignement operator
456 *fresult = *obj.fresult;
458 fChi2max = obj.fChi2max;
459 **fvettid = **obj.fvettid;
460 fflagvert = obj.fflagvert;
461 Int_t imax=200,jmax=450;
462 frl = new AliITSRad(imax,jmax);
464 fFieldFactor = obj.fFieldFactor;
467 fNlad[i] = obj.fNlad[i];
468 fNdet[i] = obj.fNdet[i];
469 fAvrad[i] = obj.fAvrad[i];
470 fDetx[i] = obj.fDetx[i];
471 fDetz[i] = obj.fDetz[i];
473 fzmin = new Double_t*[6];
474 fzmax = new Double_t*[6];
475 Int_t im1, im2, im2max;
476 for(im1=0; im1<6; im1++) {
478 fzmin[im1] = new Double_t[im2max]; fzmax[im1] = new Double_t[im2max];
480 fphimin = new Double_t*[6]; fphimax = new Double_t*[6];
481 for(im1=0;im1<6;im1++) {
483 fphimin[im1] = new Double_t[im2max];
484 fphimax[im1] = new Double_t[im2max];
487 fphidet = new Double_t*[6];
488 for(im1=0; im1<6; im1++) {
490 fphidet[im1] = new Double_t[im2max];
492 for(im1=0; im1<6; im1++) {
494 for(im2=0; im2<im2max; im2++) {
495 fzmin[im1][im2]=obj.fzmin[im1][im2];
496 fzmax[im1][im2]=obj.fzmax[im1][im2];
499 for(im1=0; im1<6; im1++) {
501 for(im2=0; im2<im2max; im2++) {
502 fphimin[im1][im2]=obj.fphimin[im1][im2];
503 fphimax[im1][im2]=obj.fphimax[im1][im2];
504 fphidet[im1][im2]=obj.fphidet[im1][im2];
508 AliITSgeom *g1 = fITS->GetITSgeom();
509 Int_t NumOfModules = g1->GetIndexMax();
511 fRecCylR = new Float_t *[NumOfModules];
512 fRecCylPhi = new Float_t *[NumOfModules];
513 fRecCylZ = new Float_t *[NumOfModules];
515 fRecCylR = new Double_t *[NumOfModules];
516 fRecCylPhi = new Double_t *[NumOfModules];
517 fRecCylZ = new Double_t *[NumOfModules];
518 fNRecPoints = new Int_t[NumOfModules];
519 for(Int_t module=0; module<NumOfModules; module++) {
520 Int_t nRecPoints=fNRecPoints[module]=obj.fNRecPoints[module];
522 fRecCylR[module] = new Float_t[nRecPoints];
523 fRecCylPhi[module] = new Float_t[nRecPoints];
524 fRecCylZ[module] = new Float_t[nRecPoints];
526 fRecCylR[module] = new Double_t[nRecPoints];
527 fRecCylPhi[module] = new Double_t[nRecPoints];
528 fRecCylZ[module] = new Double_t[nRecPoints];
530 for(ind=0; ind<nRecPoints; ind++) {
531 fRecCylR[module][ind]=obj.fRecCylR[module][ind];
532 fRecCylPhi[module][ind]=obj.fRecCylPhi[module][ind];
533 fRecCylZ[module][ind]=obj.fRecCylZ[module][ind];
540 //______________________________________________________________________
541 void AliITSTrackerV1::DoTracking(Int_t evNumber,Int_t minTr,Int_t maxTr,
542 TFile *file, Bool_t realmass) {
543 // Origin A. Badala' and G.S. Pappalardo:
544 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
545 // The method needs the event number, the minimum and maximum order
546 // number of TPC tracks that
547 // are to be tracked trough the ITS, and the file where the recpoints
549 // The method can be called by a macro. It preforms the tracking for
550 // all good TPC tracks
552 printf("begin DoTracking - file %p\n",file);
554 gAlice->GetEvent(evNumber); //modificato per gestire hbt
556 AliKalmanTrack *kkprov;
557 //kkprov->SetConvConst(100/0.299792458/0.2/fFieldFactor);
558 kkprov->SetConvConst(1000/0.299792458/gAlice->Field()->SolenoidField());
559 // cout<<" field = "<<gAlice->Field()->SolenoidField()<<endl;
562 TFile *cf=TFile::Open("AliTPCclusters.root");
563 AliTPCParam *digp= (AliTPCParam*)cf->Get("75x40_100x60_150x60");
564 if (!digp) { cerr<<"TPC parameters have not been found !\n"; getchar();}
567 AliTPCtracker *tracker = new AliTPCtracker(digp,evNumber);
570 tracker->LoadInnerSectors();
571 tracker->LoadOuterSectors();
574 TFile *tf=TFile::Open("AliTPCtracksSorted.root");
576 cerr<<"Can't open AliTPCtracksSorted.root !\n";
579 TObjArray tracks(200000);
581 sprintf(tname,"TreeT_TPC_%d",evNumber);
583 TTree *tracktree=(TTree*)tf->Get(tname);
584 if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n";}
585 TBranch *tbranch=tracktree->GetBranch("tracks");
586 Int_t nentr=(Int_t)tracktree->GetEntries();
589 AliTPCtrack *ioTrackTPC=0;
590 for (kk=0; kk<nentr; kk++) {
591 ioTrackTPC=new AliTPCtrack;
592 tbranch->SetAddress(&ioTrackTPC);
593 tracktree->GetEvent(kk);
594 tracker->CookLabel(ioTrackTPC,0.1);
595 tracks.AddLast(ioTrackTPC);
600 Int_t nt = tracks.GetEntriesFast();
601 cerr<<"Number of found tracks "<<nt<<endl;
604 TTree *tr=gAlice->TreeR();
605 Int_t nent=(Int_t)tr->GetEntries();
606 frecPoints = fITS->RecPoints();
610 Int_t *np = new Int_t[nent];
611 fvettid = new Int_t* [nent];
614 for (mod=0; mod<nent; mod++) {
616 fITS->ResetRecPoints();
617 gAlice->TreeR()->GetEvent(mod);
618 numbpoints = frecPoints->GetEntries();
619 totalpoints+=numbpoints;
620 np[mod] = numbpoints;
621 //cout<<" mod = "<<mod<<" numbpoints = "<<numbpoints<<"\n";getchar();
622 fvettid[mod] = new Int_t[numbpoints];
624 for (ii=0;ii<numbpoints; ii++) *(fvettid[mod]+ii)=0;
627 AliTPCtrack *track=0;
629 if(minTr < 0) {minTr = 0; maxTr = nt-1;}
633 TTree tracktree1("TreeT","Tree with ITS tracks");
634 AliITSIOTrack *ioTrack=0;
635 tracktree1.Branch("ITStracks","AliITSIOTrack",&ioTrack,32000,0);
637 TDatabasePDG * db = new TDatabasePDG;
640 for (j=minTr; j<=maxTr; j++) {
641 track=(AliTPCtrack*)tracks.UncheckedAt(j);
642 if (!track) continue;
644 /// mass definition ////////////////////////
645 Double_t mass=0.13956995;
646 Int_t pcode=211; // a pion by default
648 Int_t TPClabel=TMath::Abs( track->GetLabel() );
649 TParticle *p = (TParticle*)gAlice->Particle(TPClabel);
650 pcode=p->GetPdgCode();
651 // Int_t mothercode=p->GetFirstMother();
652 //if(mothercode>0 ) numofsecondaries++; else numofprimaries++;
654 //if(!pcode) pcode=211;
655 if(TMath::Abs(pcode)<20443) mass=db->GetParticle(pcode)->Mass();
658 ///////////////////////////////////////////////
660 ////// propagation to the end of TPC //////////////
662 track->PropagateTo(xk, 28.94, 1.204e-3,mass); //Ne
664 track->PropagateTo(xk, 44.77, 1.71,mass); //Tedlar
666 track->PropagateTo(xk, 44.86, 1.45,mass); //kevlar
668 track->PropagateTo(xk, 41.28, 0.029,mass); //Nomex
670 track->PropagateTo(xk,36.2,1.98e-3,mass); //C02
672 track->PropagateTo(xk, 24.01, 2.7,mass); //Al
674 track->PropagateTo(xk, 44.77, 1.71,mass); //Tedlar
676 track->PropagateTo(xk, 44.86, 1.45,mass); //kevlar
678 track->PropagateTo(xk, 41.28, 0.029,mass); //Nomex
679 ////////////////////////////////////////////////////////////////////
681 // new propagation to the end of TPC
683 // track->PropagateTo(xk,0.,0.); //Ne if it's still there //attenzione funziona solo se modifica in TPC
684 // Double_t xk=77.415;
685 track->PropagateTo(xk, 28.94, 1.204e-3);
687 track->PropagateTo(xk, 44.77,1.71); //Tedlar
689 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
691 track->PropagateTo(xk, 41.28, 0.029);//Nomex
693 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
695 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
698 // track->PropagateTo(xk,0.,0.); //C02
699 track->PropagateTo(xk,36.2,1.98e-3); //C02 //attenzione funziona solo se modifica in TPC
702 track->PropagateTo(xk, 24.01, 2.7); //Al
704 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
706 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
708 track->PropagateTo(xk, 41.28, 0.029); //Nomex
710 track->PropagateTo(xk, 44.86, 1.45); //Kevlar
712 track->PropagateTo(xk, 44.77, 1.71); //Tedlar
714 track->PropagateTo(xk, 24.01, 2.7); //Al
716 ////////////////////////////////////////////////////////////////////////////////////////////////////////
717 //AliITSTrackV1 trackITS(*track);
718 AliITSTrackV1 trackITS(*track, fFieldFactor);
719 //cout<<" fFieldFactor = "<<fFieldFactor<<"\n";
720 trackITS.PutMass(mass); //new to add mass to track
721 if(fresult){ delete fresult; fresult=0;}
722 fresult = new AliITSTrackV1(trackITS);
724 AliITSTrackV1 primaryTrack(trackITS);
725 vgeant=(*fresult).GetVertex();
727 // Definition of dv and zv for vertex constraint
728 Double_t sigmaDv=0.0050; Double_t sigmaZv=0.010;
729 //Double_t sigmaDv=0.0015; Double_t sigmaZv=0.0015;
730 Double_t uniform= gRandom->Uniform();
732 if(uniform<=0.5) signdv=-1.;
736 Double_t vr=TMath::Sqrt(vgeant(0)*vgeant(0)+ vgeant(1)*vgeant(1));
737 Double_t dv=gRandom->Gaus(signdv*vr,(Float_t)sigmaDv);
738 Double_t zv=gRandom->Gaus(vgeant(2),(Float_t)sigmaZv);
739 //cout<<" Dv e Zv = "<<dv<<" "<<zv<<"\n";
742 trackITS.SetsigmaDv(sigmaDv);
743 trackITS.SetsigmaZv(sigmaZv);
744 (*fresult).SetDv(dv);
745 (*fresult).SetZv(zv);
746 (*fresult).SetsigmaDv(sigmaDv);
747 (*fresult).SetsigmaZv(sigmaZv);
748 primaryTrack.SetDv(dv);
749 primaryTrack.SetZv(zv);
750 primaryTrack.SetsigmaDv(sigmaDv);
751 primaryTrack.SetsigmaZv(sigmaZv);
752 primaryTrack.PrimaryTrack(frl);
753 TVector d2=primaryTrack.Getd2();
754 TVector tgl2=primaryTrack.Gettgl2();
755 TVector dtgl=primaryTrack.Getdtgl();
756 trackITS.Setd2(d2); trackITS.Settgl2(tgl2);
757 trackITS.Setdtgl(dtgl);
758 (*fresult).Setd2(d2); (*fresult).Settgl2(tgl2);
759 (*fresult).Setdtgl(dtgl);
761 trackITS.SetVertex(vertex); trackITS.SetErrorVertex(ervertex);
762 (*result).SetVertex(vertex); (*result).SetErrorVertex(ervertex);
764 TList *list= new TList();
766 list->AddLast(&trackITS);
768 fPtref=TMath::Abs( (trackITS).GetPt() );
769 //cout<<" fPtref = " <<fPtref<<"\n";
770 if(fPtref>1.0) fChi2max=40.;
771 if(fPtref<=1.0) fChi2max=20.;
772 if(fPtref<0.4 ) fChi2max=100.;
773 if(fPtref<0.2 ) fChi2max=40.;
774 // if(fPtref<0.4 ) fChi2max=30.;
775 // if(fPtref<0.2 ) fChi2max=20.;
776 //if(fPtref<0.2 ) fChi2max=10.;
777 //if(fPtref<0.1 ) fChi2max=5.;
778 //cout << "\n Pt = " << fPtref <<"\n"; //stampa
779 RecursiveTracking(list);
784 TVector vecTotLabRef(18);
786 for(lay=5; lay>=0; lay--) {
787 TVector vecLabRef(3);
788 vecLabRef=(*fresult).GetLabTrack(lay);
789 Float_t clustZ=(*fresult).GetZclusterTrack( lay);
791 Int_t lpp=(Int_t)vecLabRef(k);
793 TParticle *p=(TParticle*) gAlice->Particle(lpp);
794 Int_t pcode=p->GetPdgCode();
795 if(pcode==11) vecLabRef(k)=p->GetFirstMother();
797 itot++; vecTotLabRef(itot)=vecLabRef(k);
798 if(vecLabRef(k)==0. && clustZ == -1.) vecTotLabRef(itot) =-3.;
803 (*fresult).Search(vecTotLabRef, labref, freq);
805 //if(freq < 6) labref=-labref; // cinque - sei
806 if(freq < 5) labref=-labref; // cinque - sei
807 (*fresult).SetLabel(labref);
809 // cout<<" progressive track number = "<<j<<"\r";
811 Int_t numOfCluster=(*fresult).GetNumClust();
812 //cout<<" progressive track number = "<<j<<"\n"; // stampa
813 Long_t labITS=(*fresult).GetLabel();
814 //cout << " ITS track label = " << labITS << "\n"; // stampa
815 Int_t lab=track->GetLabel();
816 //cout << " TPC track label = " << lab <<"\n"; // stampa
817 //propagation to vertex
820 if((*fresult).DoNotCross(rbeam)) continue; //no intersection with beampipe
821 (*fresult).Propagation(rbeam);
822 Double_t c00,c10,c11,c20,c21,c22,c30,c31,c32,c33,c40,c41,c42,c43,c44;
823 (*fresult).GetCElements(c00,
827 c40,c41,c42,c43,c44);
829 Double_t pt=TMath::Abs((*fresult).GetPt());
830 Double_t dr=(*fresult).GetD();
831 Double_t z=(*fresult).GetZ();
832 Double_t tgl=(*fresult).GetTgl();
833 Double_t c=(*fresult).GetC();
835 Double_t dz=z-(tgl/cy)*TMath::ASin((*fresult).Arga(rbeam));
837 // cout<<" dr e dz alla fine = "<<dr<<" "<<dz<<"\n"; getchar();
838 Double_t phi=(*fresult).Getphi();
839 Double_t phivertex = phi - TMath::ASin((*fresult).ArgA(rbeam));
840 Double_t duepi=2.*TMath::Pi();
841 if(phivertex>duepi) phivertex-=duepi;
842 if(phivertex<0.) phivertex+=duepi;
843 /////////////////////////////////////////////////////////////
844 Int_t idmodule,idpoint;
845 if(numOfCluster >=5) { // cinque - sei
846 //if(numOfCluster ==6) { // cinque - sei
847 AliITSIOTrack outTrack;
849 ioTrack->SetStatePhi(phi);
850 ioTrack->SetStateZ(z);
851 ioTrack->SetStateD(dr);
852 ioTrack->SetStateTgl(tgl);
853 ioTrack->SetStateC(c);
854 Double_t radius=(*fresult).Getrtrack();
855 ioTrack->SetRadius(radius);
857 if(c>0.) charge=-1; else charge=1;
858 ioTrack->SetCharge(charge);
859 ioTrack->SetCovMatrix(c00,
863 c40,c41,c42,c43,c44);
864 Double_t px=pt*TMath::Cos(phivertex);
865 Double_t py=pt*TMath::Sin(phivertex);
867 Double_t xtrack=dr*TMath::Sin(phivertex);
868 Double_t ytrack=dr*TMath::Cos(phivertex);
869 Double_t ztrack=dz+vgeant(2);
873 ioTrack->SetX(xtrack);
874 ioTrack->SetY(ytrack);
875 ioTrack->SetZ(ztrack);
876 ioTrack->SetLabel(labITS);
877 ioTrack->SetTPCLabel(lab);
880 for(il=0;il<6; il++){
881 ioTrack->SetIdPoint(il,(*fresult).GetIdPoint(il));
882 ioTrack->SetIdModule(il,(*fresult).GetIdModule(il));
885 for (il=0;il<6;il++) {
886 idpoint=(*fresult).GetIdPoint(il);
887 idmodule=(*fresult).GetIdModule(il);
888 if(idmodule>0.) *(fvettid[idmodule]+idpoint)=1;
890 ioTrack->SetIdPoint(il,idpoint);
891 ioTrack->SetIdModule(il,idmodule);
893 } // end if on numOfCluster
894 //gObjectTable->Print(); // stampa memoria
895 } // end for (int j=minTr; j<=maxTr; j++)
897 static Bool_t first=kTRUE;
900 tfile=new TFile("itstracks.root","RECREATE");
901 //cout<<"I have opened itstracks.root file "<<endl;
907 sprintf(hname,"TreeT%d",evNumber);
908 tracktree1.Write(hname);
910 TTree *fAli=gAlice->TreeK();
912 if (fAli) fileAli =fAli->GetCurrentFile();
914 ////////////////////////////////////////////////////////////////////
916 printf("delete vectors\n");
918 if(fvettid) delete [] fvettid;
919 if(fresult) {delete fresult; fresult=0;}
921 //______________________________________________________________________
922 void AliITSTrackerV1::RecursiveTracking(TList *trackITSlist) {
923 // This function perform the recursive tracking in ITS detectors
924 // reference is a pointer to the final best track
925 // Origin A. Badala' and G.S. Pappalardo:
926 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
927 // The authors thank Mariana Bondila to have help them to resolve some
928 // problems. July-2000
930 //Rlayer[0]=4.; Rlayer[1]=7.; Rlayer[2]=14.9;
931 // Rlayer[3]=23.8; Rlayer[4]=39.1; Rlayer[5]=43.6; //vecchio
933 //////////////////////
934 Float_t sigmaphil[6], sigmazl[6];
935 sigmaphil[0]=1.44e-6/(fAvrad[0]*fAvrad[0]);
936 sigmaphil[1]=1.44e-6/(fAvrad[1]*fAvrad[1]);
937 sigmaphil[2]=1.444e-5/(fAvrad[2]*fAvrad[2]);
938 sigmaphil[3]=1.444e-5/(fAvrad[3]*fAvrad[3]);
939 sigmaphil[4]=4e-6/(fAvrad[4]*fAvrad[4]);
940 sigmaphil[5]=4e-6/(fAvrad[5]*fAvrad[5]);
947 ///////////////////////////////////////////////////////////
949 AliITSgeom *g1 = fITS->GetITSgeom();
950 AliITSRecPoint *recp;
951 for(index =0; index<trackITSlist->GetSize(); index++) {
952 AliITSTrackV1 *trackITS = (AliITSTrackV1 *) trackITSlist->At(index);
953 if((*trackITS).GetLayer()==7) fresult->SetChi2(10.223e140);
954 // cout <<" Layer inizio = "<<(*trackITS).GetLayer()<<"\n";
955 // cout<<"fvtrack =" <<"\n";
956 // cout << (*trackITS)(0) << " "<<(*trackITS)(1)<<" "
957 // <<(*trackITS)(2)<<" "<<(*trackITS)(3)<<" "
958 // <<(*trackITS)(4)<<"\n";
959 // cout<< " rtrack = "<<(*trackITS).Getrtrack()<<"\n";
960 // cout<< " Pt = "<<(*trackITS).GetPt()<<"\n";
962 Double_t chi2Now, chi2Ref;
963 Float_t numClustRef = fresult->GetNumClust();
964 if((*trackITS).GetLayer()==1 ) {
965 chi2Now = trackITS->GetChi2();
966 Float_t numClustNow = trackITS->GetNumClust();
967 if(trackITS->GetNumClust())
968 chi2Now /= (Double_t)trackITS->GetNumClust();
969 chi2Ref = fresult->GetChi2();
970 if(fresult->GetNumClust())
971 chi2Ref /= (Double_t)fresult->GetNumClust();
972 //cout<<" chi2Now and chi2Ref = "<<chi2Now<<" "<<chi2Ref<<"\n";
973 if( numClustNow > numClustRef ) {*fresult = *trackITS;}
974 if((numClustNow == numClustRef )&&
975 (chi2Now < chi2Ref)) {
976 *fresult = *trackITS;
981 if(trackITS->Getfnoclust()>=2) continue;
982 Float_t numClustNow = trackITS->GetNumClust();
984 chi2Now = trackITS->GetChi2();
986 if(numClustNow<numClustRef && chi2Now>fresult->GetChi2()) continue;
987 //cout<<" chi2Now = "<<chi2Now<<"\n";
989 chi2Now/=numClustNow;
990 if(fPtref > 1.0 && chi2Now > 30.) continue;
991 if((fPtref >= 0.6 && fPtref<=1.0) && chi2Now > 40.) continue;
992 // if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 40.) continue;
993 // if(fPtref <= 0.2 && chi2Now > 8.) continue;
994 if((fPtref <= 0.6 && fPtref>0.2)&& chi2Now > 30.) continue;
995 if(fPtref <= 0.2 && chi2Now > 7.) continue;
996 /////////////////////////////
999 Int_t layerInit = (*trackITS).GetLayer();
1000 Int_t layernew = layerInit - 2;// -1 for new layer, -1 for matrix index
1002 Int_t ladp, ladm, detp,detm,ladinters,detinters;
1003 Int_t layerfin=layerInit-1;
1004 // cout<<"Prima di intersection \n";
1005 Int_t outinters=Intersection(*trackITS,layerfin,ladinters,detinters);
1006 // cout<<" outinters = "<<outinters<<"\n";
1007 // cout<<" Layer ladder detector intersection ="
1008 // <<layerfin<<" "<<ladinters<<" "<<detinters<<"\n";
1009 // cout << " phiinters zinters = "<<(*trackITS)(0)
1010 // << " "<<(*trackITS)(1)<<"\n"; getchar();
1011 if(outinters==-1) continue;
1014 TVector toucLad(9), toucDet(9);
1015 Int_t lycur=layerfin;
1018 if(ladm <= 0) ladm=fNlad[layerfin-1];
1019 if(ladp > fNlad[layerfin-1]) ladp=1;
1024 toucLad(0)=ladinters; toucLad(1)=ladm; toucLad(2)=ladp;
1025 toucLad(3)=ladinters; toucLad(4)=ladm; toucLad(5)=ladp;
1026 toucLad(6)=ladinters; toucLad(7)=ladm; toucLad(8)=ladp;
1027 toucDet(0)=detinters; toucDet(1)=detinters; toucDet(2)=detinters;
1028 if(detm > 0 && detp <= fNdet[layerfin-1]) {
1030 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1031 toucDet(6)=detp; toucDet(7)=detp; toucDet(8)=detp;
1033 if(detm > 0 && detp > fNdet[layerfin-1]) {
1035 toucDet(3)=detm; toucDet(4)=detm; toucDet(5)=detm;
1037 if(detm <= 0 && detp <= fNdet[layerfin-1]) {
1039 toucDet(3)=detp; toucDet(4)=detp; toucDet(5)=detp;
1042 Float_t epsphi=5.0, epsz=5.0;
1043 if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1044 // new definition of idetot e toucLad e toucDet to be
1045 // transformed in a method
1046 // these values could be modified
1047 Float_t pigre=TMath::Pi();
1048 Float_t rangephi=5., rangez=5.;
1049 if(layerfin==1 || layerfin ==2){
1050 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1051 (*trackITS).GetSigmaphi());
1052 rangez = 40.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1053 (*trackITS).GetSigmaZ());
1055 if(layerfin==3 || layerfin ==4){
1056 //rangephi=30.*fepsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1057 // (*trackITS).GetSigmaphi());
1058 //rangez = 40.*fepsz*TMath::Sqrt(sigmazl[layerfin-1]+
1059 // (*trackITS).GetSigmaZ());
1060 rangephi=40.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1061 (*trackITS).GetSigmaphi());
1062 rangez = 50.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1063 (*trackITS).GetSigmaZ());
1065 if(layerfin==5 || layerfin ==6){
1066 rangephi=20.*epsphi*TMath::Sqrt(sigmaphil[layerfin-1]+
1067 (*trackITS).GetSigmaphi());
1068 rangez =5.*epsz*TMath::Sqrt(sigmazl[layerfin-1]+
1069 (*trackITS).GetSigmaZ());
1071 Float_t phinters, zinters;
1072 phinters=(*trackITS).Getphi();
1073 zinters=(*trackITS).GetZ();
1074 Float_t distz = 0.0;
1075 Float_t phicm, phicp, distphim, distphip;
1077 if(phinters>fphimax[layerfin-1][ladm-1]) phicm=phinters-2*pigre; //corretto il 20-11-2001
1078 distphim=TMath::Abs(phicm-fphimax[layerfin-1][ladm-1]); //corretto il 20-11-2001
1080 //cout<<" fNlad[layerfin-1] e ladp = "<<fNlad[layerfin-1]<<" "<<ladp<<endl;
1081 if(phinters>fphimin[layerfin-1][ladp-1]) phicp=phinters-2.*pigre; //corretto il 20-11-2001
1082 distphip=TMath::Abs(phicp-fphimin[layerfin-1][ladp-1]); //corretto il 20-11-2001
1086 toucLad(0)=ladinters; toucDet(0)=detinters;
1087 if(detm>0) distz=TMath::Abs(zinters-fzmax[layerfin-1][detm-1]);
1088 if(detm>0 && rangez>=distz){
1090 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detm;
1091 if(rangephi>=distphim){
1093 toucLad(idetot-1)=ladm;
1094 toucDet(idetot-1)=detinters;
1096 toucLad(idetot-1)=ladm;
1097 toucDet(idetot-1)=detm;
1099 if(rangephi>=distphip){
1101 toucLad(idetot-1)=ladp;
1102 toucDet(idetot-1)=detinters;
1104 toucLad(idetot-1)=ladp;
1105 toucDet(idetot-1)=detm;
1108 if(detp<=fNdet[layerfin-1])
1109 distz=TMath::Abs(zinters-fzmin[layerfin-1][detp-1]);
1110 if(detp<=fNdet[layerfin-1] && rangez>=distz){
1112 idetot++; toucLad(idetot-1)=ladinters; toucDet(idetot-1)=detp;
1113 if(rangephi>=distphim){
1114 idetot++; toucLad(idetot-1)=ladm; toucDet(idetot-1)=detp;
1117 toucLad(idetot-1)=ladm;
1118 toucDet(idetot-1)=detinters;
1121 if(rangephi>=distphip){
1123 toucLad(idetot-1)=ladp;
1124 toucDet(idetot-1)=detp;
1127 toucLad(idetot-1)=ladp;
1128 toucDet(idetot-1)=detinters;
1131 } //end detm<fNdet[.......
1133 if(flagzmin == 0 && flagzmax==0){
1134 if(rangephi>=distphim){
1136 toucLad(idetot-1)=ladm;
1137 toucDet(idetot-1)=detinters;
1139 if(rangephi>=distphip){
1141 toucLad(idetot-1)=ladp;
1142 toucDet(idetot-1)=detinters;
1145 ////////////////////////////////////////////////////////////
1147 for (iriv=0; iriv<idetot; iriv++) { //for on detectors
1148 ///////////////////////////////////////////////////////
1149 /*** Rec points sorted by module *****/
1150 /**************************************/
1152 indexmod = g1->GetModuleIndex(lycur,(Int_t)toucLad(iriv),
1153 (Int_t)toucDet(iriv));
1154 fITS->ResetRecPoints();
1155 gAlice->TreeR()->GetEvent(indexmod);
1156 Int_t npoints=frecPoints->GetEntries();
1159 for(indnew=0; indnew<npoints; indnew++){
1160 if (*(fvettid[indexmod]+indnew)==0)
1161 recp =(AliITSRecPoint*)frecPoints->UncheckedAt(indnew);
1164 TVector cluster(3),vecclust(9);
1165 //vecclust(6)=vecclust(7)=vecclust(8)=-1.;
1167 // now vecclust is with cylindrical cohordinates
1168 vecclust(0)=(Float_t)fRecCylR[indexmod][indnew];
1169 vecclust(1)=(Float_t)fRecCylPhi[indexmod][indnew];
1170 vecclust(2)=(Float_t)fRecCylZ[indexmod][indnew];
1171 vecclust(3) = (Float_t)recp->fTracks[0];
1172 vecclust(4) = (Float_t)indnew;
1173 vecclust(5) = (Float_t)indexmod;
1174 vecclust(6) = (Float_t)recp->fTracks[0];
1175 vecclust(7) = (Float_t)recp->fTracks[1];
1176 vecclust(8) = (Float_t)recp->fTracks[2];
1177 sigma[0] = (Double_t) recp->GetSigmaX2();
1178 sigma[1] = (Double_t) recp->GetSigmaZ2();
1180 cluster(0)=fRecCylR[indexmod][indnew];
1181 cluster(1)=fRecCylPhi[indexmod][indnew];
1182 cluster(2)=fRecCylZ[indexmod][indnew];
1184 // cout<<" layer = "<<play<<"\n";
1185 // cout<<" cluster prima = "<<vecclust(0)<<" "
1186 // <<vecclust(1)<<" "
1187 // <<vecclust(2)<<"\n"; getchar();
1189 Float_t sigmatotphi, sigmatotz;
1190 // Float_t epsphi=5.0, epsz=5.0;
1191 //if(fPtref<0.2) {epsphi=3.; epsz=3.;}
1192 Double_t rTrack=(*trackITS).Getrtrack();
1193 Double_t sigmaphi=sigma[0]/(rTrack*rTrack);
1194 sigmatotphi=epsphi*TMath::Sqrt(sigmaphi +
1195 (*trackITS).GetSigmaphi());
1196 sigmatotz=epsz*TMath::Sqrt(sigma[1] +
1197 (*trackITS).GetSigmaZ());
1198 //cout<<"cluster e sigmatotphi e track = "<<cluster(0)
1199 // <<" "<<cluster(1)<<" "<<sigmatotphi<<" "
1200 // <<vecclust(3)<<"\n";
1201 //if(vecclust(3)==481) getchar();
1202 if(cluster(1)<6. && (*trackITS).Getphi()>6.)
1203 cluster(1)=cluster(1)+(2.*TMath::Pi());
1204 if(cluster(1)>6. && (*trackITS).Getphi()<6.)
1205 cluster(1)=cluster(1)-(2.*TMath::Pi());
1206 if(TMath::Abs(cluster(1)-(*trackITS).Getphi())>sigmatotphi)
1208 // cout<<" supero sigmaphi \n";
1209 AliITSTrackV1 *newTrack = new AliITSTrackV1((*trackITS));
1210 (*newTrack).SetLayer((*trackITS).GetLayer()-1);
1211 if (TMath::Abs(rTrack-cluster(0))/rTrack>1e-6)
1212 (*newTrack).Correct(Double_t(cluster(0)));
1213 //cout<<" cluster(2) e(*newTrack).GetZ()="<<cluster(2)<<" "
1214 // << (*newTrack).GetZ()<<"\n";
1215 if(TMath::Abs(cluster(2)-(*newTrack).GetZ()) > sigmatotz){
1219 Double_t sigmanew[2];
1220 sigmanew[0]= sigmaphi;
1221 sigmanew[1]=sigma[1];
1225 // Double_t chi2pred=newTrack->GetPredChi2(m,sigmanew);
1226 // cout<<" chi2pred = "<<chi2pred<<"\n";
1227 // if(chi2pred>fChi2max) continue; //aggiunto il 30-7-2001
1228 if(iriv == 0) flaghit=1;
1229 (*newTrack).AddMS(frl); // add the multiple scattering
1230 //matrix to the covariance matrix
1231 (*newTrack).AddEL(frl,1.,0);
1234 KalmanFilterVert(newTrack,cluster,sigmanew);
1235 //KalmanFilterVert(newTrack,cluster,sigmanew,chi2pred);
1237 KalmanFilter(newTrack,cluster,sigmanew);
1239 (*newTrack).PutCluster(layernew, vecclust);
1240 newTrack->AddClustInTrack();
1241 listoftrack.AddLast(newTrack);
1243 } // end of for on detectors (iriv)
1244 }//end if(outinters==0)
1246 if(flaghit==0 || outinters==-2) {
1247 AliITSTrackV1 *newTrack = new AliITSTrackV1(*trackITS);
1248 (*newTrack).Setfnoclust();
1249 (*newTrack).SetLayer((*trackITS).GetLayer()-1);
1250 (*newTrack).AddMS(frl); // add the multiple scattering matrix
1251 // to the covariance matrix
1252 (*newTrack).AddEL(frl,1.,0);
1253 listoftrack.AddLast(newTrack);
1256 //gObjectTable->Print(); // stampa memoria
1258 RecursiveTracking(&listoftrack);
1259 listoftrack.Delete();
1260 } // end of for on tracks (index)
1262 //gObjectTable->Print(); // stampa memoria
1264 //______________________________________________________________________
1265 Int_t AliITSTrackerV1::Intersection(AliITSTrackV1 &track,Int_t layer,
1266 Int_t &ladder,Int_t &detector) {
1267 // Origin A. Badala' and G.S. Pappalardo
1268 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1269 // Found the intersection and the detector
1271 Double_t rk=fAvrad[layer-1];
1272 if(track.DoNotCross(rk)){ /*cout<< " Do not cross \n";*/ return -1;}
1273 track.Propagation(rk);
1274 Double_t zinters=track.GetZ();
1275 Double_t phinters=track.Getphi();
1276 //cout<<"zinters = "<<zinters<<" phinters = "<<phinters<<"\n";
1280 TVector distZCenter(2);
1284 for(iD = 1; iD<= fNdet[layer-1]; iD++) {
1285 if(zinters > fzmin[layer-1][iD-1] && zinters <= fzmax[layer-1][iD-1]) {
1287 cout<< " Errore su iz in Intersection \n";
1290 listDet(iz)= iD; distZCenter(iz)=TMath::Abs(zinters-det(2));
1296 if(iz==0) {/* cout<< " No detector along Z \n";*/ return -2;}
1297 detector=Int_t (listDet(0));
1298 if(iz>1 && (distZCenter(0)>distZCenter(1))) detector=Int_t (listDet(1));
1301 TVector distPhiCenter(2);
1303 Double_t pigre=TMath::Pi();
1305 for(iLd = 1; iLd<= fNlad[layer-1]; iLd++) {
1306 Double_t phimin=fphimin[layer-1][iLd-1];
1307 Double_t phimax=fphimax[layer-1][iLd-1];
1308 Double_t phidet=fphidet[layer-1][iLd-1];
1309 Double_t phiconfr=phinters;
1311 //if(phimin <5.5) {cout<<" Error in Intersection for phi \n";
1314 if(phinters>(1.5*pigre)) phiconfr=phinters-(2.*pigre);
1315 if(phidet>(1.5*pigre)) phidet-=(2.*pigre);
1317 if(phiconfr>phimin && phiconfr<= phimax) {
1319 cout<< " Errore su ip in Intersection \n"; getchar();
1322 distPhiCenter(ip)=TMath::Abs(phiconfr-phidet); ip++;
1326 if(ip==0) { cout<< " No detector along phi \n"; getchar();}
1327 ladder=Int_t (listLad(0));
1328 if(ip>1 && (distPhiCenter(0)>distPhiCenter(1))) ladder=Int_t (listLad(1));
1331 //______________________________________________________________________
1332 void AliITSTrackerV1::KalmanFilter(AliITSTrackV1 *newTrack,TVector &cluster,
1334 //Origin A. Badala' and G.S. Pappalardo:
1335 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1336 // Kalman filter without vertex constraint
1337 ////// Evaluation of the measurement vector ////////////////////////
1339 Double_t rk,phik,zk;
1340 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1342 //////////////////////// Evaluation of the error matrix V /////////
1343 Double_t v00=sigma[0];
1344 Double_t v11=sigma[1];
1345 ////////////////////////////////////////////////////////////////////
1346 Double_t cin00,cin10,cin20,cin30,cin40,cin11,cin21,cin31,cin41,cin22,
1347 cin32,cin42,cin33,cin43,cin44;
1349 newTrack->GetCElements(cin00,
1352 cin30,cin31,cin32,cin33,
1353 cin40,cin41,cin42,cin43,cin44); //get C matrix
1354 Double_t rold00=cin00+v00;
1355 Double_t rold10=cin10;
1356 Double_t rold11=cin11+v11;
1357 ////////////////////// R matrix inversion /////////////////////////
1358 Double_t det=rold00*rold11-rold10*rold10;
1359 Double_t r00=rold11/det;
1360 Double_t r10=-rold10/det;
1361 Double_t r11=rold00/det;
1362 ////////////////////////////////////////////////////////////////////
1363 Double_t k00=cin00*r00+cin10*r10;
1364 Double_t k01=cin00*r10+cin10*r11;
1365 Double_t k10=cin10*r00+cin11*r10;
1366 Double_t k11=cin10*r10+cin11*r11;
1367 Double_t k20=cin20*r00+cin21*r10;
1368 Double_t k21=cin20*r10+cin21*r11;
1369 Double_t k30=cin30*r00+cin31*r10;
1370 Double_t k31=cin30*r10+cin31*r11;
1371 Double_t k40=cin40*r00+cin41*r10;
1372 Double_t k41=cin40*r10+cin41*r11;
1373 Double_t x0,x1,x2,x3,x4;
1374 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1375 Double_t savex0=x0, savex1=x1;
1376 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1);
1377 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1);
1378 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1);
1379 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1);
1380 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1);
1381 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1382 c00=cin00-k00*cin00-k01*cin10;
1383 c10=cin10-k00*cin10-k01*cin11;
1384 c20=cin20-k00*cin20-k01*cin21;
1385 c30=cin30-k00*cin30-k01*cin31;
1386 c40=cin40-k00*cin40-k01*cin41;
1387 c11=cin11-k10*cin10-k11*cin11;
1388 c21=cin21-k10*cin20-k11*cin21;
1389 c31=cin31-k10*cin30-k11*cin31;
1390 c41=cin41-k10*cin40-k11*cin41;
1391 c22=cin22-k20*cin20-k21*cin21;
1392 c32=cin32-k20*cin30-k21*cin31;
1393 c42=cin42-k20*cin40-k21*cin41;
1394 c33=cin33-k30*cin30-k31*cin31;
1395 c43=cin43-k30*cin40-k31*cin41;
1396 c44=cin44-k40*cin40-k41*cin41;
1397 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1398 newTrack->PutCElements(c00,
1402 c40,c41,c42,c43,c44); // put in track the
1404 Double_t vmcold00=v00-c00;
1405 Double_t vmcold10=-c10;
1406 Double_t vmcold11=v11-c11;
1407 ////////////////////// Matrix vmc inversion ///////////////////////
1408 det=vmcold00*vmcold11-vmcold10*vmcold10;
1409 Double_t vmc00=vmcold11/det;
1410 Double_t vmc10=-vmcold10/det;
1411 Double_t vmc11=vmcold00/det;
1412 ////////////////////////////////////////////////////////////////////
1413 Double_t chi2=(m[0]-x0)*( vmc00*(m[0]-x0) + 2.*vmc10*(m[1]-x1) ) +
1414 (m[1]-x1)*vmc11*(m[1]-x1);
1415 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1417 //----------------------------------------------------------------------
1418 //void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1419 // TVector &cluster,Double_t sigma[2]){
1420 void AliITSTrackerV1::KalmanFilterVert(AliITSTrackV1 *newTrack,
1421 TVector &cluster,Double_t sigma[2]
1422 /*, Double_t chi2pred*/){
1423 //Origin A. Badala' and G.S. Pappalardo:
1424 // e-mail Angela.Badala@ct.infn.it, Giuseppe.S.Pappalardo@ct.infn.it
1425 // Kalman filter with vertex constraint
1426 ///////////////////// Evaluation of the measurement vector m ////////
1428 Double_t rk,phik,zk;
1429 rk=cluster(0); phik=cluster(1); zk=cluster(2);
1431 Double_t cc=(*newTrack).GetC();
1432 Double_t zv=(*newTrack).GetZv();
1433 Double_t dv=(*newTrack).GetDv();
1435 Double_t tgl= (zk-zv)*cy/TMath::ASin(cy*rk);
1437 /////////////////////// Evaluation of the error matrix V //////////
1438 Int_t layer=newTrack->GetLayer();
1439 Double_t v00=sigma[0];
1440 Double_t v11=sigma[1];
1441 Double_t v31=sigma[1]/rk;
1442 Double_t sigmaDv=newTrack->GetsigmaDv();
1443 Double_t v22=sigmaDv*sigmaDv + newTrack->Getd2(layer-1);
1444 Double_t v32=newTrack->Getdtgl(layer-1);
1445 Double_t sigmaZv=newTrack->GetsigmaZv();
1446 Double_t v33=(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+newTrack->Gettgl2(layer-1);
1447 //////////////////////////////////////////////////////////////////
1448 Double_t cin00,cin10,cin11,cin20,cin21,cin22,
1449 cin30,cin31,cin32,cin33,cin40,cin41,cin42,cin43,cin44;
1450 newTrack->GetCElements(cin00,
1453 cin30,cin31,cin32,cin33,
1454 cin40,cin41,cin42,cin43,cin44); //get C matrix
1462 r[3][1]=cin31+sigma[1]/rk;
1463 r[2][2]=cin22+sigmaDv*sigmaDv+newTrack->Getd2(layer-1);
1464 r[3][2]=cin32+newTrack->Getdtgl(layer-1);
1465 r[3][3]=cin33+(sigma[1]+sigmaZv*sigmaZv)/(rk*rk)+
1466 newTrack->Gettgl2(layer-1);
1467 r[0][1]=r[1][0]; r[0][2]=r[2][0]; r[0][3]=r[3][0];
1468 r[1][2]=r[2][1]; r[1][3]=r[3][1]; r[2][3]=r[3][2];
1469 ///////////////////// Matrix R inversion //////////////////////////
1473 Int_t ll[kn],mm[kn];
1476 for(k=0; k<kn; k++) {
1480 for(j=k; j<kn ; j++) {
1481 for (i=j; i<kn; i++) {
1482 if(TMath::Abs(big) < TMath::Abs(r[i][j]) ) {
1492 for(i=0; i<kn; i++) {
1501 for(j=0; j<kn; j++) {
1510 cout << "Singular matrix\n";
1512 for(i=0; i<kn; i++) {
1513 if(i == k) { continue; }
1514 r[i][k]=r[i][k]/(-big);
1517 for(i=0; i<kn; i++) {
1519 for(j=0; j<kn; j++) {
1520 if(i == k || j == k) continue;
1521 r[i][j]=hold*r[k][j]+r[i][j];
1525 for(j=0; j<kn; j++) {
1526 if(j == k) continue;
1527 r[k][j]=r[k][j]/big;
1535 for(k=kn-1; k>=0; k--) {
1538 for (j=0; j<kn; j++) {
1546 for (i=0; i<kn; i++) {
1553 ////////////////////////////////////////////////////////////////////
1554 Double_t k00=cin00*r[0][0]+cin10*r[1][0]+cin20*r[2][0]+cin30*r[3][0];
1555 Double_t k01=cin00*r[1][0]+cin10*r[1][1]+cin20*r[2][1]+cin30*r[3][1];
1556 Double_t k02=cin00*r[2][0]+cin10*r[2][1]+cin20*r[2][2]+cin30*r[3][2];
1557 Double_t k03=cin00*r[3][0]+cin10*r[3][1]+cin20*r[3][2]+cin30*r[3][3];
1558 Double_t k10=cin10*r[0][0]+cin11*r[1][0]+cin21*r[2][0]+cin31*r[3][0];
1559 Double_t k11=cin10*r[1][0]+cin11*r[1][1]+cin21*r[2][1]+cin31*r[3][1];
1560 Double_t k12=cin10*r[2][0]+cin11*r[2][1]+cin21*r[2][2]+cin31*r[3][2];
1561 Double_t k13=cin10*r[3][0]+cin11*r[3][1]+cin21*r[3][2]+cin31*r[3][3];
1562 Double_t k20=cin20*r[0][0]+cin21*r[1][0]+cin22*r[2][0]+cin32*r[3][0];
1563 Double_t k21=cin20*r[1][0]+cin21*r[1][1]+cin22*r[2][1]+cin32*r[3][1];
1564 Double_t k22=cin20*r[2][0]+cin21*r[2][1]+cin22*r[2][2]+cin32*r[3][2];
1565 Double_t k23=cin20*r[3][0]+cin21*r[3][1]+cin22*r[3][2]+cin32*r[3][3];
1566 Double_t k30=cin30*r[0][0]+cin31*r[1][0]+cin32*r[2][0]+cin33*r[3][0];
1567 Double_t k31=cin30*r[1][0]+cin31*r[1][1]+cin32*r[2][1]+cin33*r[3][1];
1568 Double_t k32=cin30*r[2][0]+cin31*r[2][1]+cin32*r[2][2]+cin33*r[3][2];
1569 Double_t k33=cin30*r[3][0]+cin31*r[3][1]+cin32*r[3][2]+cin33*r[3][3];
1570 Double_t k40=cin40*r[0][0]+cin41*r[1][0]+cin42*r[2][0]+cin43*r[3][0];
1571 Double_t k41=cin40*r[1][0]+cin41*r[1][1]+cin42*r[2][1]+cin43*r[3][1];
1572 Double_t k42=cin40*r[2][0]+cin41*r[2][1]+cin42*r[2][2]+cin43*r[3][2];
1573 Double_t k43=cin40*r[3][0]+cin41*r[3][1]+cin42*r[3][2]+cin43*r[3][3];
1575 Double_t x0,x1,x2,x3,x4;
1576 newTrack->GetXElements(x0,x1,x2,x3,x4); // get the state vector
1577 Double_t savex0=x0, savex1=x1, savex2=x2, savex3=x3;
1578 x0+=k00*(m[0]-savex0)+k01*(m[1]-savex1)+k02*(m[2]-savex2)+
1580 x1+=k10*(m[0]-savex0)+k11*(m[1]-savex1)+k12*(m[2]-savex2)+
1582 x2+=k20*(m[0]-savex0)+k21*(m[1]-savex1)+k22*(m[2]-savex2)+
1584 x3+=k30*(m[0]-savex0)+k31*(m[1]-savex1)+k32*(m[2]-savex2)+
1586 x4+=k40*(m[0]-savex0)+k41*(m[1]-savex1)+k42*(m[2]-savex2)+
1588 Double_t c00,c10,c20,c30,c40,c11,c21,c31,c41,c22,c32,c42,c33,c43,c44;
1589 c00=cin00-k00*cin00-k01*cin10-k02*cin20-k03*cin30;
1590 c10=cin10-k00*cin10-k01*cin11-k02*cin21-k03*cin31;
1591 c20=cin20-k00*cin20-k01*cin21-k02*cin22-k03*cin32;
1592 c30=cin30-k00*cin30-k01*cin31-k02*cin32-k03*cin33;
1593 c40=cin40-k00*cin40-k01*cin41-k02*cin42-k03*cin43;
1594 c11=cin11-k10*cin10-k11*cin11-k12*cin21-k13*cin31;
1595 c21=cin21-k10*cin20-k11*cin21-k12*cin22-k13*cin32;
1596 c31=cin31-k10*cin30-k11*cin31-k12*cin32-k13*cin33;
1597 c41=cin41-k10*cin40-k11*cin41-k12*cin42-k13*cin43;
1598 c22=cin22-k20*cin20-k21*cin21-k22*cin22-k23*cin32;
1599 c32=cin32-k20*cin30-k21*cin31-k22*cin32-k23*cin33;
1600 c42=cin42-k20*cin40-k21*cin41-k22*cin42-k23*cin43;
1601 c33=cin33-k30*cin30-k31*cin31-k32*cin32-k33*cin33;
1602 c43=cin43-k30*cin40-k31*cin41-k32*cin42-k33*cin43;
1603 c44=cin44-k40*cin40-k41*cin41-k42*cin42-k43*cin43;
1605 newTrack->PutXElements(x0,x1,x2,x3,x4); // put the new state vector
1606 newTrack->PutCElements(c00,
1610 c40,c41,c42,c43,c44); // put in track the
1613 vmc[0][0]=v00-c00; vmc[1][0]=-c10; vmc[2][0]=-c20; vmc[3][0]=-c30;
1614 vmc[1][1]=v11-c11; vmc[2][1]=-c21; vmc[3][1]=v31-c31;
1615 vmc[2][2]=v22-c22; vmc[3][2]=v32-c32;
1617 vmc[0][1]=vmc[1][0]; vmc[0][2]=vmc[2][0]; vmc[0][3]=vmc[3][0];
1618 vmc[1][2]=vmc[2][1]; vmc[1][3]=vmc[3][1];
1619 vmc[2][3]=vmc[3][2];
1620 /////////////////////// vmc matrix inversion ///////////////////////
1622 for(k=0; k<kn; k++) {
1626 for(j=k; j<kn ; j++) {
1627 for (i=j; i<kn; i++) {
1628 if(TMath::Abs(big) < TMath::Abs(vmc[i][j]) ) {
1638 for(i=0; i<kn; i++) {
1640 vmc[k][i]=vmc[j][i];
1647 for(j=0; j<kn; j++) {
1649 vmc[j][k]=vmc[j][i];
1656 cout << "Singular matrix\n";
1658 for(i=0; i<kn; i++) {
1659 if(i == k) continue;
1660 vmc[i][k]=vmc[i][k]/(-big);
1663 for(i=0; i<kn; i++) {
1665 for(j=0; j<kn; j++) {
1666 if(i == k || j == k) continue;
1667 vmc[i][j]=hold*vmc[k][j]+vmc[i][j];
1671 for(j=0; j<kn; j++) {
1672 if(j == k) continue;
1673 vmc[k][j]=vmc[k][j]/big;
1681 for(k=kn-1; k>=0; k--) {
1684 for (j=0; j<kn; j++) {
1686 vmc[j][k]=-vmc[j][i];
1692 for (i=0; i<kn; i++) {
1694 vmc[k][i]=-vmc[j][i];
1699 ////////////////////////////////////////////////////////////////////
1700 Double_t chi2=(m[0]-x0)*( vmc[0][0]*(m[0]-x0) + 2.*vmc[1][0]*(m[1]-x1) +
1701 2.*vmc[2][0]*(m[2]-x2)+ 2.*vmc[3][0]*(m[3]-x3) )+
1702 (m[1]-x1)* ( vmc[1][1]*(m[1]-x1) + 2.*vmc[2][1]*(m[2]-x2)+
1703 2.*vmc[3][1]*(m[3]-x3) ) +
1704 (m[2]-x2)* ( vmc[2][2]*(m[2]-x2)+ 2.*vmc[3][2]*(m[3]-x3) ) +
1705 (m[3]-x3)*vmc[3][3]*(m[3]-x3);
1706 //cout<<" chi2 kalman = "<<chi2<<"\n"; getchar();
1707 newTrack->SetChi2(newTrack->GetChi2()+chi2);
1708 // newTrack->SetChi2(newTrack->GetChi2()+chi2pred);